source: trunk/LMDZ.MARS/util/expandstartfi.F90 @ 2739

Last change on this file since 2739 was 1469, checked in by tblmd, 9 years ago

TB: bug fixed in expandstartfi at south pole for subsurface variables

File size: 17.8 KB
Line 
1program expandstartfi
2
3! program which takes a GCM startfi.nc file and puts it on the
4! corresponding lonxlat grid (so it can be displayed using Grads, Ferret, etc.)
5! usage:
6! expandstartfi  [infile.nc] [outfile.nc]
7!     if outfile is not specified it is built as "infile_ex.nc"
8!     if infile is not specified, "startfi.nc" is used as default
9
10use netcdf
11implicit none
12
13! Input and output files:
14character(len=256) :: infile="startfi.nc"      ! default input file
15character(len=256) :: outfile="startfi_ex.nc"  ! default output file
16
17! NetCDF stuff
18integer :: status ! NetCDF return code
19integer :: inid,outid ! NetCDF file IDs
20integer :: varid ! to store the ID of a variable
21integer :: datashape(4) ! to store dimension IDs of a given dataset
22integer :: corner(3),edges(3) ! to read data with a time axis
23character(len=90) :: varname ! name of a variable
24character(len=90) :: varatt ! name of attribute of a variable
25character(len=90) :: varattcontent ! content of the attribute
26! Input file dimension IDs:
27integer :: physical_points_dimid
28integer :: subsurface_layers_dimid
29integer :: nlayer_plus_1_dimid
30integer :: number_of_advected_fields_dimid
31integer :: time_dimid
32integer :: nbindims ! number of dimensions in input file
33integer :: nbinvars ! number of variables in input file
34integer :: invarid ! to store ID of an input variable
35!Input file variables
36integer :: physical_points
37integer :: subsurface_layers
38integer :: nlayer_plus_1
39integer :: number_of_advected_fields
40integer :: timelen
41real,allocatable :: surf_field(:) ! to store a 1D field of physical_points elements
42real,allocatable :: subsurf_field(:,:) ! to store subsurface (2D field)
43
44! Output file dimensions:
45integer :: latlen
46integer :: lonlen
47! Output file variables:
48real,allocatable :: longitude(:)
49real,allocatable :: latitude(:)
50real,allocatable :: depth(:)
51real,allocatable :: out_surf_field(:,:)
52real,allocatable :: out_subsurf_field(:,:,:)
53! Output file dimension IDs
54integer :: lon_dimid
55integer :: lat_dimid
56integer :: depth_dimid
57! IDs of Output file variables
58integer :: lon_varid
59integer :: lat_varid
60integer :: depth_varid
61
62integer :: i,j,k,ig0,ivar
63integer :: nbdim,nbatt,shape(4)
64integer :: nbarg ! # of program arguments
65character(len=256) :: arg ! to store a program argument
66real :: pi ! 3.14159...
67
68pi=2.*asin(1.)
69
70! 0. Preliminary stuff, check arguments (input and output files)
71! get number of arguments this program was called with
72nbarg=command_argument_count()
73
74if (nbarg.ge.1) then
75  call get_command_argument(1,arg) ! get argument # 1
76  infile=trim(arg)
77  if (nbarg.eq.2) then
78    call get_command_argument(2,arg) ! get argument # 2
79    outfile=trim(arg)
80  else
81   ! build outfile from infile (replace ".nc" with "_ex.nc"
82   outfile=trim(infile(1:index(infile,".nc",back=.true.)-1))//"_ex.nc"
83  endif
84  if (nbarg.ge.3) then
85    write(*,*) ' Warning: Too many arguments...'
86    write(*,*) '         will only use the first 2 '
87  endif
88endif
89
90! 1. open input file
91status=nf90_open(infile,NF90_NOWRITE,inid)
92if (status.ne.nf90_noerr) then
93  write(*,*)"Failed to open datafile ",trim(infile)
94  write(*,*)trim(nf90_strerror(status))
95 stop
96endif
97write(*,*) "Reading input file: ",trim(infile)
98
99! 1.2 Identify/load dimensions in input file
100status=nf90_inq_dimid(inid,"physical_points",physical_points_dimid)
101if (status.ne.nf90_noerr) then
102  write(*,*)"Failed to find physical_points dimension"
103  write(*,*)trim(nf90_strerror(status))
104 stop
105endif
106status=nf90_inquire_dimension(inid,physical_points_dimid,len=physical_points)
107if (status.ne.nf90_noerr) then
108  write(*,*)"Failed to find physical_points value"
109  write(*,*)trim(nf90_strerror(status))
110 stop
111else
112  write(*,*) " physical_points = ",physical_points
113endif
114
115status=nf90_inq_dimid(inid,"subsurface_layers",subsurface_layers_dimid)
116if (status.ne.nf90_noerr) then
117  write(*,*)"Failed to find subsurface_layers dimension"
118  write(*,*)trim(nf90_strerror(status))
119 stop
120endif
121status=nf90_inquire_dimension(inid,subsurface_layers_dimid,len=subsurface_layers)
122if (status.ne.nf90_noerr) then
123  write(*,*)"Failed to find subsurface_layers value"
124  write(*,*)trim(nf90_strerror(status))
125 stop
126else
127  write(*,*) " subsurface_layers = ",subsurface_layers
128endif
129
130status=nf90_inq_dimid(inid,"nlayer_plus_1",nlayer_plus_1_dimid)
131if (status.ne.nf90_noerr) then
132  write(*,*)"Failed to find nlayer_plus_1 dimension"
133  write(*,*)trim(nf90_strerror(status))
134 stop
135endif
136status=nf90_inquire_dimension(inid,nlayer_plus_1_dimid,len=nlayer_plus_1)
137if (status.ne.nf90_noerr) then
138  write(*,*)"Failed to find nlayer_plus_1 value"
139  write(*,*)trim(nf90_strerror(status))
140 stop
141else
142  write(*,*) " nlayer_plus_1 = ",nlayer_plus_1
143endif
144
145status=nf90_inq_dimid(inid,"number_of_advected_fields",number_of_advected_fields_dimid)
146if (status.ne.nf90_noerr) then
147  write(*,*)"Failed to find number_of_advected_fields dimension"
148  write(*,*)trim(nf90_strerror(status))
149 stop
150endif
151status=nf90_inquire_dimension(inid,number_of_advected_fields_dimid,len=number_of_advected_fields)
152if (status.ne.nf90_noerr) then
153  write(*,*)"Failed to find number_of_advected_fields value"
154  write(*,*)trim(nf90_strerror(status))
155 stop
156else
157  write(*,*) " number_of_advected_fields = ",number_of_advected_fields
158endif
159
160status=nf90_inq_dimid(inid,"Time",time_dimid)
161if (status.ne.nf90_noerr) then
162  write(*,*)"Failed to find Time dimension"
163  write(*,*)trim(nf90_strerror(status))
164  timelen = 0
165else
166  status=nf90_inquire_dimension(inid,time_dimid,len=timelen)
167  if (status.ne.nf90_noerr) then
168    write(*,*)"Failed to read Time dimension"
169    write(*,*)trim(nf90_strerror(status))
170   stop
171  else
172    write(*,*) " time length = ",timelen
173  endif
174endif
175
176! 1.3 Allocate memory for input fields
177allocate(surf_field(physical_points))
178allocate(subsurf_field(physical_points,subsurface_layers))
179
180! 2.1. Create output file
181status=nf90_create(outfile,NF90_CLOBBER,outid)
182if (status.ne.nf90_noerr) then
183  write(*,*) "Failed to create output file: ",trim(outfile)
184  write(*,*)trim(nf90_strerror(status))
185  stop
186endif
187write(*,*) " "
188write(*,*) "Created output file: ",trim(outfile)
189
190! 2.2. Build dimensions for output file
191
192! 2.2.1 Compute lonlen and latlen from  physical_points
193! load "longitude(physical_points)" from input file
194status=nf90_inq_varid(inid,"longitude",varid)
195if (status.ne.nf90_noerr) then
196  write(*,*) "Failed to find longitude ID"
197  write(*,*)trim(nf90_strerror(status))
198  stop
199endif
200! read longitude
201status=nf90_get_var(inid,varid,surf_field)
202if (status.ne.nf90_noerr) then
203  write(*,*) "Failed to load longitude"
204  write(*,*)trim(nf90_strerror(status))
205  stop
206endif
207
208! count the number of points before longitude(i) wraps around
209i=3
210lonlen=1
211!write(*,*) "longitude(2)=",surf_field(2)
212do while (surf_field(i).ne.surf_field(2))
213!write(*,*) "i:",i,"surf_field(i)=",surf_field(i)
214 i=i+1
215 lonlen=lonlen+1
216enddo
217! and add 1 because we want a redundant lon=180 point
218lonlen=lonlen+1
219write(*,*) " => lonlen=",lonlen
220
221allocate(longitude(lonlen))
222! fill longitude(1:lonlen)
223longitude(1:lonlen-1)=surf_field(2:lonlen)
224longitude(lonlen)=-longitude(1) ! redundant +Pi/2
225! convert to degrees
226longitude(:)=longitude(:)*(180./pi)
227
228! knowing lonlen, compute latlen
229latlen=2+(physical_points-2)/(lonlen-1)
230write(*,*) " => latlen=",latlen
231
232allocate(latitude(latlen))
233! get field "latitude(physical_points)" from infile
234status=nf90_inq_varid(inid,"latitude",varid)
235if (status.ne.nf90_noerr) then
236  write(*,*) "Failed to find latitude ID"
237  write(*,*)trim(nf90_strerror(status))
238  stop
239endif
240! read latitude
241status=nf90_get_var(inid,varid,surf_field)
242if (status.ne.nf90_noerr) then
243  write(*,*) "Failed to load latitude"
244  write(*,*)trim(nf90_strerror(status))
245  stop
246endif
247
248latitude(1)=surf_field(1)
249!write(*,*) "latitude(1)=",latitude(1)
250do i=2,latlen-1
251  latitude(i)=surf_field((i-1)*(lonlen-1))
252!  write(*,*) "i:",i,"latitude(i)=",latitude(i)
253enddo
254latitude(latlen)=surf_field(physical_points)
255!write(*,*) "latitude(latlen)=",latitude(latlen)
256! convert to degrees
257latitude(:)=latitude(:)*(180./pi)
258
259! depth:
260allocate(depth(subsurface_layers))
261! load "soildepth(subsurface_layers)" from input file
262status=nf90_inq_varid(inid,"soildepth",varid)
263if (status.ne.nf90_noerr) then
264  write(*,*) "Failed to find soildepth ID"
265  write(*,*)trim(nf90_strerror(status))
266  stop
267endif
268! read longitude
269status=nf90_get_var(inid,varid,depth)
270if (status.ne.nf90_noerr) then
271  write(*,*) "Failed to load soildepth"
272  write(*,*)trim(nf90_strerror(status))
273  stop
274endif
275
276! 2.2.2 define dimensions to output file
277! longitude:
278status=nf90_def_dim(outid,"longitude",lonlen,lon_dimid)
279if (status.ne.nf90_noerr) then
280  write(*,*) "Failed creating longitude dimension"
281  write(*,*)trim(nf90_strerror(status))
282  stop
283endif
284! latitude:
285status=nf90_def_dim(outid,"latitude",latlen,lat_dimid)
286if (status.ne.nf90_noerr) then
287  write(*,*) "Failed creating longitude dimension"
288  write(*,*)trim(nf90_strerror(status))
289  stop
290endif
291! depth:
292status=nf90_def_dim(outid,"depth",subsurface_layers,depth_dimid)
293if (status.ne.nf90_noerr) then
294  write(*,*) "Failed creating depth dimension"
295  write(*,*)trim(nf90_strerror(status))
296  stop
297endif
298
299!2.3. write variables to output file
300! 2.3.1. define 1D variables
301! longitude
302datashape(1)=lon_dimid
303status=nf90_def_var(outid,"longitude",NF90_REAL,lon_dimid,lon_varid)
304if (status.ne.nf90_noerr) then
305  write(*,*) "Failed creating longitude variable"
306  write(*,*)trim(nf90_strerror(status))
307  stop
308endif
309! longitude attributes
310status=nf90_put_att(outid,lon_varid,"long_name","East longitude")
311status=nf90_put_att(outid,lon_varid,"units","degrees_east")
312
313!latitude
314datashape(2)=lat_dimid
315status=nf90_def_var(outid,"latitude",NF90_REAL,lat_dimid,lat_varid)
316if (status.ne.nf90_noerr) then
317  write(*,*) "Failed creating latitude variable"
318  write(*,*)trim(nf90_strerror(status))
319  stop
320endif
321! latitude attributes
322status=nf90_put_att(outid,lat_varid,"long_name","North latitude")
323status=nf90_put_att(outid,lat_varid,"units","degrees_north")
324
325! depth
326status=nf90_def_var(outid,"depth",NF90_REAL,depth_dimid,depth_varid)
327if (status.ne.nf90_noerr) then
328  write(*,*) "Failed creating depth variable"
329  write(*,*)trim(nf90_strerror(status))
330  stop
331endif
332!depth atributes
333status=nf90_put_att(outid,depth_varid,"long_name","Soil mid-layer depth")
334status=nf90_put_att(outid,depth_varid,"units","m")
335status=nf90_put_att(outid,depth_varid,"positive","down")
336
337! 2.3.2 write 1D variable
338
339! swich out of NetCDF define mode
340status=nf90_enddef(outid)
341if (status.ne.nf90_noerr) then
342  write(*,*) "Failed to swich out of define mode"
343  write(*,*)trim(nf90_strerror(status))
344  stop
345endif
346
347! write longitude
348status=nf90_put_var(outid,lon_varid,longitude)
349if (status.ne.nf90_noerr) then
350  write(*,*) "Failed writing longitude"
351  write(*,*)trim(nf90_strerror(status))
352  stop
353endif
354
355! write latitude
356status=nf90_put_var(outid,lat_varid,latitude)
357if (status.ne.nf90_noerr) then
358  write(*,*) "Failed writing latitude"
359  write(*,*)trim(nf90_strerror(status))
360  stop
361endif
362
363! write depth
364status=nf90_put_var(outid,depth_varid,depth)
365if (status.ne.nf90_noerr) then
366  write(*,*) "Failed writing depth"
367  write(*,*)trim(nf90_strerror(status))
368  stop
369endif
370
371! 2.3. 2D (surface) variables
372! First find out how many variables there are in input file
373status=nf90_inquire(inid,nbindims,nbinvars)
374if (status.ne.nf90_noerr) then
375  write(*,*) "Failed obtaining nbindims and nbinvars from input file"
376  write(*,*)trim(nf90_strerror(status))
377  stop
378endif
379
380allocate(out_surf_field(lonlen,latlen))
381
382shape(:) = 0
383do ivar=1,nbinvars ! loop on all input variables
384  ! find out what dimensions are linked to this variable
385  status=nf90_inquire_variable(inid,ivar,name=varname,ndims=nbdim,&
386                               dimids=shape,natts=nbatt)
387  if (((nbdim==1).and.(shape(1)==physical_points_dimid))&
388  .or.((nbdim==2).and.(shape(1)==physical_points_dimid)&
389                 .and.(shape(2)==time_dimid))) then
390 
391    corner(1) = 1
392    corner(2) = timelen
393    edges(1)  = physical_points
394    edges(2)  = 1
395   
396    ! skip "longitude" and "latitude"
397    if (trim(varname)=="longitude") cycle
398    if (trim(varname)=="latitude") cycle
399   
400    write(*,*) " processing: ",trim(varname)
401
402    ! load input data:
403    status=nf90_inq_varid(inid,varname,invarid)
404    status=nf90_get_var(inid,invarid,surf_field,corner,edges)
405   
406    ! switch output file to to define mode
407    status=nf90_redef(outid)
408    if (status.ne.nf90_noerr) then
409      write(*,*) "Failed to swich to define mode"
410      write(*,*)trim(nf90_strerror(status))
411      stop
412    endif
413    !define the variable
414    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
415                         (/lon_dimid,lat_dimid/),varid)
416    if (status.ne.nf90_noerr) then
417      write(*,*) "Failed creating variable ",trim(varname)
418      write(*,*)trim(nf90_strerror(status))
419      stop
420    endif
421
422    ! variable attributes
423    do k=1,nbatt
424      status=nf90_inq_attname(inid,invarid,k,varatt)
425      if (status.ne.nf90_noerr) then
426        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
427        write(*,*)trim(nf90_strerror(status))
428      stop
429      endif
430      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
431      if (status.ne.nf90_noerr) then
432        write(*,*) "Failed loading attribute ",trim(varatt)
433        write(*,*)trim(nf90_strerror(status))
434      stop
435      endif
436
437      ! write the attribute and its value to output
438      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
439      if (status.ne.nf90_noerr) then
440        write(*,*) "Failed writing attribute ",trim(varatt)
441        write(*,*)trim(nf90_strerror(status))
442      stop
443      endif
444    enddo ! do k=1,nbatt
445
446    ! swich out of NetCDF define mode
447    status=nf90_enddef(outid)
448    if (status.ne.nf90_noerr) then
449      write(*,*) "Failed to swich out of define mode"
450      write(*,*)trim(nf90_strerror(status))
451      stop
452    endif
453   
454    ! "convert" from surf_field(:) to out_surf_field(:,:)
455    do i=1,lonlen
456      out_surf_field(i,1)=surf_field(1) ! North Pole
457      out_surf_field(i,latlen)=surf_field(physical_points) ! South Pole
458    enddo
459    do j=2,latlen-1
460      ig0=1+(j-2)*(lonlen-1)
461      do i=1,lonlen-1
462        out_surf_field(i,j)=surf_field(ig0+i)
463      enddo
464      out_surf_field(lonlen,j)=out_surf_field(1,j) ! redundant lon=180
465    enddo
466   
467    ! write the variable
468    status=nf90_put_var(outid,varid,out_surf_field)
469    if (status.ne.nf90_noerr) then
470      write(*,*) "Failed to write ",trim(varname)
471      write(*,*)trim(nf90_strerror(status))
472      stop
473    endif
474  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid))
475enddo ! of do i=1,nbinvars
476
477! 2.4. 3D (subsurface) variables
478allocate(out_subsurf_field(lonlen,latlen,subsurface_layers))
479
480do ivar=1,nbinvars ! loop on all input variables
481  ! find out what dimensions are linked to this variable
482  status=nf90_inquire_variable(inid,ivar,name=varname,ndims=nbdim,&
483                               dimids=shape,natts=nbatt)
484  if (((nbdim==2).and.(shape(1)==physical_points_dimid)&
485                 .and.(shape(2)==subsurface_layers_dimid))&
486  .or.((nbdim==3).and.(shape(1)==physical_points_dimid)&
487                 .and.(shape(2)==subsurface_layers_dimid)&
488                 .and.(shape(3)==time_dimid))) then
489   
490    corner(1) = 1
491    corner(2) = 1
492    corner(3) = timelen
493    edges(1)  = physical_points
494    edges(2)  = subsurface_layers
495    edges(3)  = 1
496   
497    write(*,*) " processing: ",trim(varname)
498
499    ! load input data:
500    status=nf90_inq_varid(inid,varname,invarid)
501    status=nf90_get_var(inid,invarid,subsurf_field,corner,edges)
502   
503    ! switch output file to to define mode
504    status=nf90_redef(outid)
505    if (status.ne.nf90_noerr) then
506      write(*,*) "Failed to swich to define mode"
507      write(*,*)trim(nf90_strerror(status))
508      stop
509    endif
510    !define the variable
511    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
512                         (/lon_dimid,lat_dimid,depth_dimid/),varid)
513    if (status.ne.nf90_noerr) then
514      write(*,*) "Failed creating variable ",trim(varname)
515      write(*,*)trim(nf90_strerror(status))
516      stop
517    endif
518
519    ! variable attributes
520    do k=1,nbatt
521      status=nf90_inq_attname(inid,invarid,k,varatt)
522      if (status.ne.nf90_noerr) then
523        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
524        write(*,*)trim(nf90_strerror(status))
525      stop
526      endif
527      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
528      if (status.ne.nf90_noerr) then
529        write(*,*) "Failed loading attribute ",trim(varatt)
530        write(*,*)trim(nf90_strerror(status))
531      stop
532      endif
533
534      ! write the attribute and its value to output
535      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
536      if (status.ne.nf90_noerr) then
537        write(*,*) "Failed writing attribute ",trim(varatt)
538        write(*,*)trim(nf90_strerror(status))
539      stop
540      endif
541    enddo ! do k=1,nbatt
542
543    ! swich out of NetCDF define mode
544    status=nf90_enddef(outid)
545    if (status.ne.nf90_noerr) then
546      write(*,*) "Failed to swich out of define mode"
547      write(*,*)trim(nf90_strerror(status))
548      stop
549    endif
550   
551    ! "convert" from subsurf_field(:,:) to out_subsurf_field(:,:,:)
552    do i=1,lonlen
553      out_subsurf_field(i,1,:)=subsurf_field(1,:) ! North Pole
554      out_subsurf_field(i,latlen,:)=subsurf_field(physical_points,:) ! South Pole
555    enddo
556    do j=2,latlen-1
557      ig0=1+(j-2)*(lonlen-1)
558      do i=1,lonlen-1
559        out_subsurf_field(i,j,:)=subsurf_field(ig0+i,:)
560      enddo
561      out_subsurf_field(lonlen,j,:)=out_subsurf_field(1,j,:) ! redundant lon=180
562    enddo
563   
564    ! write the variable
565    status=nf90_put_var(outid,varid,out_subsurf_field)
566    if (status.ne.nf90_noerr) then
567      write(*,*) "Failed to write ",trim(varname)
568      write(*,*)trim(nf90_strerror(status))
569      stop
570    endif
571  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid)...)
572enddo ! of do i=1,nbinvars
573
574
575! 3. Finish things and cleanup
576! Close output file
577status=nf90_close(outid)
578if (status.ne.nf90_noerr) then
579  write(*,*)"Failed to close file: ",trim(outfile)
580  write(*,*)trim(nf90_strerror(status))
581  stop
582endif
583
584end program expandstartfi
Note: See TracBrowser for help on using the repository browser.