source: trunk/LMDZ.GENERIC/utilities/extrapol_icefield.F90 @ 3537

Last change on this file since 3537 was 1513, checked in by mturbet, 9 years ago

New routine to extrapolate ice fields

File size: 14.5 KB
Line 
1PROGRAM extrapol_icefield
2
3!--------------------------------------------------------------------------------------------------
4! This program is a tool to accelerate the calculation of ice fieds evolution.
5!
6! It uses data files (diagfi.nc) to extrapolate surface
7! physical fields (ice fields typically) in time.
8!
9! 1. We load data file(s) 'diagfi.nc' and get dimensions (longitude,latitude,altitude,time).
10! 2. We get a surface field from the loaded 'diagfi.nc' file.
11! 3. We make the extrapolation.
12! 4. We load a start file 'startfi.nc' and copy it into a new start file 'startfi_extrapolated.nc'.
13! 5. We modify the 'startfi_extrapolated.nc' according to the extrapolation calculations.
14!
15! --> 'startfi_extrapolated.nc' is the interpolated new start file that can be used to run
16!
17!
18! Author : M. Turbet (2016) [Adapted from E. Millour previous work]
19!
20!
21!--------------------------------------------------------------------------------------------------
22
23use netcdf
24 
25implicit none
26
27integer :: j,ig,t,flag
28
29! NetCDF variables
30integer :: status_d,inid,varid_d
31integer :: lat_dimid,lon_dimid,alt_dimid,time_dimid
32integer :: status_s,fileid,varid_s
33integer :: physical_points_dimid
34
35character(len=128) :: diag_file=""
36character(len=128) :: start_file_in=""
37character(len=128) :: start_file_out="startfi_extrapolated.nc"
38character(len=128) :: varname_d=""
39character(len=128) :: varname_s=""
40character(len=128) :: varname_area="area"
41 
42real,dimension(:),allocatable :: longitude,latitude,altitude,time ! # dimensions
43integer lonlen,latlen,altlen,timelen ! # of grid points along longitude/latitude/altitude/time
44
45integer :: physical_points ! number of atmospheric columns (created)
46integer :: physical_points_s ! number of atmospheric columns (extracted from startfi file)
47
48real,dimension(:,:,:),allocatable :: field_3D ! LON x LAT x TIME field
49real,dimension(:,:),allocatable :: field_phys_3D ! Physical_points x Time field
50real,dimension(:),allocatable :: field_phys_3D_reduced ! Physical_points field
51real,allocatable :: surf_field_3D(:) ! Physical_points field
52real,allocatable :: area(:) ! to store the 1D field of the area
53
54real :: field_ini
55real :: field_fin
56
57integer :: nb_field
58integer :: datashape(3)
59
60integer :: n_years
61
62integer :: extrapolation_mode
63
64write(*,*) ""
65write(*,*) "Welcome in extrapol_icefield routine !"
66write(*,*) "This tool was designed to extrapolate ice field evolution."
67write(*,*) ""
68
69
70!===============================================================================
71!
72! I. INPUT FILES
73!
74!===============================================================================
75
76write(*,*) ""
77write(*,*) "Enter nectdf (diagfi.nc file type) data file name:"
78read(*,'(a128)') diag_file
79write(*,*) ""
80 
81write(*,*) ""
82write(*,*) "Enter nectdf (startfi.nc file type) input file name:"
83read(*,'(a128)') start_file_in
84write(*,*) ""
85
86write(*,*) "Output file name is now ", start_file_out
87write(*,*) ""
88
89! We copy the startfi (input) file and work now with the startfi (output) file.
90CALL system ("cp "//trim(adjustl(start_file_in))//" "//trim(adjustl(start_file_out)))
91
92 
93!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
94!
95! I.1. Open diagfi.nc file
96!
97!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
98
99
100status_d=nf90_open(diag_file,NF90_NOWRITE,inid)
101if (status_d.ne.nf90_noerr) then
102  write(*,*)"Failed to open datafile ",trim(diag_file)
103  write(*,*)trim(nf90_strerror(status_d))
104 stop
105endif
106
107
108!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
109!
110! I.2. Get grids for dimensions longitude,latitude,altitude and time
111!
112!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
113
114
115! Get latitudes.
116status_d=nf90_inq_dimid(inid,"latitude",lat_dimid)
117if (status_d.ne.nf90_noerr) then
118  write(*,*)"Failed to find latitude dimension"
119  write(*,*)trim(nf90_strerror(status_d))
120 stop
121endif
122status_d=nf90_inquire_dimension(inid,lat_dimid,len=latlen)
123if (status_d.ne.nf90_noerr) then
124  write(*,*)"Failed to find latitude length"
125  write(*,*)trim(nf90_strerror(status_d))
126endif
127allocate(latitude(latlen))
128status_d=nf90_inq_varid(inid,"latitude",varid_d)
129if (status_d.ne.nf90_noerr) then
130  write(*,*) "Failed to find latitude ID"
131  write(*,*)trim(nf90_strerror(status_d))
132  stop
133endif
134status_d=nf90_get_var(inid,varid_d,latitude)
135if (status_d.ne.nf90_noerr) then
136  write(*,*) "Failed to load latitude"
137  write(*,*)trim(nf90_strerror(status_d))
138  stop
139endif
140
141! Get longitudes.
142status_d=nf90_inq_dimid(inid,"longitude",lon_dimid)
143if (status_d.ne.nf90_noerr) then
144  write(*,*)"Failed to find longitude dimension"
145  write(*,*)trim(nf90_strerror(status_d))
146 stop
147endif
148status_d=nf90_inquire_dimension(inid,lon_dimid,len=lonlen)
149if (status_d.ne.nf90_noerr) then
150  write(*,*)"Failed to find longitude length"
151  write(*,*)trim(nf90_strerror(status_d))
152endif
153allocate(longitude(lonlen))
154status_d=nf90_inq_varid(inid,"longitude",varid_d)
155if (status_d.ne.nf90_noerr) then
156  write(*,*) "Failed to find longitude ID"
157  write(*,*)trim(nf90_strerror(status_d))
158  stop
159endif
160status_d=nf90_get_var(inid,varid_d,longitude)
161if (status_d.ne.nf90_noerr) then
162  write(*,*) "Failed to load longitude"
163  write(*,*)trim(nf90_strerror(status_d))
164  stop
165endif
166
167! Get times.
168status_d=nf90_inq_dimid(inid,"Time",time_dimid)
169if (status_d.ne.nf90_noerr) then
170  write(*,*)"Failed to find Time dimension"
171  write(*,*)trim(nf90_strerror(status_d))
172 stop
173endif
174status_d=nf90_inquire_dimension(inid,time_dimid,len=timelen)
175if (status_d.ne.nf90_noerr) then
176  write(*,*)"Failed to find Time length"
177  write(*,*)trim(nf90_strerror(status_d))
178endif
179allocate(time(timelen))
180status_d=nf90_inq_varid(inid,"Time",varid_d)
181if (status_d.ne.nf90_noerr) then
182  write(*,*) "Failed to find Time ID"
183  write(*,*)trim(nf90_strerror(status_d))
184  stop
185endif
186status_d=nf90_get_var(inid,varid_d,time)
187if (status_d.ne.nf90_noerr) then
188  write(*,*) "Failed to load Time"
189  write(*,*)trim(nf90_strerror(status_d))
190  stop
191endif
192
193! Get altitudes.
194status_d=nf90_inq_dimid(inid,"altitude",alt_dimid)
195if (status_d.ne.nf90_noerr) then
196  write(*,*)"Failed to find altitude dimension"
197  write(*,*)trim(nf90_strerror(status_d))
198 stop
199endif
200status_d=nf90_inquire_dimension(inid,alt_dimid,len=altlen)
201if (status_d.ne.nf90_noerr) then
202  write(*,*)"Failed to find altitude length"
203  write(*,*)trim(nf90_strerror(status_d))
204endif
205allocate(altitude(altlen))
206status_d=nf90_inq_varid(inid,"altitude",varid_d)
207if (status_d.ne.nf90_noerr) then
208  write(*,*) "Failed to find altitude ID"
209  write(*,*)trim(nf90_strerror(status_d))
210  stop
211endif
212status_d=nf90_get_var(inid,varid_d,altitude)
213if (status_d.ne.nf90_noerr) then
214  write(*,*) "Failed to load altitude"
215  write(*,*)trim(nf90_strerror(status_d))
216  stop
217endif
218
219
220!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
221!
222! I.3. Get the data field.
223!
224!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
225
226
227 write(*,*) ""
228 write(*,*) "Enter surface tracer field (in diagfi.nc) to extrapolate"
229 read(*,'(a128)') varname_d
230 write(*,*) ""
231 
232 write(*,*) ""
233 write(*,*) "Enter corresponding surface tracer field (in startfi.nc) to extrapolate"
234 read(*,'(a128)') varname_s
235 write(*,*) ""
236
237
238status_d=nf90_inq_varid(inid,trim(varname_d),varid_d)
239if (status_d.ne.nf90_noerr) then
240  write(*,*) "Failed to find variable ",trim(varname_d)," in ",trim(diag_file)
241  write(*,*)trim(nf90_strerror(status_d))
242  write(*,*) " Might as well stop here."
243  stop
244endif
245
246! Sanity checks on the variable.
247status_d=nf90_inquire_variable(inid,varid_d,ndims=nb_field,dimids=datashape)
248if (status_d.ne.nf90_noerr) then
249  write(*,*) "Failed to obtain information on variable ",trim(varname_d)
250  write(*,*)trim(nf90_strerror(status_d))
251  write(*,*) " Might as well stop here."
252  stop
253else
254  ! check that it is a 3D variable.
255  if (nb_field.ne.3) then
256    write(*,*) "Error, expected a 3D (lon-lat-time) variable!"
257    stop
258  endif
259  ! check that its dimensions are indeed lon,lat,time.
260  if (datashape(1).ne.lon_dimid) then
261    write(*,*) "Error, expected first dimension to be longitude!"
262    write(*,*) "Please check your diagfi.nc file"
263    stop
264  endif
265  if (datashape(2).ne.lat_dimid) then
266    write(*,*) "Error, expected second dimension to be latitude!"
267    write(*,*) "Please check your diagfi.nc file"
268    stop
269  endif
270  if (datashape(3).ne.time_dimid) then
271    write(*,*) "Error, expected third dimension to be time!"
272    write(*,*) "Please check your diagfi.nc file"
273    stop
274  endif
275endif
276
277allocate(field_3D(lonlen,latlen,timelen))
278
279! Load the data field.
280status_d=nf90_get_var(inid,varid_d,field_3D)
281if (status_d.ne.nf90_noerr) then
282  write(*,*) "Failed to load ",trim(varname_d)
283  write(*,*)trim(nf90_strerror(status_d))
284  stop
285else
286  write(*,*) "Loaded ",trim(varname_d)
287endif
288
289
290!===============================================================================
291!
292! II. MANIPULATION OF DATA FIELDS
293!
294!===============================================================================
295
296
297!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
298!
299! II.1. Transfer the 2D data onto the "physics" 1-dimensional grid.
300!
301!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
302
303! Calculate the physical 1D grid.
304physical_points=2+(latlen-2)*(lonlen-1)
305write(*,*) " lonlen=",lonlen," latlen=",latlen
306write(*,*) " physical_points=",physical_points
307
308allocate(field_phys_3D(physical_points,timelen))
309allocate(field_phys_3D_reduced(physical_points))
310
311! North pole :
312field_phys_3D(1,1:timelen)=field_3D(1,1,1:timelen)
313! South pole :
314field_phys_3D(physical_points,1:timelen)=field_3D(lonlen,latlen,1:timelen)
315! Rest of the world :
316do j=2,latlen-1
317  ig=2+(j-2)*(lonlen-1)
318  field_phys_3D(ig:ig+(lonlen-2),1:timelen)=field_3D(1:lonlen-1,j,1:timelen)
319enddo
320
321field_phys_3D_reduced(1:physical_points)=0.
322
323!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
324!
325! II.2. Extrapolation of data
326!
327!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
328
329
330 write(*,*) ""
331 write(*,*) "Extrapolation/Speed Factor ?"
332 read(*,*) n_years
333 write(*,*) ""
334 
335 write(*,*) ""
336 write(*,*) "Extrapolation Modes :"
337 write(*,*) "Mode 1 : We extrapolate using only the first value and the last value."
338 write(*,*) "Mode 2 : We extrapolate using the first value and the mean value."
339 write(*,*) ""
340 write(*,*) "Write the number of the mode you want :"
341 read(*,*) extrapolation_mode
342 write(*,*) ""
343
344
345if (extrapolation_mode .eq. 1) then
346
347   do ig = 1, physical_points
348      flag = 0
349      do t = 1,timelen
350         if (field_phys_3D(ig,t) .eq. 0) flag = 1 ! Here, we remove the "seasonal" ice deposit.
351      enddo
352      if (flag .eq. 0) then
353         field_phys_3D_reduced(ig)=(field_phys_3D(ig,timelen)-field_phys_3D(ig,1))*n_years
354         if(field_phys_3D_reduced(ig) .lt. 0) field_phys_3D_reduced(ig) = 0   
355      else
356         field_phys_3D_reduced(ig)=0. ! Here, we remove the "seasonal" ice deposit.
357      endif
358   enddo
359   
360else if (extrapolation_mode .eq. 2) then
361
362   do ig = 1, physical_points
363      flag = 0
364      do t = 1,timelen
365         if (field_phys_3D(ig,t) .eq. 0) flag = 1 ! Here, we remove the "seasonal" ice deposit.
366      enddo
367      if (flag .eq. 0) then
368         do t = 1,timelen
369            field_phys_3D_reduced(ig)=field_phys_3D_reduced(ig)+field_phys_3D(ig,t)
370         enddo
371         field_phys_3D_reduced(ig)=field_phys_3D_reduced(ig)/timelen
372         field_phys_3D_reduced(ig)=(field_phys_3D_reduced(ig)-field_phys_3D(ig,1))*2*n_years
373         if(field_phys_3D_reduced(ig) .lt. 0) field_phys_3D_reduced(ig) = 0   
374      else
375         field_phys_3D_reduced(ig)=0.
376      endif
377   enddo
378else
379    write(*,*) "Please check your extrapolation mode because "
380    write(*,*) "it does not match the possibilities !"
381    stop
382endif
383
384
385!===============================================================================
386!
387! III. CREATION OF THE EXTRAPOLATED STARTFI FILE
388!
389!===============================================================================
390
391
392! Open output file in read/write mode.
393status_s=nf90_open(start_file_out,NF90_WRITE,fileid)
394if (status_s.ne.nf90_noerr) then
395  write(*,*)"Failed to open output file ",trim(start_file_out)
396  write(*,*)trim(nf90_strerror(status_s))
397 stop
398endif
399write(*,*) "Reading output file: ",trim(start_file_out)
400
401! Get value of physical_points.
402status_s=nf90_inq_dimid(fileid,"physical_points",physical_points_dimid)
403if (status_s.ne.nf90_noerr) then
404  write(*,*)"Failed to find physical_points dimension"
405  write(*,*)trim(nf90_strerror(status_s))
406 stop
407endif
408
409status_s=nf90_inquire_dimension(fileid,physical_points_dimid,len=physical_points_s)
410if (status_s.ne.nf90_noerr) then
411  write(*,*)"Failed to find physical_points value"
412  write(*,*)trim(nf90_strerror(status_s))
413 stop
414else
415  write(*,*) " physical_points = ",physical_points_s
416endif
417
418if (physical_points_s .ne. physical_points) then
419   write(*,*) "the number of physical points in diagfi and startfi are not equal!"
420   stop
421endif
422
423! Load surface field
424allocate(surf_field_3D(physical_points))
425allocate(area(physical_points))
426
427! Get id of the area
428status_s=nf90_inq_varid(fileid,varname_area,varid_s)
429if (status_s.ne.nf90_noerr) then
430  write(*,*)"Failed to find ",trim(varname_area)," varid_s"
431  write(*,*)"Check that the area is included in your files"
432  write(*,*)trim(nf90_strerror(status_s))
433  stop
434endif
435
436! Read variable from file
437status_s=nf90_get_var(fileid,varid_s,area)
438if (status_s.ne.nf90_noerr) then
439  write(*,*)"Failed to load ",trim(varname_area)
440  write(*,*)trim(nf90_strerror(status_s))
441  stop
442else
443  write(*,*)"Loaded ",trim(varname_s)
444endif
445
446! Get id of variable
447status_s=nf90_inq_varid(fileid,varname_s,varid_s)
448if (status_s.ne.nf90_noerr) then
449  write(*,*)"Failed to find ",trim(varname_s)," varid_s"
450  write(*,*)trim(nf90_strerror(status_s))
451  stop
452endif
453
454! Read variable from file
455status_s=nf90_get_var(fileid,varid_s,surf_field_3D)
456if (status_s.ne.nf90_noerr) then
457  write(*,*)"Failed to load ",trim(varname_s)
458  write(*,*)trim(nf90_strerror(status_s))
459  stop
460else
461  write(*,*)"Loaded ",trim(varname_s)
462endif
463
464
465! Here, we calculate the final ice field.
466! In particular, we do some calculations to conserve the total amount of ice.
467
468field_ini=0.
469do ig=1,physical_points
470   field_ini = field_ini + area(ig)*surf_field_3D(ig)
471enddo
472
473surf_field_3D(1:physical_points)= surf_field_3D(1:physical_points) &
474   + field_phys_3D_reduced(1:physical_points)
475
476field_fin=0.
477do ig=1,physical_points
478   field_fin = field_fin + area(ig)*surf_field_3D(ig)
479enddo
480
481do ig=1, physical_points
482   surf_field_3D(ig) = surf_field_3D(ig)*field_ini/field_fin
483enddo
484
485
486! 3. Write the data to the startfi file
487status_s=nf90_put_var(fileid,varid_s,surf_field_3D)
488if (status_s.ne.nf90_noerr) then
489  write(*,*)"Failed to write ",trim(varname_s)
490  write(*,*)trim(nf90_strerror(status_s))
491  stop
492else
493  write(*,*)"Wrote ",trim(varname_s)
494endif
495
496! 4. Close file
497status_s=nf90_close(fileid)
498
499END PROGRAM extrapol_icefield
Note: See TracBrowser for help on using the repository browser.