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

Last change on this file since 3468 was 2990, checked in by emillour, 17 months ago

Mars PCM:
Minor fix in utility "expandstartfi": correctly handle the case of
polar mesh areas.
EM

File size: 27.3 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(4),edges(4) ! 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 :: nslope_dimid
31integer :: number_of_advected_fields_dimid
32integer :: time_dimid
33integer :: nbindims ! number of dimensions in input file
34integer :: nbinvars ! number of variables in input file
35integer :: invarid ! to store ID of an input variable
36!Input file variables
37integer :: physical_points
38integer :: subsurface_layers
39integer :: nlayer_plus_1
40integer :: nslope
41integer :: number_of_advected_fields
42integer :: timelen
43real,allocatable :: surf_field(:) ! to store a 1D field of physical_points elements
44real,allocatable :: subsurf_field(:,:) ! to store subsurface (2D field)
45real,allocatable :: subslope_field(:,:)
46real,allocatable :: nlayer_plus_1_field(:,:)
47real,allocatable :: subsurf_subslope_field(:,:,:)
48
49! Output file dimensions:
50integer :: latlen
51integer :: lonlen
52! Output file variables:
53real,allocatable :: longitude(:)
54real,allocatable :: latitude(:)
55real,allocatable :: depth(:)
56real,allocatable :: out_surf_field(:,:)
57real,allocatable :: out_subsurf_field(:,:,:)
58real,allocatable :: out_subslope_field(:,:,:)
59real,allocatable :: out_nlayer_plus_1_field(:,:,:)
60real,allocatable :: out_subsurf_subslope_field(:,:,:,:)
61! Output file dimension IDs
62integer :: lon_dimid
63integer :: lat_dimid
64integer :: depth_dimid
65integer :: nslope_out_dimid
66integer :: nlayer_plus_1_out_dimid
67! IDs of Output file variables
68integer :: lon_varid
69integer :: lat_varid
70integer :: depth_varid
71
72integer :: i,j,k,ig0,ivar
73integer :: nbdim,nbatt,shape(4)
74integer :: nbarg ! # of program arguments
75character(len=256) :: arg ! to store a program argument
76real :: pi ! 3.14159...
77
78pi=2.*asin(1.)
79
80! 0. Preliminary stuff, check arguments (input and output files)
81! get number of arguments this program was called with
82nbarg=command_argument_count()
83
84if (nbarg.ge.1) then
85  call get_command_argument(1,arg) ! get argument # 1
86  infile=trim(arg)
87  if (nbarg.eq.2) then
88    call get_command_argument(2,arg) ! get argument # 2
89    outfile=trim(arg)
90  else
91   ! build outfile from infile (replace ".nc" with "_ex.nc"
92   outfile=trim(infile(1:index(infile,".nc",back=.true.)-1))//"_ex.nc"
93  endif
94  if (nbarg.ge.3) then
95    write(*,*) ' Warning: Too many arguments...'
96    write(*,*) '         will only use the first 2 '
97  endif
98endif
99
100! 1. open input file
101status=nf90_open(infile,NF90_NOWRITE,inid)
102if (status.ne.nf90_noerr) then
103  write(*,*)"Failed to open datafile ",trim(infile)
104  write(*,*)trim(nf90_strerror(status))
105 stop
106endif
107write(*,*) "Reading input file: ",trim(infile)
108
109! 1.2 Identify/load dimensions in input file
110status=nf90_inq_dimid(inid,"physical_points",physical_points_dimid)
111if (status.ne.nf90_noerr) then
112  write(*,*)"Failed to find physical_points dimension"
113  write(*,*)trim(nf90_strerror(status))
114 stop
115endif
116status=nf90_inquire_dimension(inid,physical_points_dimid,len=physical_points)
117if (status.ne.nf90_noerr) then
118  write(*,*)"Failed to find physical_points value"
119  write(*,*)trim(nf90_strerror(status))
120 stop
121else
122  write(*,*) " physical_points = ",physical_points
123endif
124
125status=nf90_inq_dimid(inid,"subsurface_layers",subsurface_layers_dimid)
126if (status.ne.nf90_noerr) then
127  write(*,*)"Failed to find subsurface_layers dimension"
128  write(*,*)trim(nf90_strerror(status))
129 stop
130endif
131status=nf90_inquire_dimension(inid,subsurface_layers_dimid,len=subsurface_layers)
132if (status.ne.nf90_noerr) then
133  write(*,*)"Failed to find subsurface_layers value"
134  write(*,*)trim(nf90_strerror(status))
135 stop
136else
137  write(*,*) " subsurface_layers = ",subsurface_layers
138endif
139
140status=nf90_inq_dimid(inid,"nlayer_plus_1",nlayer_plus_1_dimid)
141if (status.ne.nf90_noerr) then
142  write(*,*)"Failed to find nlayer_plus_1 dimension"
143  write(*,*)trim(nf90_strerror(status))
144 stop
145endif
146status=nf90_inquire_dimension(inid,nlayer_plus_1_dimid,len=nlayer_plus_1)
147if (status.ne.nf90_noerr) then
148  write(*,*)"Failed to find nlayer_plus_1 value"
149  write(*,*)trim(nf90_strerror(status))
150 stop
151else
152  write(*,*) " nlayer_plus_1 = ",nlayer_plus_1
153endif
154
155status=nf90_inq_dimid(inid,"nslope",nslope_dimid)
156if (status.ne.nf90_noerr) then
157  write(*,*)"Failed to find slope dimension, old startfi file"
158  nslope=0
159  write(*,*)trim(nf90_strerror(status))
160else
161status=nf90_inquire_dimension(inid,nslope_dimid,len=nslope)
162if (status.ne.nf90_noerr) then
163  write(*,*)"Failed to find nslope value"
164  write(*,*)trim(nf90_strerror(status))
165 stop
166else
167  write(*,*) " nslope = ",nslope
168endif
169endif
170
171status=nf90_inq_dimid(inid,"Time",time_dimid)
172if (status.ne.nf90_noerr) then
173  write(*,*)"Failed to find Time dimension"
174  write(*,*)trim(nf90_strerror(status))
175  timelen = 0
176else
177  status=nf90_inquire_dimension(inid,time_dimid,len=timelen)
178  if (status.ne.nf90_noerr) then
179    write(*,*)"Failed to read Time dimension"
180    write(*,*)trim(nf90_strerror(status))
181   stop
182  else
183    write(*,*) " time length = ",timelen
184  endif
185endif
186
187! 1.3 Allocate memory for input fields
188allocate(surf_field(physical_points))
189allocate(subsurf_field(physical_points,subsurface_layers))
190allocate(subslope_field(physical_points,nslope))
191allocate(nlayer_plus_1_field(physical_points,nlayer_plus_1))
192allocate(subsurf_subslope_field(physical_points,subsurface_layers,nslope))
193
194! 2.1. Create output file
195status=nf90_create(outfile,NF90_CLOBBER,outid)
196if (status.ne.nf90_noerr) then
197  write(*,*) "Failed to create output file: ",trim(outfile)
198  write(*,*)trim(nf90_strerror(status))
199  stop
200endif
201write(*,*) " "
202write(*,*) "Created output file: ",trim(outfile)
203
204! 2.2. Build dimensions for output file
205
206! 2.2.1 Compute lonlen and latlen from  physical_points
207! load "longitude(physical_points)" from input file
208status=nf90_inq_varid(inid,"longitude",varid)
209if (status.ne.nf90_noerr) then
210  write(*,*) "Failed to find longitude ID"
211  write(*,*)trim(nf90_strerror(status))
212  stop
213endif
214! read longitude
215status=nf90_get_var(inid,varid,surf_field)
216if (status.ne.nf90_noerr) then
217  write(*,*) "Failed to load longitude"
218  write(*,*)trim(nf90_strerror(status))
219  stop
220endif
221
222! count the number of points before longitude(i) wraps around
223i=3
224lonlen=1
225!write(*,*) "longitude(2)=",surf_field(2)
226do while (surf_field(i).ne.surf_field(2))
227!write(*,*) "i:",i,"surf_field(i)=",surf_field(i)
228 i=i+1
229 lonlen=lonlen+1
230enddo
231! and add 1 because we want a redundant lon=180 point
232lonlen=lonlen+1
233write(*,*) " => lonlen=",lonlen
234
235allocate(longitude(lonlen))
236! fill longitude(1:lonlen)
237longitude(1:lonlen-1)=surf_field(2:lonlen)
238longitude(lonlen)=-longitude(1) ! redundant +Pi/2
239! convert to degrees
240longitude(:)=longitude(:)*(180./pi)
241
242! knowing lonlen, compute latlen
243latlen=2+(physical_points-2)/(lonlen-1)
244write(*,*) " => latlen=",latlen
245
246allocate(latitude(latlen))
247! get field "latitude(physical_points)" from infile
248status=nf90_inq_varid(inid,"latitude",varid)
249if (status.ne.nf90_noerr) then
250  write(*,*) "Failed to find latitude ID"
251  write(*,*)trim(nf90_strerror(status))
252  stop
253endif
254! read latitude
255status=nf90_get_var(inid,varid,surf_field)
256if (status.ne.nf90_noerr) then
257  write(*,*) "Failed to load latitude"
258  write(*,*)trim(nf90_strerror(status))
259  stop
260endif
261
262latitude(1)=surf_field(1)
263!write(*,*) "latitude(1)=",latitude(1)
264do i=2,latlen-1
265  latitude(i)=surf_field((i-1)*(lonlen-1))
266!  write(*,*) "i:",i,"latitude(i)=",latitude(i)
267enddo
268latitude(latlen)=surf_field(physical_points)
269!write(*,*) "latitude(latlen)=",latitude(latlen)
270! convert to degrees
271latitude(:)=latitude(:)*(180./pi)
272
273! depth:
274allocate(depth(subsurface_layers))
275! load "soildepth(subsurface_layers)" from input file
276status=nf90_inq_varid(inid,"soildepth",varid)
277if (status.ne.nf90_noerr) then
278  write(*,*) "Failed to find soildepth ID"
279  write(*,*)trim(nf90_strerror(status))
280  stop
281endif
282! read longitude
283status=nf90_get_var(inid,varid,depth)
284if (status.ne.nf90_noerr) then
285  write(*,*) "Failed to load soildepth"
286  write(*,*)trim(nf90_strerror(status))
287  stop
288endif
289
290! 2.2.2 define dimensions to output file
291! longitude:
292status=nf90_def_dim(outid,"longitude",lonlen,lon_dimid)
293if (status.ne.nf90_noerr) then
294  write(*,*) "Failed creating longitude dimension"
295  write(*,*)trim(nf90_strerror(status))
296  stop
297endif
298! latitude:
299status=nf90_def_dim(outid,"latitude",latlen,lat_dimid)
300if (status.ne.nf90_noerr) then
301  write(*,*) "Failed creating longitude dimension"
302  write(*,*)trim(nf90_strerror(status))
303  stop
304endif
305! depth:
306status=nf90_def_dim(outid,"depth",subsurface_layers,depth_dimid)
307if (status.ne.nf90_noerr) then
308  write(*,*) "Failed creating depth dimension"
309  write(*,*)trim(nf90_strerror(status))
310  stop
311endif
312! nslope:
313status=nf90_def_dim(outid,"nslope",nslope,nslope_out_dimid)
314if (status.ne.nf90_noerr) then
315  write(*,*) "Failed creating nslope dimension"
316  write(*,*)trim(nf90_strerror(status))
317  stop
318endif
319! nslope:
320status=nf90_def_dim(outid,"nlayer_plus_1",nlayer_plus_1,nlayer_plus_1_out_dimid)
321if (status.ne.nf90_noerr) then
322  write(*,*) "Failed creating nslope dimension"
323  write(*,*)trim(nf90_strerror(status))
324  stop
325endif
326
327!2.3. write variables to output file
328! 2.3.1. define 1D variables
329! longitude
330datashape(1)=lon_dimid
331status=nf90_def_var(outid,"longitude",NF90_REAL,lon_dimid,lon_varid)
332if (status.ne.nf90_noerr) then
333  write(*,*) "Failed creating longitude variable"
334  write(*,*)trim(nf90_strerror(status))
335  stop
336endif
337! longitude attributes
338status=nf90_put_att(outid,lon_varid,"long_name","East longitude")
339status=nf90_put_att(outid,lon_varid,"units","degrees_east")
340
341!latitude
342datashape(2)=lat_dimid
343status=nf90_def_var(outid,"latitude",NF90_REAL,lat_dimid,lat_varid)
344if (status.ne.nf90_noerr) then
345  write(*,*) "Failed creating latitude variable"
346  write(*,*)trim(nf90_strerror(status))
347  stop
348endif
349! latitude attributes
350status=nf90_put_att(outid,lat_varid,"long_name","North latitude")
351status=nf90_put_att(outid,lat_varid,"units","degrees_north")
352
353! depth
354status=nf90_def_var(outid,"depth",NF90_REAL,depth_dimid,depth_varid)
355if (status.ne.nf90_noerr) then
356  write(*,*) "Failed creating depth variable"
357  write(*,*)trim(nf90_strerror(status))
358  stop
359endif
360!depth atributes
361status=nf90_put_att(outid,depth_varid,"long_name","Soil mid-layer depth")
362status=nf90_put_att(outid,depth_varid,"units","m")
363status=nf90_put_att(outid,depth_varid,"positive","down")
364
365! 2.3.2 write 1D variable
366
367! swich out of NetCDF define mode
368status=nf90_enddef(outid)
369if (status.ne.nf90_noerr) then
370  write(*,*) "Failed to swich out of define mode"
371  write(*,*)trim(nf90_strerror(status))
372  stop
373endif
374
375! write longitude
376status=nf90_put_var(outid,lon_varid,longitude)
377if (status.ne.nf90_noerr) then
378  write(*,*) "Failed writing longitude"
379  write(*,*)trim(nf90_strerror(status))
380  stop
381endif
382
383! write latitude
384status=nf90_put_var(outid,lat_varid,latitude)
385if (status.ne.nf90_noerr) then
386  write(*,*) "Failed writing latitude"
387  write(*,*)trim(nf90_strerror(status))
388  stop
389endif
390
391! write depth
392status=nf90_put_var(outid,depth_varid,depth)
393if (status.ne.nf90_noerr) then
394  write(*,*) "Failed writing depth"
395  write(*,*)trim(nf90_strerror(status))
396  stop
397endif
398
399! 2.3. 2D (surface) variables
400! First find out how many variables there are in input file
401status=nf90_inquire(inid,nbindims,nbinvars)
402if (status.ne.nf90_noerr) then
403  write(*,*) "Failed obtaining nbindims and nbinvars from input file"
404  write(*,*)trim(nf90_strerror(status))
405  stop
406endif
407
408allocate(out_surf_field(lonlen,latlen))
409allocate(out_subsurf_field(lonlen,latlen,subsurface_layers))
410allocate(out_subslope_field(lonlen,latlen,nslope))
411allocate(out_nlayer_plus_1_field(lonlen,latlen,nlayer_plus_1))
412allocate(out_subsurf_subslope_field(lonlen,latlen,subsurface_layers,nslope))
413
414do ivar=1,nbinvars ! loop on all input variables
415  shape(:) = 0
416  ! find out what dimensions are linked to this variable
417  status=nf90_inquire_variable(inid,ivar,name=varname,ndims=nbdim,&
418                               dimids=shape,natts=nbatt)
419
420  if (((nbdim==1).and.(shape(1)==physical_points_dimid))&
421  .or.((nbdim==2).and.(shape(1)==physical_points_dimid)&
422                 .and.(shape(2)==time_dimid))) then
423 
424    corner(1) = 1
425    corner(2) = timelen
426    edges(1)  = physical_points
427    edges(2)  = 1
428   
429    ! skip "longitude" and "latitude"
430    if (trim(varname)=="longitude") cycle
431    if (trim(varname)=="latitude") cycle
432   
433    write(*,*) " processing: ",trim(varname)
434
435    ! load input data:
436    status=nf90_inq_varid(inid,varname,invarid)
437    status=nf90_get_var(inid,invarid,surf_field,corner,edges)
438   
439    ! switch output file to to define mode
440    status=nf90_redef(outid)
441    if (status.ne.nf90_noerr) then
442      write(*,*) "Failed to swich to define mode"
443      write(*,*)trim(nf90_strerror(status))
444      stop
445    endif
446    !define the variable
447    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
448                         (/lon_dimid,lat_dimid/),varid)
449    if (status.ne.nf90_noerr) then
450      write(*,*) "Failed creating variable ",trim(varname)
451      write(*,*)trim(nf90_strerror(status))
452      stop
453    endif
454
455    ! variable attributes
456    do k=1,nbatt
457      status=nf90_inq_attname(inid,invarid,k,varatt)
458      if (status.ne.nf90_noerr) then
459        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
460        write(*,*)trim(nf90_strerror(status))
461      stop
462      endif
463      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
464      if (status.ne.nf90_noerr) then
465        write(*,*) "Failed loading attribute ",trim(varatt)
466        write(*,*)trim(nf90_strerror(status))
467      stop
468      endif
469
470      ! write the attribute and its value to output
471      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
472      if (status.ne.nf90_noerr) then
473        write(*,*) "Failed writing attribute ",trim(varatt)
474        write(*,*)trim(nf90_strerror(status))
475      stop
476      endif
477    enddo ! do k=1,nbatt
478
479    ! swich out of NetCDF define mode
480    status=nf90_enddef(outid)
481    if (status.ne.nf90_noerr) then
482      write(*,*) "Failed to swich out of define mode"
483      write(*,*)trim(nf90_strerror(status))
484      stop
485    endif
486   
487    ! "convert" from surf_field(:) to out_surf_field(:,:)
488    if (trim(varname)=="area") then
489     ! Very specific handling of mesh area at the poles
490     ! which need to be split along all longitudes
491     do i=1,lonlen
492      out_surf_field(i,1)=surf_field(1)/(lonlen-1)
493      out_surf_field(i,latlen)=surf_field(physical_points)/(lonlen-1)
494     enddo   
495    else ! All other variables
496     do i=1,lonlen
497      out_surf_field(i,1)=surf_field(1) ! North Pole
498      out_surf_field(i,latlen)=surf_field(physical_points) ! South Pole
499     enddo
500    endif ! of if (trim(varname)=="area")
501    do j=2,latlen-1
502      ig0=1+(j-2)*(lonlen-1)
503      do i=1,lonlen-1
504        out_surf_field(i,j)=surf_field(ig0+i)
505      enddo
506      out_surf_field(lonlen,j)=out_surf_field(1,j) ! redundant lon=180
507    enddo
508   
509    ! write the variable
510    status=nf90_put_var(outid,varid,out_surf_field)
511    if (status.ne.nf90_noerr) then
512      write(*,*) "Failed to write ",trim(varname)
513      write(*,*)trim(nf90_strerror(status))
514      stop
515    endif
516  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid))
517
518  if ((nbdim==2).and.(shape(1)==physical_points_dimid)&
519                 .and.(shape(2)==subsurface_layers_dimid)) then
520   
521    corner(1) = 1
522    corner(2) = 1
523    corner(3) = timelen
524    edges(1)  = physical_points
525    edges(2)  = subsurface_layers
526    edges(3)  = 1
527   
528    write(*,*) " processing: ",trim(varname)
529
530    ! load input data:
531    status=nf90_inq_varid(inid,varname,invarid)
532    status=nf90_get_var(inid,invarid,subsurf_field,corner,edges)
533   
534    ! switch output file to to define mode
535    status=nf90_redef(outid)
536    if (status.ne.nf90_noerr) then
537      write(*,*) "Failed to swich to define mode"
538      write(*,*)trim(nf90_strerror(status))
539      stop
540    endif
541    !define the variable
542    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
543                         (/lon_dimid,lat_dimid,depth_dimid/),varid)
544    if (status.ne.nf90_noerr) then
545      write(*,*) "Failed creating variable ",trim(varname)
546      write(*,*)trim(nf90_strerror(status))
547      stop
548    endif
549
550    ! variable attributes
551    do k=1,nbatt
552      status=nf90_inq_attname(inid,invarid,k,varatt)
553      if (status.ne.nf90_noerr) then
554        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
555        write(*,*)trim(nf90_strerror(status))
556      stop
557      endif
558      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
559      if (status.ne.nf90_noerr) then
560        write(*,*) "Failed loading attribute ",trim(varatt)
561        write(*,*)trim(nf90_strerror(status))
562      stop
563      endif
564
565      ! write the attribute and its value to output
566      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
567      if (status.ne.nf90_noerr) then
568        write(*,*) "Failed writing attribute ",trim(varatt)
569        write(*,*)trim(nf90_strerror(status))
570      stop
571      endif
572    enddo ! do k=1,nbatt
573
574    ! swich out of NetCDF define mode
575    status=nf90_enddef(outid)
576    if (status.ne.nf90_noerr) then
577      write(*,*) "Failed to swich out of define mode"
578      write(*,*)trim(nf90_strerror(status))
579      stop
580    endif
581   
582    ! "convert" from subsurf_field(:,:) to out_subsurf_field(:,:,:)
583    do i=1,lonlen
584      out_subsurf_field(i,1,:)=subsurf_field(1,:) ! North Pole
585      out_subsurf_field(i,latlen,:)=subsurf_field(physical_points,:) ! South Pole
586    enddo
587    do j=2,latlen-1
588      ig0=1+(j-2)*(lonlen-1)
589      do i=1,lonlen-1
590        out_subsurf_field(i,j,:)=subsurf_field(ig0+i,:)
591      enddo
592      out_subsurf_field(lonlen,j,:)=out_subsurf_field(1,j,:) ! redundant lon=180
593    enddo
594   
595    ! write the variable
596    status=nf90_put_var(outid,varid,out_subsurf_field)
597    if (status.ne.nf90_noerr) then
598      write(*,*) "Failed to write ",trim(varname)
599      write(*,*)trim(nf90_strerror(status))
600      stop
601    endif
602  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid)...)
603
604
605
606
607  if ((nbdim==3).and.(shape(1)==physical_points_dimid)&
608                 .and.(shape(2)==nslope_dimid)) then
609   
610    corner(1) = 1
611    corner(2) = 1
612    corner(3) = timelen
613    edges(1)  = physical_points
614    edges(2)  = nslope
615    edges(3)  = 1
616   
617    write(*,*) " processing: ",trim(varname)
618
619    ! load input data:
620    status=nf90_inq_varid(inid,varname,invarid)
621    status=nf90_get_var(inid,invarid,subslope_field,corner,edges)
622   
623    ! switch output file to to define mode
624    status=nf90_redef(outid)
625    if (status.ne.nf90_noerr) then
626      write(*,*) "Failed to swich to define mode"
627      write(*,*)trim(nf90_strerror(status))
628      stop
629    endif
630    !define the variable
631    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
632                         (/lon_dimid,lat_dimid,nslope_out_dimid/),varid)
633    if (status.ne.nf90_noerr) then
634      write(*,*) "Failed creating variable ",trim(varname)
635      write(*,*)trim(nf90_strerror(status))
636      stop
637    endif
638
639    ! variable attributes
640    do k=1,nbatt
641      status=nf90_inq_attname(inid,invarid,k,varatt)
642      if (status.ne.nf90_noerr) then
643        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
644        write(*,*)trim(nf90_strerror(status))
645      stop
646      endif
647      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
648      if (status.ne.nf90_noerr) then
649        write(*,*) "Failed loading attribute ",trim(varatt)
650        write(*,*)trim(nf90_strerror(status))
651      stop
652      endif
653
654      ! write the attribute and its value to output
655      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
656      if (status.ne.nf90_noerr) then
657        write(*,*) "Failed writing attribute ",trim(varatt)
658        write(*,*)trim(nf90_strerror(status))
659      stop
660      endif
661    enddo ! do k=1,nbatt
662
663    ! swich out of NetCDF define mode
664    status=nf90_enddef(outid)
665    if (status.ne.nf90_noerr) then
666      write(*,*) "Failed to swich out of define mode"
667      write(*,*)trim(nf90_strerror(status))
668      stop
669    endif
670   
671    ! "convert" from subsurf_field(:,:) to out_subslope_field(:,:,:)
672    do i=1,lonlen
673      out_subslope_field(i,1,:)=subslope_field(1,:) ! North Pole
674      out_subslope_field(i,latlen,:)=subslope_field(physical_points,:) ! South Pole
675    enddo
676    do j=2,latlen-1
677      ig0=1+(j-2)*(lonlen-1)
678      do i=1,lonlen-1
679        out_subslope_field(i,j,:)=subslope_field(ig0+i,:)
680      enddo
681      out_subslope_field(lonlen,j,:)=out_subslope_field(1,j,:) ! redundant lon=180
682    enddo
683   
684    ! write the variable
685    status=nf90_put_var(outid,varid,out_subslope_field)
686    if (status.ne.nf90_noerr) then
687      write(*,*) "Failed to write ",trim(varname)
688      write(*,*)trim(nf90_strerror(status))
689      stop
690    endif
691  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid)...)
692
693
694
695
696
697  if ((nbdim==3).and.(shape(1)==physical_points_dimid)&
698                 .and.(shape(2)==nlayer_plus_1_dimid)) then
699   
700    corner(1) = 1
701    corner(2) = 1
702    corner(3) = timelen
703    edges(1)  = physical_points
704    edges(2)  = nlayer_plus_1
705    edges(3)  = 1
706   
707    write(*,*) " processing: ",trim(varname)
708
709    ! load input data:
710    status=nf90_inq_varid(inid,varname,invarid)
711    status=nf90_get_var(inid,invarid,nlayer_plus_1_field,corner,edges)
712   
713    ! switch output file to to define mode
714    status=nf90_redef(outid)
715    if (status.ne.nf90_noerr) then
716      write(*,*) "Failed to swich to define mode"
717      write(*,*)trim(nf90_strerror(status))
718      stop
719    endif
720    !define the variable
721    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
722                         (/lon_dimid,lat_dimid,nlayer_plus_1_out_dimid/),varid)
723    if (status.ne.nf90_noerr) then
724      write(*,*) "Failed creating variable ",trim(varname)
725      write(*,*)trim(nf90_strerror(status))
726      stop
727    endif
728
729    ! variable attributes
730    do k=1,nbatt
731      status=nf90_inq_attname(inid,invarid,k,varatt)
732      if (status.ne.nf90_noerr) then
733        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
734        write(*,*)trim(nf90_strerror(status))
735      stop
736      endif
737      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
738      if (status.ne.nf90_noerr) then
739        write(*,*) "Failed loading attribute ",trim(varatt)
740        write(*,*)trim(nf90_strerror(status))
741      stop
742      endif
743
744      ! write the attribute and its value to output
745      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
746      if (status.ne.nf90_noerr) then
747        write(*,*) "Failed writing attribute ",trim(varatt)
748        write(*,*)trim(nf90_strerror(status))
749      stop
750      endif
751    enddo ! do k=1,nbatt
752
753    ! swich out of NetCDF define mode
754    status=nf90_enddef(outid)
755    if (status.ne.nf90_noerr) then
756      write(*,*) "Failed to swich out of define mode"
757      write(*,*)trim(nf90_strerror(status))
758      stop
759    endif
760   
761    ! "convert" from subsurf_field(:,:) to out_nlayer_plus_1_field(:,:,:)
762    do i=1,lonlen
763      out_nlayer_plus_1_field(i,1,:)=nlayer_plus_1_field(1,:) ! North Pole
764      out_nlayer_plus_1_field(i,latlen,:)=nlayer_plus_1_field(physical_points,:) ! South Pole
765    enddo
766    do j=2,latlen-1
767      ig0=1+(j-2)*(lonlen-1)
768      do i=1,lonlen-1
769        out_nlayer_plus_1_field(i,j,:)=nlayer_plus_1_field(ig0+i,:)
770      enddo
771      out_nlayer_plus_1_field(lonlen,j,:)=out_nlayer_plus_1_field(1,j,:) ! redundant lon=180
772    enddo
773   
774    ! write the variable
775    status=nf90_put_var(outid,varid,out_nlayer_plus_1_field)
776    if (status.ne.nf90_noerr) then
777      write(*,*) "Failed to write ",trim(varname)
778      write(*,*)trim(nf90_strerror(status))
779      stop
780    endif
781  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid)...)
782
783
784
785
786  if (nbdim==4) then
787   
788    corner(1) = 1
789    corner(2) = 1
790    corner(3) = 1
791    corner(4) = timelen
792    edges(1)  = physical_points
793    edges(2)  = subsurface_layers
794    edges(3)  = nslope
795    edges(4)  = 1
796   
797    write(*,*) " processing: ",trim(varname)
798
799    ! load input data:
800    status=nf90_inq_varid(inid,varname,invarid)
801    status=nf90_get_var(inid,invarid,subsurf_subslope_field,corner,edges)
802   
803    ! switch output file to to define mode
804    status=nf90_redef(outid)
805    if (status.ne.nf90_noerr) then
806      write(*,*) "Failed to swich to define mode"
807      write(*,*)trim(nf90_strerror(status))
808      stop
809    endif
810    !define the variable
811    print *, "subsurface_layers_dimid", subsurface_layers_dimid
812    print *, "nslope_out_dimid", nslope_out_dimid
813    status=nf90_def_var(outid,trim(varname),NF90_REAL,&
814                         (/lon_dimid,lat_dimid,depth_dimid,nslope_out_dimid/),varid)
815    if (status.ne.nf90_noerr) then
816      write(*,*) "Failed creating variable ",trim(varname)
817      write(*,*)trim(nf90_strerror(status))
818      stop
819    endif
820
821    ! variable attributes
822    do k=1,nbatt
823      status=nf90_inq_attname(inid,invarid,k,varatt)
824      if (status.ne.nf90_noerr) then
825        write(*,*) "Failed getting attribute number",k," for ",trim(varname)
826        write(*,*)trim(nf90_strerror(status))
827      stop
828      endif
829      status=nf90_get_att(inid,invarid,trim(varatt),varattcontent)
830      if (status.ne.nf90_noerr) then
831        write(*,*) "Failed loading attribute ",trim(varatt)
832        write(*,*)trim(nf90_strerror(status))
833      stop
834      endif
835
836      ! write the attribute and its value to output
837      status=nf90_put_att(outid,varid,trim(varatt),trim(varattcontent))
838      if (status.ne.nf90_noerr) then
839        write(*,*) "Failed writing attribute ",trim(varatt)
840        write(*,*)trim(nf90_strerror(status))
841      stop
842      endif
843    enddo ! do k=1,nbatt
844
845    ! swich out of NetCDF define mode
846    status=nf90_enddef(outid)
847    if (status.ne.nf90_noerr) then
848      write(*,*) "Failed to swich out of define mode"
849      write(*,*)trim(nf90_strerror(status))
850      stop
851    endif
852   
853    ! "convert" from subsurf_field(:,:) to out_subsurf_subslope_field(:,:,:)
854    do i=1,lonlen
855      out_subsurf_subslope_field(i,1,:,:)=subsurf_subslope_field(1,:,:) ! North Pole
856      out_subsurf_subslope_field(i,latlen,:,:)=subsurf_subslope_field(physical_points,:,:) ! South Pole
857    enddo
858    do j=2,latlen-1
859      ig0=1+(j-2)*(lonlen-1)
860      do i=1,lonlen-1
861        out_subsurf_subslope_field(i,j,:,:)=subsurf_subslope_field(ig0+i,:,:)
862      enddo
863      out_subsurf_subslope_field(lonlen,j,:,:)=out_subsurf_subslope_field(1,j,:,:) ! redundant lon=180
864    enddo
865   
866    ! write the variable
867    status=nf90_put_var(outid,varid,out_subsurf_subslope_field)
868    if (status.ne.nf90_noerr) then
869      write(*,*) "Failed to write ",trim(varname)
870      write(*,*)trim(nf90_strerror(status))
871      stop
872    endif
873  endif ! of if ((nbdim==1).and.(shape(1)==physical_points_dimid)...)
874
875
876
877enddo ! of do i=1,nbinvars
878
879! 3. Finish things and cleanup
880! Close output file
881status=nf90_close(outid)
882if (status.ne.nf90_noerr) then
883  write(*,*)"Failed to close file: ",trim(outfile)
884  write(*,*)trim(nf90_strerror(status))
885  stop
886endif
887
888  write(*,*)"Done writing file ",trim(outfile)
889
890end program expandstartfi
Note: See TracBrowser for help on using the repository browser.