source: trunk/MESOSCALE/LMD_MM_MARS/SRC/PREP_MARS/ukdiag/newdiag.F90 @ 1601

Last change on this file since 1601 was 1440, checked in by aslmd, 10 years ago

MESOSCALE. changes inspired by Liam Steele. mcmodel in makemeso and mpifort for larger domains. cleaning the PREP_MARS repository. added Liam's ukdiag tool to convert UK dynamical core outputs for use in mesoscale simulations.

  • Property svn:executable set to *
File size: 13.5 KB
Line 
1! ====================================================================
2!
3! This routine reads in a diagfi.nc file output from the UK-LMD MGCM
4! and converts it into a format suitable for producing initial and
5! boundary conditions for the mesoscale model. Specifically, the
6! arrays are interpolated onto latitudes ranging from +90 --> -90
7! and longitudes ranging from -180 --> +180. The surface geopotential
8! is also added if it doesn't already exist, and a check is made to
9! see if the starting sol is present in the controle array.
10!
11! Arrays written on regolith layers are now interpolated back on to
12! soil layers for use in the mesoscale model (sometimes different
13! layers had to be specified for the regolith scheme (compared to
14! the soil layers) to ensure numerical stability at high obliquities).
15!
16! The new file created is called diagfi_new.nc, which can be linked
17! to the PREP_MARS directory as input_diagfi.nc
18!
19! Liam Steele June 2014, March 2015
20!
21! ====================================================================
22
23program newdiag
24implicit none
25
26include "netcdf.inc"
27
28integer :: i,j,k,ierr,nid,nid2,nid3
29integer :: ndims,nvars,ngatts,unlimdimid
30integer :: lonid,latid,sigid,soilid,soildim,regid,timeid,contid,phisid
31integer :: nlatnew,nlonnew,nsoil,nreg,npos
32integer :: vtype,vatts,vdims,solstart
33integer, dimension(4) :: dimids
34integer, dimension(:), allocatable :: id,dimid,dimarr,loc,corners,edges
35real, dimension(:), allocatable :: lat,lon,sig,soil,regol,time,controle
36real, dimension(:), allocatable :: newlat,newlon,weight
37real, dimension(:,:), allocatable :: geopot,geopot_new,array2,array2_new
38real, dimension(:,:,:), allocatable :: array3,array3_new
39real, dimension(:,:,:,:), allocatable :: array4,array4_new
40real, dimension(:,:,:,:), allocatable :: new_array
41character*20, dimension(:), allocatable :: dimname
42character*20 :: vname
43character*2 :: printid,res
44logical :: wphis,skipv
45
46! Open file for reading
47ierr = nf_open("diagfi.nc",nf_nowrite,nid)
48if (ierr /= nf_noerr) then
49  print*, 'Oh dear: cannot find input_diagfi.nc'
50  stop
51endif
52
53! Check properties of diagfi
54ierr = nf_inq(nid, ndims, nvars, ngatts, unlimdimid)
55print*, '--------------'
56print*, 'File contains:'
57write(printid,'(i2)'),ndims
58print*, printid//' dimensions'
59write(printid,'(i2)'),nvars
60print*, printid//' variables'
61print*, '--------------'
62allocate(dimarr(ndims))
63allocate(dimname(ndims))
64
65! Get dimensions
66do i = 1, ndims
67  ierr = nf_inq_dim(nid,i,dimname(i),dimarr(i))
68  if (ierr /= nf_noerr) then
69    print*, 'Whoops: problem reading dimensions in diagfi.nc'
70    stop
71  endif
72  if (trim(dimname(i)) == 'latitude' .or. trim(dimname(i)) == 'lat') then
73    dimname(i) = 'latitude'
74    allocate(lat(dimarr(i)))
75    ierr = nf_get_var_real(nid,i,lat)
76    nlatnew = dimarr(i)+1
77  else if (trim(dimname(i)) == 'longitude' .or. trim(dimname(i)) == 'lon') then
78    dimname(i) = 'longitude'
79    allocate(lon(dimarr(i)))
80    ierr = nf_get_var_real(nid,i,lon)
81    nlonnew = dimarr(i)+1
82  else if (trim(dimname(i)) == 'sigma') then
83    allocate(sig(dimarr(i)))
84    ierr = nf_get_var_real(nid,i,sig)
85  else if (trim(dimname(i)) == 'soildepth' .or. trim(dimname(i)) == 'soil') then
86    dimname(i) = 'soildepth'
87    allocate(soil(dimarr(i)))
88    ierr = nf_get_var_real(nid,i,soil)
89    nsoil = dimarr(i)
90    soildim = i
91  else if (trim(dimname(i)) == 'regdepth') then
92    dimname(i) = 'regdepth'
93    allocate(regol(dimarr(i)))
94    ierr = nf_get_var_real(nid,i,regol)
95    nreg = dimarr(i)
96  else if (trim(dimname(i)) == 'Time' .or. trim(dimname(i)) == 'time') then
97    dimname(i) = 'Time'
98    allocate(time(dimarr(i)))
99    ierr = nf_get_var_real(nid,i,time)
100  else if (trim(dimname(i)) == 'lentable') then
101    allocate(controle(dimarr(i)))
102    ierr = nf_get_var_real(nid,i,controle)
103  endif
104enddo
105
106! Check for starting sol in diagfi
107if (controle(4) == 0) then
108  print*, 'No starting sol sumber in contole array. Please enter [1-669]:'
109  read*, solstart
110  do while (solstart > 669 .or. solstart < 0)
111    if (solstart > 669) print*, 'Sol > 669. Have another go...'
112    if (solstart < 0) print*, 'A negative sol is just silly. Try again...'
113    read*, solstart
114  enddo
115  controle(4) = solstart
116else if (controle(4) > 669) then
117  controle(4) = mod(controle(4),669.)
118endif
119
120! Add soil mid-layer depths if not already present
121if (soil(1) == 0 .or. soil(1) == 1) then
122  print*, 'Adding the following soil levels (depths in m):'
123  do i = 1, nsoil
124    soil(i) = 2.e-4*(2.**(i-1.5))
125    print*, i, soil(i)
126  enddo
127endif
128
129! Make new lat/lon arrays
130allocate(newlat(nlatnew))
131allocate(newlon(nlonnew))
132do i = 1, nlatnew
133  newlat(i) = 90. + (i-1)*(lat(2)-lat(1))
134enddo
135do i = 1, nlonnew
136  newlon(i) = -180. + (i-1)*(lon(2)-lon(1))
137enddo
138
139! Set weights and locations for latitudinal interpolation (assumes linear spacing between lats)
140allocate(weight(nlatnew))
141allocate(loc(nlatnew))
142weight(:) = 0.5
143weight(1) = 1.5
144weight(nlatnew) = -0.5
145do i = 1, nlatnew
146  loc(i) = i-1
147enddo
148loc(1) = 1
149loc(nlatnew) = nlatnew-2
150
151! Create new netcdf file to write data to (doesn't matter about writing descriptions and units)
152ierr = nf_create('diagfi_new.nc',nf_clobber,nid2)
153
154! Define dimensions
155do i = 1, ndims
156  write(printid,'(i2)'),i
157  print*, printid//' writing '//trim(dimname(i))//' to file'
158  if (trim(dimname(i)) == 'longitude') then
159    ierr = nf_def_dim(nid2,'longitude',nlonnew,lonid)
160    ierr = nf_def_var(nid2,'longitude',nf_float,1,lonid,lonid)
161  else if (trim(dimname(i)) == 'latitude') then
162    ierr = nf_def_dim(nid2,'latitude',nlatnew,latid)
163    ierr = nf_def_var(nid2,'latitude',nf_float,1,latid,latid)
164  else  if (trim(dimname(i)) == 'sigma') then
165    ierr = nf_def_dim(nid2,'sigma',size(sig),sigid)
166    ierr = nf_def_var(nid2,'sigma',nf_float,1,sigid,sigid)
167  else if (trim(dimname(i)) == 'soildepth') then
168    ierr = nf_def_dim(nid2,'soildepth',size(soil),soilid)
169    ierr = nf_def_var(nid2,'soildepth',nf_float,1,soilid,soilid)
170  else if (trim(dimname(i)) == 'regdepth') then
171    ierr = nf_def_dim(nid2,'regdepth',size(regol),regid)
172    ierr = nf_def_var(nid2,'regdepth',nf_float,1,regid,regid)
173  else if (trim(dimname(i)) == 'Time') then
174    ierr = nf_def_dim(nid2,'Time',size(time),timeid)
175    ierr = nf_def_var(nid2,'Time',nf_float,1,timeid,timeid)
176  else if (trim(dimname(i)) == 'lentable') then
177    ierr = nf_def_dim(nid2,'lentable',size(controle),contid)
178    ierr = nf_def_var(nid2,'controle',nf_float,1,contid,contid)
179  endif
180  if (ierr /= nf_noerr) then
181    print*, "Yikes! Problem defining dimensions in new_diagfi.nc"
182    stop
183  endif
184enddo
185
186ierr = nf_enddef(nid2)
187
188! Write dimension arrays
189ierr = nf_put_var_real(nid2,lonid,newlon)
190ierr = nf_put_var_real(nid2,latid,newlat)
191ierr = nf_put_var_real(nid2,sigid,sig)
192ierr = nf_put_var_real(nid2,soilid,soil)
193ierr = nf_put_var_real(nid2,regid,regol)
194ierr = nf_put_var_real(nid2,timeid,time)
195ierr = nf_put_var_real(nid2,contid,controle)
196
197! Read each variable in turn from original file, convert to new lat/lon and write to new file
198npos = ndims+1
199wphis = .true.
200do i = ndims+1, nvars
201
202  ! Read variable dimensions etc.
203  ierr = nf_inq_var(nid, i, vname, vtype, vdims, dimids, vatts)
204  write(printid,'(i2)'),i
205  if (trim(vname) == 'phisinit') wphis = .false.
206  allocate(id(vdims))
207  allocate(corners(vdims))
208  allocate(edges(vdims))
209  allocate(dimid(vdims))
210 
211  ! Require lon/lat to be first two dimensions (everything should satisfy this condition, but you never know!)
212  skipv = .true.
213  if (vdims >= 2) then
214    if (dimname(dimids(1)) == 'longitude' .and. dimname(dimids(2)) == 'latitude')  skipv = .false.
215    if (dimname(dimids(1)) == 'latitude'  .and. dimname(dimids(2)) == 'longitude') skipv = .false.
216    if (skipv) then
217      print*, " > "//trim(vname)//": lat and lon aren't first two dimensions. Skipping"
218      cycle
219    endif
220  endif
221 
222  ! Change certain names to match LMD (more might need added later)
223  if (trim(vname) == 'h2ommr')    vname = 'q02'
224  if (trim(vname) == 'icemmr')    vname = 'q01'
225  if (trim(vname) == 'h2o_ice_s') vname = 'qsurf01'
226 
227  print*, printid//' writing '//trim(vname)//' to file'
228 
229  ! Define arrays for specifying location in diagfi
230  do k = 1, vdims
231    id(k) = dimarr(dimids(k))
232    dimid(k) = dimids(k)
233    corners(k) = 1
234    edges(k) = dimarr(dimids(k))
235  enddo
236  edges(lonid) = nlonnew
237  edges(latid) = nlatnew
238  edges(vdims) = 1
239 
240  ! Allocate and read old arrays
241  if (vdims == 2) then
242    allocate(array2(id(1),id(2)))
243    ierr = nf_get_var_real(nid,i,array2)
244  else if (vdims == 3) then
245    allocate(array3(id(1),id(2),id(3)))
246    ierr = nf_get_var_real(nid,i,array3)
247  else if (vdims == 4) then
248    allocate(array4(id(1),id(2),id(3),id(4)))
249    ierr = nf_get_var_real(nid,i,array4)
250  endif
251 
252  ! Check if regolith depth is the vertical dimension (if so, will convert fields to soil depths)
253  do k = 1, vdims
254    if (dimname(dimids(k)).eq.'regdepth') then
255   
256      if (vdims == 3) then
257        allocate(new_array(id(1),id(2),dimarr(soildim),1))
258        call interpol_soil(soil,regol,nsoil,nreg,(/id(1:vdims),1/),array3,new_array)
259        deallocate(array3)
260        allocate(array3(id(1),id(2),dimarr(soildim)))
261        array3 = new_array(1:id(1),1:id(2),1:dimarr(soildim),1)
262        deallocate(new_array)
263      else if (vdims == 4) then
264        allocate(new_array(id(1),id(2),dimarr(soildim),id(4)))
265        call interpol_soil(soil,regol,nsoil,nreg,id(1:vdims),array4,new_array)
266        deallocate(array4)
267        allocate(array4(id(1),id(2),dimarr(soildim),id(4)))
268        array4 = new_array
269        deallocate(new_array)
270      endif   
271     
272      dimids(k) = soildim
273      dimid(k) = soildim
274      id(k) = dimarr(dimids(k))
275      edges(k) = dimarr(dimids(k))
276
277      exit
278    endif
279  enddo
280 
281  ! Allocate new arrays
282  id(lonid) = nlonnew
283  id(latid) = nlatnew
284  if (vdims == 2) allocate(array2_new(id(1),id(2)))
285  if (vdims == 3) allocate(array3_new(id(1),id(2),id(3)))
286  if (vdims == 4) allocate(array4_new(id(1),id(2),id(3),id(4)))
287 
288  ! Transform onto new lat-lon grid (loops assume array dimensions ordered lon, lat, sig, time)
289  if (vdims == 2) then
290    do k = 1, nlatnew
291      array2_new(1:nlonnew-1,k) = weight(k)*array2(:,loc(k)) + (1-weight(k))*array2(:,loc(k)+1)
292    enddo
293    array2_new(nlonnew,:) = array2_new(1,:)
294    array2_new(:,1) = sum(array2_new(:,2))/nlonnew
295    array2_new(:,nlatnew) = sum(array2_new(:,nlatnew-1))/nlonnew
296    if (minval(array2) > 0.0) then
297      where (array2_new < 0.0) array2_new = 0.0
298    endif
299  else if (vdims == 3) then
300    do k = 1, nlatnew
301      array3_new(1:nlonnew-1,k,:) = weight(k)*array3(:,loc(k),:) + (1-weight(k))*array3(:,loc(k)+1,:)
302    enddo
303    array3_new(nlonnew,:,:) = array3_new(1,:,:)
304    do k = 1, id(3)
305      array3_new(:,1,k) = sum(array3_new(:,2,k))/nlonnew
306      array3_new(:,nlatnew,k) = sum(array3_new(:,nlatnew-1,k))/nlonnew
307    enddo
308    if (minval(array3) > 0.0) then
309      where (array3_new < 0.0) array3_new = 0.0
310    endif
311  else if (vdims == 4) then
312    do k = 1, nlatnew
313      array4_new(1:nlonnew-1,k,:,:) = weight(k)*array4(:,loc(k),:,:) + (1-weight(k))*array4(:,loc(k)+1,:,:)
314    enddo
315    array4_new(nlonnew,:,:,:) = array4_new(1,:,:,:)
316    do k = 1, id(4)
317      do j = 1, id(3)
318        array4_new(:,1,j,k) = sum(array4_new(:,2,j,k))/nlonnew
319        array4_new(:,nlatnew,j,k) = sum(array4_new(:,nlatnew-1,j,k))/nlonnew
320      enddo
321    enddo
322    if (minval(array4) > 0.0) then
323      where (array4_new < 0.0) array4_new = 0.0
324    endif
325  else
326    print*, ' > '//trim(vname)//': code does not interpolate 1D variables'
327  endif
328 
329  ! Write variable to new file
330  ierr = nf_redef(nid2)
331  ierr = nf_def_var(nid2,trim(vname),nf_float,vdims,dimid,npos)
332 
333  ierr = nf_enddef(nid2)
334  if (vdims == 2) then
335    ierr = nf_put_var_real(nid2,npos,array2_new)
336  else
337    do k = 1, id(vdims)
338      corners(vdims) = k
339      if (vdims == 3) ierr = nf_put_vara_real(nid2,npos,corners,edges,array3_new(:,:,k))
340      if (vdims == 4) ierr = nf_put_vara_real(nid2,npos,corners,edges,array4_new(:,:,:,k))
341    enddo
342    if (vdims == 3) deallocate(array3,array3_new)
343    if (vdims == 4) deallocate(array4,array4_new)
344  endif
345 
346  if (ierr /= nf_noerr) then
347    print*, "Oh no! Problem writing "//trim(vname)
348    stop
349  endif
350 
351  deallocate(id,corners,edges,dimid)
352  npos = npos + 1
353 
354enddo
355
356! Add geopotential height to file
357if (wphis) then
358  print*, 'Adding phisinit to file'
359  allocate(id(2))
360
361  if (nlatnew-1 == 8)   res = '5'
362  if (nlatnew-1 == 14)  res = '10'
363  if (nlatnew-1 == 24)  res = '21'
364  if (nlatnew-1 == 36)  res = '31'
365  if (nlatnew-1 == 48)  res = '42'
366  if (nlatnew-1 == 72)  res = '63'
367  if (nlatnew-1 == 96)  res = '85'
368  if (nlatnew-1 == 144) res = '127'
369  if (nlatnew-1 == 192) res = '170'
370  allocate(geopot(nlonnew-1,nlatnew-1))
371  allocate(geopot_new(nlonnew,nlatnew))
372  ierr = nf_open("phisinit_files/phisinit-T"//trim(res)//".nc",nf_nowrite,nid3)
373  ierr = nf_inq_varid(nid3,"phis",phisid)
374  if (ierr /= nf_noerr) then
375    print*, "Oh no! Can't find geopotential file. Please link from the ukv_phisinit folder:"
376    print*, "ln -sf $MMM/your_comp_dir/PREP_MARS/ukv_phisinit/phisinit-T"//trim(res)//".nc ."
377    stop
378  endif
379  ierr = nf_get_var_real(nid3,phisid,geopot)
380  ierr = nf_close(nid3)
381 
382  do k = 1, nlatnew
383    geopot_new(1:nlonnew-1,k) = weight(k)*geopot(:,loc(k)) + (1-weight(k))*geopot(:,loc(k)+1)
384  enddo
385  geopot_new(nlonnew,:) = geopot_new(1,:)
386  geopot_new(:,1) = sum(geopot_new(:,2))/nlonnew
387  geopot_new(:,nlatnew) = sum(geopot_new(:,nlatnew-1))/nlonnew
388 
389  id(1) = lonid
390  id(2) = latid
391  ierr = nf_redef(nid2)
392  ierr = nf_def_var(nid2,"phisinit",nf_float,2,id,npos)
393  ierr = nf_enddef(nid2)
394  ierr = nf_put_var_real(nid2,npos,geopot_new)
395  if (ierr /= nf_noerr) then
396    print*, "Oh no! Problem writing phisinit"
397    stop
398  endif
399endif
400
401print*, 'Closing files...'
402ierr = nf_close(nid)
403ierr = nf_close(nid2)
404
405end
Note: See TracBrowser for help on using the repository browser.