1 | MODULE restcom |
---|
2 | !- |
---|
3 | !$Id$ |
---|
4 | !- |
---|
5 | USE netcdf |
---|
6 | !- |
---|
7 | USE errioipsl, ONLY : ipslerr |
---|
8 | USE stringop |
---|
9 | USE calendar |
---|
10 | USE mathelp |
---|
11 | USE fliocom, ONLY : flio_dom_file,flio_dom_att |
---|
12 | !- |
---|
13 | IMPLICIT NONE |
---|
14 | !- |
---|
15 | PRIVATE |
---|
16 | !- |
---|
17 | PUBLIC :: & |
---|
18 | & restini, restget, restput, restclo, & |
---|
19 | & ioconf_setatt, ioget_vname, ioconf_expval, & |
---|
20 | & ioget_expval, ioget_vdim |
---|
21 | !- |
---|
22 | INTERFACE restput |
---|
23 | MODULE PROCEDURE & |
---|
24 | & restput_r3d, restput_r2d, restput_r1d, & |
---|
25 | & restput_opp_r2d, restput_opp_r1d |
---|
26 | END INTERFACE |
---|
27 | !- |
---|
28 | INTERFACE restget |
---|
29 | MODULE PROCEDURE & |
---|
30 | & restget_r3d,restget_r2d,restget_r1d, & |
---|
31 | & restget_opp_r2d,restget_opp_r1d |
---|
32 | END INTERFACE |
---|
33 | !- |
---|
34 | ! We do not use allocatable arrays because these sizes are safe |
---|
35 | ! and we do not know from start how many variables will be in |
---|
36 | ! the out file. |
---|
37 | !- |
---|
38 | INTEGER,PARAMETER :: & |
---|
39 | & max_var=500, max_file=50, max_dim=NF90_MAX_VAR_DIMS |
---|
40 | !- |
---|
41 | CHARACTER(LEN=9),SAVE :: calend_str='unknown' |
---|
42 | !- |
---|
43 | ! The IDs of the netCDF files are going in pairs. |
---|
44 | ! The input one (netcdf_id(?,1)) and the output one (netcdf_id(?,2)) |
---|
45 | !- |
---|
46 | INTEGER,SAVE :: nb_fi = 0 |
---|
47 | INTEGER,DIMENSION(max_file,2),SAVE :: netcdf_id = -1 |
---|
48 | !- |
---|
49 | ! Description of the content of the 'in' files and the 'out' files. |
---|
50 | ! Number of variables : nbvar_* |
---|
51 | ! Number of dimensions : nbdim_* |
---|
52 | ! ID of the time axis : tdimid_* |
---|
53 | !- |
---|
54 | INTEGER,SAVE :: nbvar_in(max_file), nbvar_out(max_file) |
---|
55 | INTEGER,SAVE :: tdimid_in(max_file), tdimid_out(max_file) |
---|
56 | !- |
---|
57 | ! Variables for one or the other file |
---|
58 | !- |
---|
59 | ! Number of dimensions in the input file : nbdim_in |
---|
60 | ! Number of variables read so far from the input file : nbvar_read |
---|
61 | ! Type of variable read from the input file : vartyp_in |
---|
62 | ! (Could be used later to test if we have a restart file) |
---|
63 | !- |
---|
64 | INTEGER,SAVE :: nbdim_in(max_file), nbvar_read(max_file) |
---|
65 | INTEGER,SAVE :: vartyp_in(max_file, max_var) |
---|
66 | !- |
---|
67 | ! Time step and time origine in the input file. |
---|
68 | !- |
---|
69 | REAL,DIMENSION(max_file),SAVE :: deltat,timeorig |
---|
70 | !- |
---|
71 | ! Description of the axes in the output file |
---|
72 | !- |
---|
73 | ! tstp_out : Index on the tie axis currently beeing written |
---|
74 | ! itau_out : Time step which is written on this index of the file |
---|
75 | !- |
---|
76 | INTEGER,DIMENSION(max_file),SAVE :: tstp_out,itau_out |
---|
77 | !- |
---|
78 | ! Description of the axes in the output file |
---|
79 | !- |
---|
80 | ! For the ?ax_infs variable the following order is used : |
---|
81 | ! ?ax_infs (if,in,1) = size of axis |
---|
82 | ! ?ax_infs (if,in,2) = id of dimension |
---|
83 | ! Number of x,y and z axes in the output file : |
---|
84 | ! ?ax_nb(if) |
---|
85 | !- |
---|
86 | INTEGER,DIMENSION(max_file,max_dim,2),SAVE :: & |
---|
87 | & xax_infs,yax_infs,zax_infs |
---|
88 | INTEGER,DIMENSION(max_file),SAVE :: & |
---|
89 | & xax_nb=0,yax_nb=0,zax_nb=0 |
---|
90 | !- |
---|
91 | ! Description of the time axes in the input and output files |
---|
92 | !- |
---|
93 | ! ID of the variable which contains the itaus : |
---|
94 | ! tind_varid_* |
---|
95 | ! ID of the variables which contains the seconds since date : |
---|
96 | ! tax_varid_* |
---|
97 | ! Size of the time axis in the input file : |
---|
98 | ! tax_size_in |
---|
99 | !- |
---|
100 | INTEGER,SAVE :: tind_varid_in(max_file), tax_varid_in(max_file), & |
---|
101 | & tind_varid_out(max_file), tax_varid_out(max_file) |
---|
102 | INTEGER,SAVE :: tax_size_in(max_file)=1 |
---|
103 | !- |
---|
104 | ! The two time axes we have in the input file : |
---|
105 | ! t_index : dates in itaus |
---|
106 | ! (thus the variable has a tstep_sec attribute) |
---|
107 | ! t_julian : Julian days of the time axis |
---|
108 | !- |
---|
109 | INTEGER,SAVE,ALLOCATABLE :: t_index(:,:) |
---|
110 | REAL,SAVE,ALLOCATABLE :: t_julian(:,:) |
---|
111 | !- |
---|
112 | ! Here we save a number of informations on the variables |
---|
113 | ! in the files we are handling |
---|
114 | !- |
---|
115 | ! Name of variables : varname_* |
---|
116 | ! ID of the variables : varid_* |
---|
117 | ! Number of dimensions of the variable : varnbdim_* |
---|
118 | ! Dimensions which are used for the variable : vardims_* |
---|
119 | ! Number of attributes for a variables : varatt_* |
---|
120 | ! A flag which markes the variables we have worked on : touched_* |
---|
121 | !- |
---|
122 | CHARACTER(LEN=20),DIMENSION(max_file,max_var),SAVE :: & |
---|
123 | & varname_in,varname_out |
---|
124 | INTEGER,DIMENSION(max_file,max_var),SAVE :: & |
---|
125 | & varid_in,varid_out,varnbdim_in,varatt_in |
---|
126 | INTEGER,DIMENSION(max_file,max_var,max_dim),SAVE :: & |
---|
127 | & vardims_in |
---|
128 | LOGICAL,DIMENSION(max_file,max_var),SAVE :: & |
---|
129 | & touched_in,touched_out |
---|
130 | !- |
---|
131 | CHARACTER(LEN=120),SAVE :: indchfun= 'scatter, fill, gather, coll' |
---|
132 | REAL,PARAMETER :: missing_val=1.e20 |
---|
133 | ! or HUGE(1.0) (maximum real number) |
---|
134 | !- |
---|
135 | ! The default value we will use for variables |
---|
136 | ! which are not present in the restart file |
---|
137 | !- |
---|
138 | REAL,SAVE :: val_exp = 999999. |
---|
139 | LOGICAL,SAVE :: lock_valexp = .FALSE. |
---|
140 | !- |
---|
141 | ! Temporary variables in which we store the attributed which are going |
---|
142 | ! to be given to a new variable which is going to be defined. |
---|
143 | !- |
---|
144 | CHARACTER(LEN=80),SAVE :: rest_units='XXXXX',rest_lname='XXXXX' |
---|
145 | !- |
---|
146 | ! For allocations |
---|
147 | !- |
---|
148 | REAL,ALLOCATABLE,DIMENSION(:),SAVE :: buff_tmp1,buff_tmp2 |
---|
149 | !- |
---|
150 | !=== |
---|
151 | CONTAINS |
---|
152 | !=== |
---|
153 | !- |
---|
154 | SUBROUTINE restini & |
---|
155 | & (fnamein,iim,jjm,lon,lat,llm,lev, & |
---|
156 | & fnameout,itau,date0,dt,fid,owrite_time_in,domain_id) |
---|
157 | !--------------------------------------------------------------------- |
---|
158 | !- This subroutine sets up all the restart process. |
---|
159 | !- It will call the subroutine which opens the input |
---|
160 | !- and output files. |
---|
161 | !- The time step (itau), date of origine (date0) and time step are |
---|
162 | !- READ from the input file. |
---|
163 | !- A file ID, which is common to the input and output file is returned |
---|
164 | !- |
---|
165 | !- If fnamein = fnameout then the same file is used for the reading |
---|
166 | !- the restart conditions and writing the new restart. |
---|
167 | !- |
---|
168 | !- A special mode can be switched in with filename='NONE'. |
---|
169 | !- This means that no restart file is present. |
---|
170 | !- Usefull for creating the first restart file |
---|
171 | !- or to get elements in a file without creating an output file. |
---|
172 | !- |
---|
173 | !- A mode needs to be written in which itau, date0 and dt |
---|
174 | !- are given to the restart process and thus |
---|
175 | !- written into the output restart file. |
---|
176 | !- |
---|
177 | !- INPUT |
---|
178 | !- |
---|
179 | !- fnamein : name of the file for the restart |
---|
180 | !- iim : Dimension in x |
---|
181 | !- jjm : Dimension in y |
---|
182 | !- lon : Longitude in the x,y domain |
---|
183 | !- lat : Latitude in the x,y domain |
---|
184 | !- llm : Dimension in the vertical |
---|
185 | !- lev : Positions of the levels |
---|
186 | !- fnameout : |
---|
187 | !- |
---|
188 | !- OUTPUT |
---|
189 | !- |
---|
190 | !- itau : Time step of the restart file and at which the model |
---|
191 | !- should restart |
---|
192 | !- date0 : Time at which itau = 0 |
---|
193 | !- dt : time step in seconds between two succesiv itaus |
---|
194 | !- fid : File identification of the restart file |
---|
195 | !- |
---|
196 | !- Optional INPUT arguments |
---|
197 | !- |
---|
198 | !- owrite_time_in : logical argument which allows to |
---|
199 | !- overwrite the time in the restart file |
---|
200 | !- domain_id : Domain identifier |
---|
201 | !--------------------------------------------------------------------- |
---|
202 | IMPLICIT NONE |
---|
203 | !- |
---|
204 | CHARACTER(LEN=*),INTENT(IN) :: fnamein,fnameout |
---|
205 | INTEGER :: iim,jjm,llm,fid,itau |
---|
206 | REAL :: lon(iim,jjm),lat(iim,jjm),lev(llm) |
---|
207 | REAL :: date0,dt |
---|
208 | LOGICAL,OPTIONAL :: owrite_time_in |
---|
209 | INTEGER,INTENT(IN),OPTIONAL :: domain_id |
---|
210 | !- |
---|
211 | ! LOCAL |
---|
212 | !- |
---|
213 | INTEGER :: ncfid |
---|
214 | REAL :: dt_tmp,date0_tmp |
---|
215 | LOGICAL :: l_fi,l_fo,l_rw |
---|
216 | LOGICAL :: overwrite_time |
---|
217 | CHARACTER(LEN=120) :: fname |
---|
218 | LOGICAL :: check = .FALSE. |
---|
219 | !--------------------------------------------------------------------- |
---|
220 | !- |
---|
221 | ! 0.0 Prepare the configuration before opening any files |
---|
222 | !- |
---|
223 | IF (.NOT.PRESENT(owrite_time_in)) THEN |
---|
224 | overwrite_time = .FALSE. |
---|
225 | ELSE |
---|
226 | overwrite_time = owrite_time_in |
---|
227 | ENDIF |
---|
228 | !- |
---|
229 | IF (check) THEN |
---|
230 | WRITE(*,*) 'restini 0.0 : ',TRIM(fnamein),' , ',TRIM(fnameout) |
---|
231 | ENDIF |
---|
232 | !- |
---|
233 | nb_fi = nb_fi+1 |
---|
234 | !- |
---|
235 | IF (nb_fi > max_file) THEN |
---|
236 | CALL ipslerr (3,'restini',& |
---|
237 | & 'Too many restart files are used. The problem can be',& |
---|
238 | & 'solved by increasing max_file in restcom.f90 ',& |
---|
239 | & 'and recompiling ioipsl.') |
---|
240 | ENDIF |
---|
241 | !- |
---|
242 | ! 0.1 Define the open flags |
---|
243 | !- |
---|
244 | l_fi = (TRIM(fnamein) /= 'NONE') |
---|
245 | l_fo = (TRIM(fnameout) /= 'NONE') |
---|
246 | IF ((.NOT.l_fi).AND.(.NOT.l_fo)) THEN |
---|
247 | CALL ipslerr (3,'restini',& |
---|
248 | & 'Input and output file names are both to NONE.',& |
---|
249 | & 'It is probably an error.','Verify your logic.') |
---|
250 | ENDIF |
---|
251 | l_rw = l_fi.AND.l_fo.AND.(TRIM(fnamein) == TRIM(fnameout)) |
---|
252 | !- |
---|
253 | IF (check) THEN |
---|
254 | WRITE(*,*) 'restini 0.1 l_fi, l_fo, l_rw : ',l_fi,l_fo,l_rw |
---|
255 | ENDIF |
---|
256 | !- |
---|
257 | ! 1.0 Open the input file. |
---|
258 | !- |
---|
259 | IF (l_fi) THEN |
---|
260 | !--- |
---|
261 | IF (check) WRITE(*,*) 'restini 1.0 : Open input file' |
---|
262 | !-- Add DOMAIN number and ".nc" suffix in file names if needed |
---|
263 | fname = fnamein |
---|
264 | CALL flio_dom_file (fname,domain_id) |
---|
265 | !-- Open the file |
---|
266 | CALL restopenin (nb_fi,fname,l_rw,iim,jjm,lon,lat,llm,lev,ncfid) |
---|
267 | netcdf_id(nb_fi,1) = ncfid |
---|
268 | !--- |
---|
269 | !-- 1.3 Extract the time information |
---|
270 | !--- |
---|
271 | CALL restsett (dt_tmp,date0_tmp,itau,overwrite_time) |
---|
272 | IF (.NOT.overwrite_time) THEN |
---|
273 | dt = dt_tmp |
---|
274 | date0 = date0_tmp |
---|
275 | ENDIF |
---|
276 | !--- |
---|
277 | ELSE |
---|
278 | !--- |
---|
279 | !-- 2.0 The case of a missing restart file is dealt with |
---|
280 | !--- |
---|
281 | IF (check) WRITE(*,*) 'restini 2.0' |
---|
282 | !--- |
---|
283 | IF ( (ALL(MINLOC(lon(:iim,:jjm)) == MAXLOC(lon(:iim,:jjm)))) & |
---|
284 | .AND.(iim > 1) ) THEN |
---|
285 | CALL ipslerr (3,'restini',& |
---|
286 | & 'For creating a restart file the longitudes of the',& |
---|
287 | & 'grid need to be provided to restini. This ',& |
---|
288 | & 'information is needed for the restart files') |
---|
289 | ENDIF |
---|
290 | IF ( (ALL(MINLOC(lat(:iim,:jjm)) == MAXLOC(lat(:iim,:jjm)))) & |
---|
291 | .AND.(jjm > 1) ) THEN |
---|
292 | CALL ipslerr (3,'restini',& |
---|
293 | & 'For creating a restart file the latitudes of the',& |
---|
294 | & 'grid need to be provided to restini. This ',& |
---|
295 | & 'information is needed for the restart files') |
---|
296 | ENDIF |
---|
297 | IF ( (ALL(MINLOC(lev(:llm)) == MAXLOC(lev(:llm)))) & |
---|
298 | .AND.(llm > 1) ) THEN |
---|
299 | CALL ipslerr (3,'restini',& |
---|
300 | & 'For creating a restart file the levels of the',& |
---|
301 | & 'grid need to be provided to restini. This',& |
---|
302 | & 'information is needed for the restart files') |
---|
303 | ENDIF |
---|
304 | !--- |
---|
305 | !-- 2.2 Allocate the time axes and write the inputed variables |
---|
306 | !--- |
---|
307 | tax_size_in(nb_fi) = 1 |
---|
308 | CALL rest_atim (check,'restini') |
---|
309 | t_index(nb_fi,1) = itau |
---|
310 | t_julian(nb_fi,1) = date0 |
---|
311 | ENDIF |
---|
312 | !- |
---|
313 | IF (l_fo.AND.(.NOT.l_rw)) THEN |
---|
314 | !-- Add DOMAIN number and ".nc" suffix in file names if needed |
---|
315 | fname = fnameout |
---|
316 | CALL flio_dom_file (fname,domain_id) |
---|
317 | !-- Open the file |
---|
318 | CALL restopenout & |
---|
319 | (nb_fi,fname,iim,jjm,lon,lat,llm,lev,dt,date0,ncfid,domain_id) |
---|
320 | netcdf_id(nb_fi,2) = ncfid |
---|
321 | ELSE IF (l_fi.AND.l_fo) THEN |
---|
322 | netcdf_id(nb_fi,2) = netcdf_id(nb_fi,1) |
---|
323 | varname_out(nb_fi,:) = varname_in(nb_fi,:) |
---|
324 | nbvar_out(nb_fi) = nbvar_in(nb_fi) |
---|
325 | tind_varid_out(nb_fi) = tind_varid_in(nb_fi) |
---|
326 | tax_varid_out(nb_fi) = tax_varid_in(nb_fi) |
---|
327 | varid_out(nb_fi,:) = varid_in(nb_fi,:) |
---|
328 | touched_out(nb_fi,:) = .TRUE. |
---|
329 | ENDIF |
---|
330 | !- |
---|
331 | ! 2.3 Set the calendar for the run. |
---|
332 | ! This should not produce any error message if |
---|
333 | ! This does not mean any change in calendar |
---|
334 | ! (to be modified in ioconf_calendar) |
---|
335 | !- |
---|
336 | IF (check) THEN |
---|
337 | WRITE(*,*) 'restini 2.3 : Configure calendar if needed : ', & |
---|
338 | calend_str |
---|
339 | ENDIF |
---|
340 | !- |
---|
341 | IF (INDEX(calend_str,'unknown') < 1) THEN |
---|
342 | CALL ioconf_calendar (calend_str) |
---|
343 | IF (check) THEN |
---|
344 | WRITE(*,*) 'restini 2.3b : new calendar : ',calend_str |
---|
345 | ENDIF |
---|
346 | ENDIF |
---|
347 | !- |
---|
348 | ! Save some data in the module |
---|
349 | !- |
---|
350 | deltat(nb_fi) = dt |
---|
351 | !- |
---|
352 | ! Prepare the variables which will be returned |
---|
353 | !- |
---|
354 | fid = nb_fi |
---|
355 | IF (check) THEN |
---|
356 | WRITE(*,*) 'SIZE of t_index :',SIZE(t_index), & |
---|
357 | SIZE(t_index,dim=1),SIZE(t_index,dim=2) |
---|
358 | WRITE(*,*) 't_index = ',t_index(fid,:) |
---|
359 | ENDIF |
---|
360 | itau = t_index(fid,1) |
---|
361 | !- |
---|
362 | IF (check) WRITE(*,*) 'restini END' |
---|
363 | !--------------------- |
---|
364 | END SUBROUTINE restini |
---|
365 | !=== |
---|
366 | SUBROUTINE restopenin & |
---|
367 | (fid,fname,l_rw,iim,jjm,lon,lat,llm,lev,ncfid) |
---|
368 | !--------------------------------------------------------------------- |
---|
369 | !- Opens the restart file and checks that it belongsd to the model. |
---|
370 | !- This means that the coordinates of the model are compared to the |
---|
371 | !- ones in the file. |
---|
372 | !- |
---|
373 | !- The number and name of variable in the file are exctracted. Also |
---|
374 | !- the time details. |
---|
375 | !--------------------------------------------------------------------- |
---|
376 | IMPLICIT NONE |
---|
377 | !- |
---|
378 | INTEGER,INTENT(IN) :: fid,iim,jjm,llm |
---|
379 | CHARACTER(LEN=*),INTENT(IN) :: fname |
---|
380 | REAL :: lon(iim,jjm),lat(iim,jjm),lev(llm) |
---|
381 | LOGICAL,INTENT(IN) :: l_rw |
---|
382 | INTEGER,INTENT(OUT) :: ncfid |
---|
383 | !- |
---|
384 | ! LOCAL |
---|
385 | !- |
---|
386 | INTEGER,DIMENSION(max_dim) :: var_dims,dimlen |
---|
387 | INTEGER :: nb_dim,nb_var,id_unl,id,iv |
---|
388 | INTEGER :: iread,jread,lread,iret |
---|
389 | INTEGER :: lon_vid,lat_vid |
---|
390 | REAL :: lon_read(iim,jjm),lat_read(iim,jjm) |
---|
391 | REAL :: lev_read(llm) |
---|
392 | REAL :: mdlon,mdlat |
---|
393 | CHARACTER(LEN=80) :: units |
---|
394 | CHARACTER(LEN=NF90_max_name),DIMENSION(max_dim) :: dimname |
---|
395 | LOGICAL :: check = .FALSE. |
---|
396 | !--------------------------------------------------------------------- |
---|
397 | !- |
---|
398 | ! If we reuse the same file for input and output |
---|
399 | ! then we open it in write mode |
---|
400 | !- |
---|
401 | IF (l_rw) THEN; id = NF90_WRITE; ELSE; id = NF90_NOWRITE; ENDIF |
---|
402 | iret = NF90_OPEN(fname,id,ncfid) |
---|
403 | IF (iret /= NF90_NOERR) THEN |
---|
404 | CALL ipslerr (3,'restopenin','Could not open file :',fname,' ') |
---|
405 | ENDIF |
---|
406 | !- |
---|
407 | IF (check) WRITE (*,*) "restopenin 0.0 ",TRIM(fname) |
---|
408 | iret = NF90_INQUIRE(ncfid,nDimensions=nb_dim, & |
---|
409 | & nVariables=nb_var,unlimitedDimId=id_unl) |
---|
410 | tdimid_in(fid) = id_unl |
---|
411 | !- |
---|
412 | IF (nb_dim > max_dim) THEN |
---|
413 | CALL ipslerr (3,'restopenin',& |
---|
414 | & 'More dimensions present in file that can be store',& |
---|
415 | & 'Please increase max_dim in the global variables ',& |
---|
416 | & 'in restcom.F90') |
---|
417 | ENDIF |
---|
418 | IF (nb_var > max_var) THEN |
---|
419 | CALL ipslerr (3,'restopenin',& |
---|
420 | & 'More variables present in file that can be store',& |
---|
421 | & 'Please increase max_var in the global variables ',& |
---|
422 | & 'in restcom.F90') |
---|
423 | ENDIF |
---|
424 | !- |
---|
425 | nbvar_in(fid) = nb_var |
---|
426 | nbdim_in(fid) = nb_dim |
---|
427 | iread = -1; jread = -1; lread = -1; |
---|
428 | DO id=1,nb_dim |
---|
429 | iret = NF90_INQUIRE_DIMENSION(ncfid,id, & |
---|
430 | & len=dimlen(id),name=dimname(id)) |
---|
431 | IF (check) THEN |
---|
432 | WRITE (*,*) "restopenin 0.0 dimname",id,TRIM(dimname(id)) |
---|
433 | ENDIF |
---|
434 | IF (TRIM(dimname(id)) == 'x') THEN |
---|
435 | iread = dimlen(id) |
---|
436 | IF (check) WRITE (*,*) "iread",iread |
---|
437 | ELSE IF (TRIM(dimname(id)) == 'y') THEN |
---|
438 | jread = dimlen(id) |
---|
439 | IF (check) WRITE (*,*) "jread",jread |
---|
440 | ELSE IF (TRIM(dimname(id)) == 'z') THEN |
---|
441 | lread = dimlen(id) |
---|
442 | IF (check) WRITE (*,*) "lread",lread |
---|
443 | ENDIF |
---|
444 | ENDDO |
---|
445 | !- |
---|
446 | IF (id_unl > 0) THEN |
---|
447 | !--- |
---|
448 | !-- 0.1 If we are going to add values to this file |
---|
449 | !-- we need to know where it ends |
---|
450 | !-- We also need to have all the dimensions in the file |
---|
451 | !--- |
---|
452 | IF (l_rw) THEN |
---|
453 | tstp_out(fid) = dimlen(id_unl) |
---|
454 | itau_out(fid) = -1 |
---|
455 | tdimid_out(fid) = tdimid_in(fid) |
---|
456 | IF (check) THEN |
---|
457 | WRITE (*,*) & |
---|
458 | & "restopenin 0.0 unlimited axis dimname", & |
---|
459 | & dimname(id_unl),tstp_out(fid) |
---|
460 | ENDIF |
---|
461 | !----- |
---|
462 | xax_nb(fid) = 0 |
---|
463 | yax_nb(fid) = 0 |
---|
464 | zax_nb(fid) = 0 |
---|
465 | !----- |
---|
466 | DO id=1,nb_dim |
---|
467 | IF (dimname(id)(1:1) == 'x') THEN |
---|
468 | xax_nb(fid) = xax_nb(fid)+1 |
---|
469 | xax_infs(fid,xax_nb(fid),1) = dimlen(id) |
---|
470 | xax_infs(fid,xax_nb(fid),2) = id |
---|
471 | ELSE IF (dimname(id)(1:1) == 'y') THEN |
---|
472 | yax_nb(fid) = yax_nb(fid)+1 |
---|
473 | yax_infs(fid,yax_nb(fid),1) = dimlen(id) |
---|
474 | yax_infs(fid,yax_nb(fid),2) = id |
---|
475 | ELSE IF (dimname(id)(1:1) == 'z') THEN |
---|
476 | zax_nb(fid) = zax_nb(fid)+1 |
---|
477 | zax_infs(fid,zax_nb(fid),1) = dimlen(id) |
---|
478 | zax_infs(fid,zax_nb(fid),2) = id |
---|
479 | ENDIF |
---|
480 | ENDDO |
---|
481 | ENDIF |
---|
482 | ELSE |
---|
483 | !--- |
---|
484 | !-- Still need to find a method for dealing with this |
---|
485 | !--- |
---|
486 | ! CALL ipslerr (3,'restopenin',& |
---|
487 | ! & ' We do not deal yet with files without time axis.',' ',' ') |
---|
488 | ENDIF |
---|
489 | !- |
---|
490 | ! 1.0 First let us check that we have the righ restart file |
---|
491 | !- |
---|
492 | IF ((iread /= iim).OR.(jread /= jjm).OR.(lread /= llm)) THEN |
---|
493 | CALL ipslerr (3,'restopenin',& |
---|
494 | & 'The grid of the restart file does not correspond',& |
---|
495 | & 'to that of the model',' ') |
---|
496 | ENDIF |
---|
497 | !- |
---|
498 | ! 2.0 Get the list of variables |
---|
499 | !- |
---|
500 | IF (check) WRITE(*,*) 'restopenin 1.2' |
---|
501 | !- |
---|
502 | lat_vid = -1 |
---|
503 | lon_vid = -1 |
---|
504 | tind_varid_in(fid) = -1 |
---|
505 | tax_varid_in(fid) = -1 |
---|
506 | !- |
---|
507 | DO iv=1,nb_var |
---|
508 | !--- |
---|
509 | varid_in(fid,iv) = iv |
---|
510 | var_dims(:) = 0 |
---|
511 | iret = NF90_INQUIRE_VARIABLE(ncfid,iv, & |
---|
512 | & name=varname_in(fid,iv),xtype=vartyp_in(fid,iv), & |
---|
513 | & ndims=varnbdim_in(fid,iv),dimids=var_dims, & |
---|
514 | & nAtts=varatt_in(fid,iv)) |
---|
515 | !--- |
---|
516 | DO id=1,varnbdim_in(fid,iv) |
---|
517 | iret = NF90_INQUIRE_DIMENSION & |
---|
518 | & (ncfid,var_dims(id),len=vardims_in(fid,iv,id)) |
---|
519 | ENDDO |
---|
520 | !--- |
---|
521 | !-- 2.1 Read the units of the variable |
---|
522 | !--- |
---|
523 | units='' |
---|
524 | iret = NF90_GET_ATT(ncfid,iv,'units',units) |
---|
525 | CALL strlowercase (units) |
---|
526 | CALL cmpblank (units) |
---|
527 | !--- |
---|
528 | !-- 2.2 Catch the time variables |
---|
529 | !--- |
---|
530 | IF (varnbdim_in(fid,iv) == 1) THEN |
---|
531 | IF ( (INDEX(units,'timesteps since') > 0) & |
---|
532 | .AND.(tind_varid_in(fid) < 0) ) THEN |
---|
533 | tind_varid_in(fid) = iv |
---|
534 | tax_size_in(fid) = vardims_in(fid,iv,1) |
---|
535 | ENDIF |
---|
536 | IF ( (INDEX(units,'seconds since') > 0) & |
---|
537 | .AND.(tax_varid_in(fid) < 0) ) THEN |
---|
538 | tax_varid_in(fid) = iv |
---|
539 | tax_size_in(fid) = vardims_in(fid,iv,1) |
---|
540 | ENDIF |
---|
541 | ENDIF |
---|
542 | !--- |
---|
543 | !-- 2.3 Catch longitude and latitude variables |
---|
544 | !--- |
---|
545 | IF (INDEX(units,'degrees_nort') >= 1) THEN |
---|
546 | lat_vid = iv |
---|
547 | ENDIF |
---|
548 | IF (INDEX(units,'degrees_east') >= 1) THEN |
---|
549 | lon_vid = iv |
---|
550 | ENDIF |
---|
551 | !--- |
---|
552 | ENDDO |
---|
553 | !- |
---|
554 | ! 2.4 None of the variables was yet read |
---|
555 | !- |
---|
556 | nbvar_read(fid) = 0 |
---|
557 | touched_in(fid,:) = .FALSE. |
---|
558 | !- |
---|
559 | ! 3.0 Reading the coordinates from the input restart file |
---|
560 | !- |
---|
561 | lon_read = missing_val |
---|
562 | lat_read = missing_val |
---|
563 | !- |
---|
564 | IF (lon_vid < 0 .OR. lat_vid < 0) THEN |
---|
565 | CALL ipslerr (3,'restopenin',& |
---|
566 | & ' No variables containing longitude or latitude were ',& |
---|
567 | & ' found in the restart file.',' ') |
---|
568 | ELSE |
---|
569 | iret = NF90_GET_VAR(ncfid,lon_vid,lon_read) |
---|
570 | iret = NF90_GET_VAR(ncfid,lat_vid,lat_read) |
---|
571 | !--- |
---|
572 | IF ( (ABS( MAXVAL(lon(:,:)) & |
---|
573 | & -MINVAL(lon(:,:))) < EPSILON(MAXVAL(lon(:,:)))) & |
---|
574 | & .AND.(ABS( MAXVAL(lat(:,:)) & |
---|
575 | & -MINVAL(lat(:,:))) < EPSILON(MAXVAL(lat(:,:)))) ) THEN |
---|
576 | !----- |
---|
577 | !---- 3.1 No longitude nor latitude are provided thus |
---|
578 | !---- they are taken from the restart file |
---|
579 | !----- |
---|
580 | lon(:,:) = lon_read(:,:) |
---|
581 | lat(:,:) = lat_read(:,:) |
---|
582 | ELSE |
---|
583 | !----- |
---|
584 | !---- 3.2 We check that the longitudes and latitudes |
---|
585 | !---- in the file and the model are the same |
---|
586 | !----- |
---|
587 | mdlon = MAXVAL(ABS(lon_read-lon)) |
---|
588 | mdlat = MAXVAL(ABS(lat_read-lat)) |
---|
589 | !----- |
---|
590 | !---- We can not test against epsilon here as the longitude |
---|
591 | !---- can be stored at another precision in the netCDF file. |
---|
592 | !---- The test here does not need to be very precise. |
---|
593 | !----- |
---|
594 | IF (mdlon > 1.e-4 .OR. mdlat > 1.e-4) THEN |
---|
595 | CALL ipslerr (3,'restopenin',& |
---|
596 | & ' The longitude or latitude found in the restart ',& |
---|
597 | & ' file are not the same as the ones used in the model.',& |
---|
598 | & ' ') |
---|
599 | ENDIF |
---|
600 | ENDIF |
---|
601 | ENDIF |
---|
602 | !------------------------ |
---|
603 | END SUBROUTINE restopenin |
---|
604 | !=== |
---|
605 | SUBROUTINE restsett (timestep,date0,itau,owrite_time_in) |
---|
606 | !--------------------------------------------------------------------- |
---|
607 | !- Here we get all the time information from the file. |
---|
608 | !- |
---|
609 | !- The time information can come in three forms : |
---|
610 | !- -global attributes which give the time origine and the |
---|
611 | !- time step is taken from the input to restinit |
---|
612 | !- -A physical time exists and thus the julian date from the |
---|
613 | !- input is used for positioning using the itau as input |
---|
614 | !- -A time-step axis exists and itau is positioned on it. |
---|
615 | !- |
---|
616 | !- What takes precedence : the model |
---|
617 | !- |
---|
618 | !- itau : Time step of the model |
---|
619 | !- |
---|
620 | !- Optional INPUT arguments |
---|
621 | !- |
---|
622 | !- owrite_time_in : logical argument which allows to |
---|
623 | !- overwrite the time in the restart file |
---|
624 | !--------------------------------------------------------------------- |
---|
625 | IMPLICIT NONE |
---|
626 | !- |
---|
627 | REAL :: date0,timestep |
---|
628 | INTEGER :: itau |
---|
629 | LOGICAL,OPTIONAL :: owrite_time_in |
---|
630 | !- |
---|
631 | ! LOCAL |
---|
632 | !- |
---|
633 | INTEGER :: ncfid,iret,it,iax,iv |
---|
634 | CHARACTER(LEN=80) :: itau_orig,tax_orig,calendar |
---|
635 | CHARACTER(LEN=9) :: tmp_cal |
---|
636 | INTEGER :: year0,month0,day0,hours0,minutes0,seci |
---|
637 | REAL :: sec0,one_day,one_year,date0_ju,ttmp |
---|
638 | CHARACTER :: strc |
---|
639 | LOGICAL :: ow_time |
---|
640 | !- |
---|
641 | LOGICAL :: check = .FALSE. |
---|
642 | !--------------------------------------------------------------------- |
---|
643 | IF (PRESENT(owrite_time_in)) THEN |
---|
644 | ow_time = owrite_time_in |
---|
645 | ELSE |
---|
646 | ow_time = .FALSE. |
---|
647 | ENDIF |
---|
648 | !- |
---|
649 | ncfid = netcdf_id(nb_fi,1) |
---|
650 | !- |
---|
651 | ! Allocate the space we need for the time axes |
---|
652 | !- |
---|
653 | CALL rest_atim (check,'restsett') |
---|
654 | !- |
---|
655 | ! Get the calendar if possible. Else it will be gregorian. |
---|
656 | !- |
---|
657 | IF (tax_size_in(nb_fi) > 0 ) THEN |
---|
658 | calendar = 'XXXXX' |
---|
659 | iret = NF90_GET_ATT(ncfid,tax_varid_in(nb_fi),'calendar',calendar) |
---|
660 | IF ( INDEX(calendar,'XXXXX') < 0 ) THEN |
---|
661 | CALL ioconf_calendar (calendar) |
---|
662 | IF (check) THEN |
---|
663 | WRITE(*,*) 'restsett : calendar of the restart ',calendar |
---|
664 | ENDIF |
---|
665 | ENDIF |
---|
666 | ENDIF |
---|
667 | CALL ioget_calendar (one_year,one_day) |
---|
668 | IF (check) THEN |
---|
669 | WRITE(*,*) 'one_year,one_day = ',one_year,one_day |
---|
670 | ENDIF |
---|
671 | !- |
---|
672 | itau_orig = 'XXXXX' |
---|
673 | tax_orig = 'XXXXX' |
---|
674 | !- |
---|
675 | ! Get the time steps of the time axis if available on the restart file |
---|
676 | !- |
---|
677 | IF (tind_varid_in(nb_fi) > 0) THEN |
---|
678 | IF (ow_time) THEN |
---|
679 | t_index(nb_fi,:) = itau |
---|
680 | IF (check) THEN |
---|
681 | WRITE(*,*) "nb_fi,t_index",nb_fi,t_index(nb_fi,:) |
---|
682 | ENDIF |
---|
683 | CALL ju2ymds (date0,year0,month0,day0,sec0) |
---|
684 | hours0 = NINT(sec0/3600) |
---|
685 | sec0 = sec0 - 3600 * hours0 |
---|
686 | minutes0 = NINT(sec0 / 60) |
---|
687 | sec0 = sec0 - 60 * minutes0 |
---|
688 | seci = NINT(sec0) |
---|
689 | strc=':' |
---|
690 | IF (check) THEN |
---|
691 | WRITE(*,*) date0 |
---|
692 | WRITE (UNIT=itau_orig,FMT='(I4.4,5(A,I2.2))') & |
---|
693 | & year0,'-',month0,'-',day0,' ',hours0,':',minutes0,':',seci |
---|
694 | WRITE(*,*) "itau_orig : ",itau_orig |
---|
695 | ENDIF |
---|
696 | ELSE |
---|
697 | iret = NF90_GET_VAR(ncfid,tind_varid_in(nb_fi),t_index(nb_fi,:)) |
---|
698 | IF (check) THEN |
---|
699 | WRITE(*,*) "restsett, time axis : ",t_index(nb_fi,:) |
---|
700 | ENDIF |
---|
701 | iret = NF90_GET_ATT(ncfid,tind_varid_in(nb_fi),'units',itau_orig) |
---|
702 | itau_orig = & |
---|
703 | & itau_orig(INDEX(itau_orig,'since')+6:LEN_TRIM(itau_orig)) |
---|
704 | iret = & |
---|
705 | & NF90_GET_ATT(ncfid,tind_varid_in(nb_fi),'tstep_sec',timestep) |
---|
706 | !----- |
---|
707 | !---- This time origin will dominate as it is linked to the time steps. |
---|
708 | !----- |
---|
709 | READ (UNIT=itau_orig,FMT='(I4.4,5(A,I2.2))') & |
---|
710 | & year0,strc,month0,strc,day0,strc, & |
---|
711 | & hours0,strc,minutes0,strc,seci |
---|
712 | sec0 = REAL(seci) |
---|
713 | sec0 = hours0*3600.+minutes0*60.+sec0 |
---|
714 | CALL ymds2ju (year0,month0,day0,sec0,date0) |
---|
715 | ENDIF |
---|
716 | ENDIF |
---|
717 | !- |
---|
718 | ! If a julian day time axis is available then we get it |
---|
719 | !- |
---|
720 | IF (tax_varid_in(nb_fi) > 0) THEN |
---|
721 | iret = NF90_GET_VAR(ncfid,tax_varid_in(nb_fi),t_julian(nb_fi,:)) |
---|
722 | iret = NF90_GET_ATT(ncfid,tax_varid_in(nb_fi),'units',tax_orig) |
---|
723 | tax_orig = tax_orig(INDEX(tax_orig,'since')+6:LEN_TRIM(tax_orig)) |
---|
724 | iret = NF90_GET_ATT(ncfid,tax_varid_in(nb_fi),'calendar',tmp_cal) |
---|
725 | IF (check) THEN |
---|
726 | WRITE(*,*) 'restsett : tmp_calendar of the restart ',tmp_cal |
---|
727 | ENDIF |
---|
728 | !--- |
---|
729 | CALL strlowercase (tmp_cal) |
---|
730 | IF (INDEX(calend_str,tmp_cal) < 0) THEN |
---|
731 | IF (INDEX(calend_str,'unknown') > 0) THEN |
---|
732 | calend_str = tmp_cal |
---|
733 | ELSE |
---|
734 | CALL ipslerr (3,'restsett', & |
---|
735 | & ' In the restart files two different calendars were found.', & |
---|
736 | & ' Please check the files you have used',' ') |
---|
737 | ENDIF |
---|
738 | ENDIF |
---|
739 | !--- |
---|
740 | !-- We need to transform that into julian days |
---|
741 | !-- to get ride of the intial date. |
---|
742 | !--- |
---|
743 | IF (check) WRITE(*,*) 'tax_orig : ',TRIM(tax_orig) |
---|
744 | READ (UNIT=tax_orig,FMT='(I4.4,5(a,I2.2))') & |
---|
745 | year0,strc,month0,strc,day0,strc, & |
---|
746 | hours0,strc,minutes0,strc,seci |
---|
747 | sec0 = REAL(seci) |
---|
748 | sec0 = hours0*3600.+minutes0*60.+sec0 |
---|
749 | CALL ymds2ju (year0,month0,day0,sec0,date0_ju) |
---|
750 | t_julian(nb_fi,:) = t_julian(nb_fi,:)/one_day+date0_ju |
---|
751 | ENDIF |
---|
752 | !- |
---|
753 | IF ( (INDEX(itau_orig,'XXXXX') > 0) & |
---|
754 | .AND.(INDEX(tax_orig,'XXXXX') < 0) ) THEN |
---|
755 | !!- Compute the t_itau from the date read and the timestep in the input |
---|
756 | ENDIF |
---|
757 | !- |
---|
758 | IF ( (INDEX(tax_orig,'XXXXX') > 0) & |
---|
759 | .AND.(INDEX(itau_orig,'XXXXX') < 0) ) THEN |
---|
760 | DO it=1,tax_size_in(nb_fi) |
---|
761 | t_julian(nb_fi,it) = itau2date(t_index(nb_fi,it),date0,timestep) |
---|
762 | ENDDO |
---|
763 | ENDIF |
---|
764 | !- |
---|
765 | ! If neither the indices or time is present then get global attributes |
---|
766 | ! This is for compatibility reasons and should not be used. |
---|
767 | !- |
---|
768 | IF ((tax_varid_in(nb_fi) < 0).AND.(tind_varid_in(nb_fi) < 0)) THEN |
---|
769 | iax = -1 |
---|
770 | DO iv=1,nbvar_in(nb_fi) |
---|
771 | IF ( (INDEX(varname_in(nb_fi,iv),'tsteps') > 0) & |
---|
772 | & .OR.(INDEX(varname_in(nb_fi,iv),'time_steps') > 0)) THEN |
---|
773 | iax = iv |
---|
774 | ENDIF |
---|
775 | ENDDO |
---|
776 | !--- |
---|
777 | IF (iax < 0) THEN |
---|
778 | CALL ipslerr (3,'restsett',& |
---|
779 | & 'No time axis was found in the restart file. Please check',& |
---|
780 | & 'that it corresponds to the convention used in restsett',& |
---|
781 | & ' ') |
---|
782 | ELSE |
---|
783 | iret = NF90_GET_VAR(ncfid,tind_varid_in(nb_fi),t_index(nb_fi,:)) |
---|
784 | iret = NF90_GET_ATT(ncfid,NF90_GLOBAL,'delta_tstep_sec',timestep) |
---|
785 | iret = NF90_GET_ATT(ncfid,NF90_GLOBAL,'year0',ttmp) |
---|
786 | year0 = NINT(ttmp) |
---|
787 | iret = NF90_GET_ATT(ncfid,NF90_GLOBAL,'month0',ttmp) |
---|
788 | month0 = NINT(ttmp) |
---|
789 | iret = NF90_GET_ATT(ncfid,NF90_GLOBAL,'day0',ttmp) |
---|
790 | day0 = NINT(ttmp) |
---|
791 | iret = NF90_GET_ATT(ncfid,NF90_GLOBAL,'sec0',sec0) |
---|
792 | !--- |
---|
793 | CALL ymds2ju (year0,month0,day0,sec0,date0) |
---|
794 | t_julian(nb_fi,1) = itau2date(t_index(nb_fi,1),date0,timestep) |
---|
795 | ENDIF |
---|
796 | ENDIF |
---|
797 | !---------------------- |
---|
798 | END SUBROUTINE restsett |
---|
799 | !=== |
---|
800 | SUBROUTINE restopenout & |
---|
801 | (fid,fname,iim,jjm, & |
---|
802 | lon,lat,llm,lev,timestep,date,ncfid,domain_id) |
---|
803 | !--------------------------------------------------------------------- |
---|
804 | !- Opens the restart file for output. |
---|
805 | !- The longitude and time variables are written. |
---|
806 | !--------------------------------------------------------------------- |
---|
807 | IMPLICIT NONE |
---|
808 | !- |
---|
809 | INTEGER,INTENT(IN) :: fid,iim,jjm,llm |
---|
810 | CHARACTER(LEN=*) :: fname |
---|
811 | REAL :: date,timestep |
---|
812 | REAL :: lon(iim,jjm),lat(iim,jjm),lev(llm) |
---|
813 | INTEGER,INTENT(OUT) :: ncfid |
---|
814 | INTEGER,INTENT(IN),OPTIONAL :: domain_id |
---|
815 | !- |
---|
816 | ! LOCAL |
---|
817 | !- |
---|
818 | INTEGER :: iret |
---|
819 | CHARACTER(LEN=70) :: str_t |
---|
820 | INTEGER :: x_id,y_id,z_id,itauid |
---|
821 | INTEGER :: nlonid,nlatid,nlevid,timeid |
---|
822 | INTEGER :: year,month,day,hours,minutes |
---|
823 | REAL :: sec |
---|
824 | CHARACTER(LEN=3),DIMENSION(12) :: & |
---|
825 | cal = (/'JAN','FEB','MAR','APR','MAY','JUN', & |
---|
826 | 'JUL','AUG','SEP','OCT','NOV','DEC'/) |
---|
827 | CHARACTER(LEN=30) :: timenow,conv |
---|
828 | LOGICAL :: check = .FALSE. |
---|
829 | !--------------------------------------------------------------------- |
---|
830 | IF (check) WRITE(*,*) "restopenout 0.0 ",TRIM(fname) |
---|
831 | !- |
---|
832 | ! If we use the same file for input and output |
---|
833 | !- we will not even call restopenout |
---|
834 | !- |
---|
835 | iret = NF90_CREATE(fname,NF90_NOCLOBBER,ncfid) |
---|
836 | IF (iret == -35) THEN |
---|
837 | CALL ipslerr (3,'restopenout',& |
---|
838 | & ' The restart file aready exists on the disc. IOIPSL ',& |
---|
839 | & ' will not overwrite it. You should remove the old one or ',& |
---|
840 | & ' generate the new one with another name') |
---|
841 | ENDIF |
---|
842 | !- |
---|
843 | iret = NF90_DEF_DIM(ncfid,'x',iim,x_id) |
---|
844 | xax_nb(fid) = xax_nb(fid)+1 |
---|
845 | xax_infs(fid,xax_nb(fid),1) = iim |
---|
846 | xax_infs(fid,xax_nb(fid),2) = x_id |
---|
847 | !- |
---|
848 | iret = NF90_DEF_DIM(ncfid,'y',jjm,y_id) |
---|
849 | yax_nb(fid) = yax_nb(fid)+1 |
---|
850 | yax_infs(fid,yax_nb(fid),1) = jjm |
---|
851 | yax_infs(fid,yax_nb(fid),2) = y_id |
---|
852 | !- |
---|
853 | iret = NF90_DEF_DIM(ncfid,'z',llm,z_id) |
---|
854 | zax_nb(fid) = zax_nb(fid)+1 |
---|
855 | zax_infs(fid,zax_nb(fid),1) = llm |
---|
856 | zax_infs(fid,zax_nb(fid),2) = z_id |
---|
857 | !- |
---|
858 | iret = NF90_DEF_DIM(ncfid,'time',NF90_UNLIMITED,tdimid_out(fid)) |
---|
859 | !- |
---|
860 | ! 1.0 Longitude |
---|
861 | !- |
---|
862 | IF (check) WRITE(*,*) "restopenout 1.0" |
---|
863 | !- |
---|
864 | iret = NF90_DEF_VAR(ncfid,"nav_lon",NF90_FLOAT,(/x_id,y_id/),nlonid) |
---|
865 | iret = NF90_PUT_ATT(ncfid,nlonid,'units',"degrees_east") |
---|
866 | iret = NF90_PUT_ATT(ncfid,nlonid,'valid_min',REAL(-180.,KIND=4)) |
---|
867 | iret = NF90_PUT_ATT(ncfid,nlonid,'valid_max',REAL( 180.,KIND=4)) |
---|
868 | iret = NF90_PUT_ATT(ncfid,nlonid,'long_name',"Longitude") |
---|
869 | !- |
---|
870 | ! 2.0 Latitude |
---|
871 | !- |
---|
872 | IF (check) WRITE(*,*) "restopenout 2.0" |
---|
873 | !- |
---|
874 | iret = NF90_DEF_VAR(ncfid,"nav_lat",NF90_FLOAT,(/x_id,y_id/),nlatid) |
---|
875 | iret = NF90_PUT_ATT(ncfid,nlatid,'units',"degrees_north") |
---|
876 | iret = NF90_PUT_ATT(ncfid,nlatid,'valid_min',REAL(-90.,KIND=4)) |
---|
877 | iret = NF90_PUT_ATT(ncfid,nlatid,'valid_max',REAL( 90.,KIND=4)) |
---|
878 | iret = NF90_PUT_ATT(ncfid,nlatid,'long_name',"Latitude") |
---|
879 | !- |
---|
880 | ! 3.0 Levels |
---|
881 | !- |
---|
882 | IF (check) WRITE(*,*) "restopenout 3.0" |
---|
883 | !- |
---|
884 | iret = NF90_DEF_VAR(ncfid,"nav_lev",NF90_FLOAT,z_id,nlevid) |
---|
885 | iret = NF90_PUT_ATT(ncfid,nlevid,'units',"model_levels") |
---|
886 | iret = NF90_PUT_ATT(ncfid,nlevid,'valid_min', & |
---|
887 | & REAL(MINVAL(lev),KIND=4)) |
---|
888 | iret = NF90_PUT_ATT(ncfid,nlevid,'valid_max', & |
---|
889 | & REAL(MAXVAL(lev),KIND=4)) |
---|
890 | iret = NF90_PUT_ATT(ncfid,nlevid,'long_name',"Model levels") |
---|
891 | !- |
---|
892 | ! 4.0 Time axis, this is the seconds since axis |
---|
893 | !- |
---|
894 | IF (check) WRITE(*,*) "restopenout 4.0" |
---|
895 | !- |
---|
896 | iret = NF90_DEF_VAR(ncfid,"time",NF90_FLOAT, & |
---|
897 | tdimid_out(fid),timeid) |
---|
898 | tax_varid_out(fid) = timeid |
---|
899 | !- |
---|
900 | timeorig(fid) = date |
---|
901 | CALL ju2ymds (date,year,month,day,sec) |
---|
902 | hours = INT(sec/(60.*60.)) |
---|
903 | minutes = INT((sec-hours*60.*60.)/60.) |
---|
904 | sec = sec-(hours*60.*60.+minutes*60.) |
---|
905 | WRITE (UNIT=str_t, & |
---|
906 | FMT='("seconds since ",I4.4,2("-",I2.2)," ",I2.2,2(":",I2.2))') & |
---|
907 | & year,month,day,hours,minutes,INT(sec) |
---|
908 | iret = NF90_PUT_ATT(ncfid,timeid,'units',TRIM(str_t)) |
---|
909 | !- |
---|
910 | CALL ioget_calendar (str_t) |
---|
911 | iret = NF90_PUT_ATT(ncfid,timeid,'calendar',TRIM(str_t)) |
---|
912 | iret = NF90_PUT_ATT(ncfid,timeid,'title','Time') |
---|
913 | iret = NF90_PUT_ATT(ncfid,timeid,'long_name','Time axis') |
---|
914 | !- |
---|
915 | WRITE(UNIT=str_t, & |
---|
916 | FMT='(" ",I4.4,"-",A3,"-",I2.2," ",I2.2,2(":",I2.2))') & |
---|
917 | & year,cal(month),day,hours,minutes,INT(sec) |
---|
918 | iret = NF90_PUT_ATT(ncfid,timeid,'time_origin',TRIM(str_t)) |
---|
919 | !- |
---|
920 | ! 5.0 Time axis, this is the time steps since axis |
---|
921 | !- |
---|
922 | IF (check) WRITE(*,*) "restopenout 5.0" |
---|
923 | !- |
---|
924 | iret = NF90_DEF_VAR(ncfid,"time_steps",NF90_INT, & |
---|
925 | & tdimid_out(fid),itauid) |
---|
926 | tind_varid_out(fid) = itauid |
---|
927 | !- |
---|
928 | CALL ju2ymds (date,year,month,day,sec) |
---|
929 | !- |
---|
930 | hours = INT(sec/(60.*60.)) |
---|
931 | minutes = INT((sec-hours*60.*60.)/60.) |
---|
932 | sec = sec-(hours*60.*60.+minutes*60.) |
---|
933 | !- |
---|
934 | WRITE (UNIT=str_t, & |
---|
935 | FMT='("timesteps since ",I4.4,2("-",I2.2)," ",I2.2,2(":",I2.2))') & |
---|
936 | & year,month,day,hours,minutes,INT(sec) |
---|
937 | !- |
---|
938 | iret = NF90_PUT_ATT(ncfid,itauid,'units',TRIM(str_t)) |
---|
939 | iret = NF90_PUT_ATT(ncfid,itauid,'title','Time steps') |
---|
940 | iret = NF90_PUT_ATT(ncfid,itauid,'tstep_sec',REAL(timestep,KIND=4)) |
---|
941 | iret = NF90_PUT_ATT(ncfid,itauid,'long_name','Time step axis') |
---|
942 | !- |
---|
943 | WRITE(UNIT=str_t, & |
---|
944 | FMT='(" ",I4.4,"-",A3,"-",I2.2," ",I2.2,2(":",I2.2))') & |
---|
945 | & year,cal(month),day,hours,minutes,INT(sec) |
---|
946 | iret = NF90_PUT_ATT(ncfid,itauid,'time_origin',TRIM(str_t)) |
---|
947 | !- |
---|
948 | ! 5.2 Write global attributes |
---|
949 | !- |
---|
950 | conv='GDT 1.2' |
---|
951 | iret = NF90_PUT_ATT(ncfid,NF90_GLOBAL,'Conventions',TRIM(conv)) |
---|
952 | iret = NF90_PUT_ATT(ncfid,NF90_GLOBAL,'file_name',TRIM(fname)) |
---|
953 | !! TO BE DONE LATER |
---|
954 | !! iret = NF90_PUT_ATT(ncfid,NF90_GLOBAL, & |
---|
955 | !! 'production',TRIM(model_name)) |
---|
956 | !! lock_modname = .TRUE. |
---|
957 | CALL ioget_timestamp (timenow) |
---|
958 | iret = NF90_PUT_ATT(ncfid,NF90_GLOBAL,'TimeStamp',TRIM(timenow)) |
---|
959 | !- |
---|
960 | ! Add DOMAIN attributes if needed |
---|
961 | !- |
---|
962 | CALL flio_dom_att (ncfid,domain_id) |
---|
963 | !- |
---|
964 | ! 6.0 The coordinates are written to the file |
---|
965 | !- |
---|
966 | iret = NF90_ENDDEF(ncfid) |
---|
967 | !- |
---|
968 | iret = NF90_PUT_VAR(ncfid,nlonid,lon) |
---|
969 | iret = NF90_PUT_VAR(ncfid,nlatid,lat) |
---|
970 | iret = NF90_PUT_VAR(ncfid,nlevid,lev) |
---|
971 | !- |
---|
972 | ! 7.0 Set a few variables related to the out file |
---|
973 | !- |
---|
974 | nbvar_out(fid) = 0 |
---|
975 | itau_out(fid) = -1 |
---|
976 | tstp_out(fid) = 0 |
---|
977 | touched_out(fid,:) = .FALSE. |
---|
978 | !- |
---|
979 | ! 7.1 The file is put back in define mode. |
---|
980 | ! This will last until itau_out >= 0 |
---|
981 | !- |
---|
982 | iret = NF90_REDEF(ncfid) |
---|
983 | !- |
---|
984 | IF (check) WRITE(*,*) "restopenout END" |
---|
985 | !------------------------- |
---|
986 | END SUBROUTINE restopenout |
---|
987 | !=== |
---|
988 | SUBROUTINE restget_opp_r1d & |
---|
989 | & (fid,vname_q,iim,jjm,llm,itau,def_beha, & |
---|
990 | & var,MY_OPERATOR,nbindex,ijndex) |
---|
991 | !--------------------------------------------------------------------- |
---|
992 | !- This subroutine serves as an interface to restget_real |
---|
993 | !- |
---|
994 | !- Should work as restput_opp_r1d but the other way around ! |
---|
995 | !--------------------------------------------------------------------- |
---|
996 | IMPLICIT NONE |
---|
997 | !- |
---|
998 | INTEGER :: fid |
---|
999 | CHARACTER(LEN=*) :: vname_q |
---|
1000 | INTEGER :: iim,jjm,llm,itau |
---|
1001 | LOGICAL def_beha |
---|
1002 | REAL :: var(:) |
---|
1003 | CHARACTER(LEN=*) :: MY_OPERATOR |
---|
1004 | INTEGER :: nbindex,ijndex(nbindex) |
---|
1005 | !- |
---|
1006 | ! LOCAL |
---|
1007 | !- |
---|
1008 | INTEGER :: req_sz,siz1 |
---|
1009 | REAL :: scal |
---|
1010 | CHARACTER(LEN=7) :: topp |
---|
1011 | LOGICAL :: check = .FALSE. |
---|
1012 | !--------------------------------------------------------------------- |
---|
1013 | !- |
---|
1014 | ! 0.0 What size should be the data in the file |
---|
1015 | !- |
---|
1016 | req_sz = 1 |
---|
1017 | IF (nbindex == iim .AND. jjm <= 1 .AND. llm <= 1) THEN |
---|
1018 | IF (xax_infs(fid,1,1) > 0) req_sz = req_sz*xax_infs(fid,1,1) |
---|
1019 | IF (yax_infs(fid,1,1) > 0) req_sz = req_sz*yax_infs(fid,1,1) |
---|
1020 | IF (zax_infs(fid,1,1) > 0) req_sz = req_sz*zax_infs(fid,1,1) |
---|
1021 | ELSE |
---|
1022 | CALL ipslerr (3,'resget_opp_r1d', & |
---|
1023 | 'Unable to performe an operation on this variable as it has',& |
---|
1024 | 'a second and third dimension',vname_q) |
---|
1025 | ENDIF |
---|
1026 | !- |
---|
1027 | ! 1.0 Allocate the temporary buffer we need |
---|
1028 | ! to put the variable in right dimension |
---|
1029 | !- |
---|
1030 | siz1 = SIZE(var) |
---|
1031 | CALL rest_alloc (1,siz1,check,'restget_opp_r1d') |
---|
1032 | CALL rest_alloc (2,req_sz,check,'restget_opp_r1d') |
---|
1033 | !- |
---|
1034 | ! 2.0 Here we get the variable from the restart file |
---|
1035 | !- |
---|
1036 | CALL restget_real & |
---|
1037 | (fid,vname_q,xax_infs(fid,1,1),yax_infs(fid,1,1), & |
---|
1038 | zax_infs(fid,1,1),itau,def_beha,buff_tmp2) |
---|
1039 | !- |
---|
1040 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1041 | ! into the variable the model expects |
---|
1042 | !- |
---|
1043 | topp = MY_OPERATOR(1:MIN(LEN_TRIM(MY_OPERATOR),7)) |
---|
1044 | !- |
---|
1045 | IF (INDEX(indchfun,topp(:LEN_TRIM(topp))) > 0) THEN |
---|
1046 | scal = missing_val |
---|
1047 | CALL mathop (topp,req_sz,buff_tmp2,missing_val, & |
---|
1048 | & nbindex,ijndex,scal,siz1,buff_tmp1) |
---|
1049 | var(:) = buff_tmp1(1:siz1) |
---|
1050 | ELSE |
---|
1051 | CALL ipslerr (3,'resget_opp_r1d', & |
---|
1052 | 'The operation you wish to do on the variable for the ',& |
---|
1053 | 'restart file is not allowed.',topp) |
---|
1054 | ENDIF |
---|
1055 | !----------------------------- |
---|
1056 | END SUBROUTINE restget_opp_r1d |
---|
1057 | !=== |
---|
1058 | SUBROUTINE restget_opp_r2d & |
---|
1059 | & (fid,vname_q,iim,jjm,llm,itau,def_beha, & |
---|
1060 | & var,MY_OPERATOR,nbindex,ijndex) |
---|
1061 | !--------------------------------------------------------------------- |
---|
1062 | !- This subroutine serves as an interface to restget_real |
---|
1063 | !- |
---|
1064 | !- Should work as restput_opp_r2d but the other way around ! |
---|
1065 | !--------------------------------------------------------------------- |
---|
1066 | IMPLICIT NONE |
---|
1067 | !- |
---|
1068 | INTEGER :: fid |
---|
1069 | CHARACTER(LEN=*) :: vname_q |
---|
1070 | INTEGER :: iim,jjm,llm,itau |
---|
1071 | LOGICAL def_beha |
---|
1072 | REAL :: var(:,:) |
---|
1073 | CHARACTER(LEN=*) :: MY_OPERATOR |
---|
1074 | INTEGER :: nbindex,ijndex(nbindex) |
---|
1075 | !- |
---|
1076 | ! LOCAL |
---|
1077 | !- |
---|
1078 | INTEGER :: jj,req_sz,ist,var_sz,siz1 |
---|
1079 | REAL :: scal |
---|
1080 | CHARACTER(LEN=7) :: topp |
---|
1081 | LOGICAL :: check = .FALSE. |
---|
1082 | !--------------------------------------------------------------------- |
---|
1083 | !- |
---|
1084 | ! 0.0 What size should be the data in the file |
---|
1085 | !- |
---|
1086 | req_sz = 1 |
---|
1087 | IF (nbindex == iim .AND. llm <= 1) THEN |
---|
1088 | IF (xax_infs(fid,1,1) > 0) req_sz = req_sz*xax_infs(fid,1,1) |
---|
1089 | IF (yax_infs(fid,1,1) > 0) req_sz = req_sz*yax_infs(fid,1,1) |
---|
1090 | ELSE |
---|
1091 | CALL ipslerr (3,'resget_opp_r2d', & |
---|
1092 | 'Unable to performe an operation on this variable as it has', & |
---|
1093 | 'a second and third dimension',vname_q) |
---|
1094 | ENDIF |
---|
1095 | !- |
---|
1096 | IF (jjm < 1) THEN |
---|
1097 | CALL ipslerr (3,'resget_opp_r2d', & |
---|
1098 | 'Please specify a second dimension which is the', & |
---|
1099 | 'layer on which the operations are performed',vname_q) |
---|
1100 | ENDIF |
---|
1101 | !- |
---|
1102 | ! 1.0 Allocate the temporary buffer we need |
---|
1103 | ! to put the variable in right dimension |
---|
1104 | !- |
---|
1105 | siz1 = SIZE(var,1) |
---|
1106 | CALL rest_alloc (1,siz1,check,'restget_opp_r2d') |
---|
1107 | CALL rest_alloc (2,req_sz*jjm,check,'restget_opp_r2d') |
---|
1108 | !- |
---|
1109 | ! 2.0 Here we get the full variable from the restart file |
---|
1110 | !- |
---|
1111 | CALL restget_real & |
---|
1112 | & (fid,vname_q,xax_infs(fid,1,1),yax_infs(fid,1,1), & |
---|
1113 | & jjm,itau,def_beha,buff_tmp2) |
---|
1114 | !- |
---|
1115 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1116 | ! into the variable the model expects |
---|
1117 | !- |
---|
1118 | topp = MY_OPERATOR(1:MIN(LEN_TRIM(MY_OPERATOR),7)) |
---|
1119 | !- |
---|
1120 | IF (INDEX(indchfun,topp(:LEN_TRIM(topp))) > 0) THEN |
---|
1121 | scal = missing_val |
---|
1122 | var_sz = siz1 |
---|
1123 | DO jj = 1,jjm |
---|
1124 | ist = (jj-1)*req_sz+1 |
---|
1125 | CALL mathop (topp,req_sz,buff_tmp2(ist:ist+req_sz-1), & |
---|
1126 | & missing_val,nbindex,ijndex,scal,var_sz,buff_tmp1) |
---|
1127 | var(:,jj) = buff_tmp1(1:siz1) |
---|
1128 | ENDDO |
---|
1129 | ELSE |
---|
1130 | CALL ipslerr (3,'resget_opp_r2d', & |
---|
1131 | 'The operation you wish to do on the variable for the ',& |
---|
1132 | 'restart file is not allowed.',topp) |
---|
1133 | ENDIF |
---|
1134 | !----------------------------- |
---|
1135 | END SUBROUTINE restget_opp_r2d |
---|
1136 | !=== |
---|
1137 | SUBROUTINE restget_r1d & |
---|
1138 | & (fid,vname_q,iim,jjm,llm,itau,def_beha,var) |
---|
1139 | !--------------------------------------------------------------------- |
---|
1140 | !- This subroutine serves as an interface to restget_real |
---|
1141 | !--------------------------------------------------------------------- |
---|
1142 | IMPLICIT NONE |
---|
1143 | !- |
---|
1144 | INTEGER :: fid |
---|
1145 | CHARACTER(LEN=*) :: vname_q |
---|
1146 | INTEGER :: iim,jjm,llm,itau |
---|
1147 | LOGICAL :: def_beha |
---|
1148 | REAL :: var(:) |
---|
1149 | !- |
---|
1150 | ! LOCAL |
---|
1151 | !- |
---|
1152 | INTEGER :: ji,jl,req_sz,var_sz,siz1 |
---|
1153 | CHARACTER(LEN=70) :: str,str2 |
---|
1154 | LOGICAL :: check = .FALSE. |
---|
1155 | !--------------------------------------------------------------------- |
---|
1156 | !- |
---|
1157 | ! 1.0 Allocate the temporary buffer we need |
---|
1158 | ! to put the variable in right dimension |
---|
1159 | !- |
---|
1160 | siz1 = SIZE(var) |
---|
1161 | var_sz = siz1 |
---|
1162 | CALL rest_alloc (1,var_sz,check,'restget_r1d') |
---|
1163 | !- |
---|
1164 | ! 2.0 Here we could check if the sizes specified agree |
---|
1165 | ! with the size of the variable provided |
---|
1166 | !- |
---|
1167 | req_sz = 1 |
---|
1168 | IF (iim > 0) req_sz = req_sz*iim |
---|
1169 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1170 | IF (llm > 0) req_sz = req_sz*llm |
---|
1171 | IF (req_sz > var_sz) THEN |
---|
1172 | WRITE(str, & |
---|
1173 | & '("Size of variable requested from file should be ",I6)') req_sz |
---|
1174 | WRITE(str2, & |
---|
1175 | & '("but the provided variable can only hold ",I6)') var_sz |
---|
1176 | CALL ipslerr (3,'restget_r1d',str,str2,' ') |
---|
1177 | ENDIF |
---|
1178 | IF (req_sz < var_sz) THEN |
---|
1179 | WRITE(str, & |
---|
1180 | & '("the size of variable requested from file is ",I6)') req_sz |
---|
1181 | WRITE(str2, & |
---|
1182 | & '("but the provided variable can hold ",I6)') var_sz |
---|
1183 | CALL ipslerr (2,'restget_r1d', & |
---|
1184 | 'There could be a problem here :',str,str2) |
---|
1185 | ENDIF |
---|
1186 | !- |
---|
1187 | CALL restget_real & |
---|
1188 | & (fid,vname_q,iim,jjm,llm,itau,def_beha,buff_tmp1) |
---|
1189 | !- |
---|
1190 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1191 | ! into the variable the model expects |
---|
1192 | !- |
---|
1193 | jl=0 |
---|
1194 | DO ji=1,siz1 |
---|
1195 | jl=jl+1 |
---|
1196 | var(ji) = buff_tmp1(jl) |
---|
1197 | ENDDO |
---|
1198 | !------------------------- |
---|
1199 | END SUBROUTINE restget_r1d |
---|
1200 | !=== |
---|
1201 | SUBROUTINE restget_r2d & |
---|
1202 | & (fid,vname_q,iim,jjm,llm,itau,def_beha,var) |
---|
1203 | !--------------------------------------------------------------------- |
---|
1204 | !- This subroutine serves as an interface to restget_real |
---|
1205 | !--------------------------------------------------------------------- |
---|
1206 | IMPLICIT NONE |
---|
1207 | !- |
---|
1208 | INTEGER :: fid |
---|
1209 | CHARACTER(LEN=*) :: vname_q |
---|
1210 | INTEGER :: iim,jjm,llm,itau |
---|
1211 | LOGICAL :: def_beha |
---|
1212 | REAL :: var(:,:) |
---|
1213 | !- |
---|
1214 | ! LOCAL |
---|
1215 | !- |
---|
1216 | INTEGER :: ji,jj,jl,req_sz,var_sz,siz1,siz2 |
---|
1217 | CHARACTER(LEN=70) :: str,str2 |
---|
1218 | LOGICAL :: check = .FALSE. |
---|
1219 | !--------------------------------------------------------------------- |
---|
1220 | !- |
---|
1221 | ! 1.0 Allocate the temporary buffer we need |
---|
1222 | ! to put the variable in right dimension |
---|
1223 | !- |
---|
1224 | siz1 = SIZE(var,1) |
---|
1225 | siz2 = SIZE(var,2) |
---|
1226 | var_sz = siz1*siz2 |
---|
1227 | CALL rest_alloc (1,var_sz,check,'restget_r2d') |
---|
1228 | !- |
---|
1229 | ! 2.0 Here we check if the sizes specified agree |
---|
1230 | ! with the size of the variable provided |
---|
1231 | !- |
---|
1232 | req_sz = 1 |
---|
1233 | IF (iim > 0) req_sz = req_sz*iim |
---|
1234 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1235 | IF (llm > 0) req_sz = req_sz*llm |
---|
1236 | IF (req_sz > var_sz) THEN |
---|
1237 | WRITE(str, & |
---|
1238 | & '("Size of variable ",A, & |
---|
1239 | & //" requested from file should be ",I6)') TRIM(vname_q),req_sz |
---|
1240 | WRITE(str2, & |
---|
1241 | & '("but the provided variable can only hold ",I6)') var_sz |
---|
1242 | CALL ipslerr (3,'restget_r2d',str,str2,' ') |
---|
1243 | ENDIF |
---|
1244 | IF (req_sz < var_sz) THEN |
---|
1245 | WRITE(str, & |
---|
1246 | & '("Size of variable ",A, & |
---|
1247 | & //" requested from file is ",I6)') TRIM(vname_q),req_sz |
---|
1248 | WRITE(str2,'("but the provided variable can hold ",I6)') var_sz |
---|
1249 | CALL ipslerr (2,'restget_r2d', & |
---|
1250 | 'There could be a problem here :',str,str2) |
---|
1251 | ENDIF |
---|
1252 | !- |
---|
1253 | CALL restget_real & |
---|
1254 | & (fid,vname_q,iim,jjm,llm,itau,def_beha,buff_tmp1) |
---|
1255 | !- |
---|
1256 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1257 | ! into the variable the model expects |
---|
1258 | !- |
---|
1259 | jl=0 |
---|
1260 | DO jj=1,siz2 |
---|
1261 | DO ji=1,siz1 |
---|
1262 | jl=jl+1 |
---|
1263 | var(ji,jj) = buff_tmp1(jl) |
---|
1264 | ENDDO |
---|
1265 | ENDDO |
---|
1266 | !------------------------- |
---|
1267 | END SUBROUTINE restget_r2d |
---|
1268 | !=== |
---|
1269 | SUBROUTINE restget_r3d & |
---|
1270 | (fid,vname_q,iim,jjm,llm,itau,def_beha,var) |
---|
1271 | !--------------------------------------------------------------------- |
---|
1272 | !- This subroutine serves as an interface to restget_real |
---|
1273 | !--------------------------------------------------------------------- |
---|
1274 | IMPLICIT NONE |
---|
1275 | !- |
---|
1276 | INTEGER :: fid |
---|
1277 | CHARACTER(LEN=*) :: vname_q |
---|
1278 | INTEGER :: iim,jjm,llm,itau |
---|
1279 | LOGICAL def_beha |
---|
1280 | REAL :: var(:,:,:) |
---|
1281 | !- |
---|
1282 | ! LOCAL |
---|
1283 | !- |
---|
1284 | INTEGER :: ji,jj,jk,jl,req_sz,var_sz,siz1,siz2,siz3 |
---|
1285 | CHARACTER(LEN=70) :: str,str2 |
---|
1286 | LOGICAL :: check = .FALSE. |
---|
1287 | !--------------------------------------------------------------------- |
---|
1288 | !- |
---|
1289 | ! 1.0 Allocate the temporary buffer we need |
---|
1290 | ! to put the variable in right dimension |
---|
1291 | !- |
---|
1292 | siz1 = SIZE(var,1) |
---|
1293 | siz2 = SIZE(var,2) |
---|
1294 | siz3 = SIZE(var,3) |
---|
1295 | var_sz = siz1*siz2*siz3 |
---|
1296 | CALL rest_alloc (1,var_sz,check,'restget_r3d') |
---|
1297 | !- |
---|
1298 | ! 2.0 Here we check if the sizes specified agree |
---|
1299 | ! with the size of the variable provided |
---|
1300 | !- |
---|
1301 | req_sz = 1 |
---|
1302 | IF (iim > 0) req_sz = req_sz*iim |
---|
1303 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1304 | IF (llm > 0) req_sz = req_sz*llm |
---|
1305 | IF (req_sz > var_sz) THEN |
---|
1306 | WRITE(str, & |
---|
1307 | & '("Size of variable ",A, & |
---|
1308 | & //" requested from file should be ",I6)') TRIM(vname_q),req_sz |
---|
1309 | WRITE(str2, & |
---|
1310 | & '("but the provided variable can only hold ",I6)') var_sz |
---|
1311 | CALL ipslerr (3,'restget_r3d',str,str2,' ') |
---|
1312 | ENDIF |
---|
1313 | IF (req_sz < var_sz) THEN |
---|
1314 | WRITE(str, & |
---|
1315 | & '("Size of variable ",A, & |
---|
1316 | & //" requested from file is ",I6)') TRIM(vname_q),req_sz |
---|
1317 | WRITE(str2,'("but the provided variable can hold ",I6)') var_sz |
---|
1318 | CALL ipslerr (2,'restget_r3d', & |
---|
1319 | 'There could be a problem here :',str,str2) |
---|
1320 | ENDIF |
---|
1321 | !- |
---|
1322 | CALL restget_real & |
---|
1323 | (fid,vname_q,iim,jjm,llm,itau,def_beha,buff_tmp1) |
---|
1324 | !- |
---|
1325 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1326 | ! into the variable the model expects |
---|
1327 | !- |
---|
1328 | jl=0 |
---|
1329 | DO jk=1,siz3 |
---|
1330 | DO jj=1,siz2 |
---|
1331 | DO ji=1,siz1 |
---|
1332 | jl=jl+1 |
---|
1333 | var(ji,jj,jk) = buff_tmp1(jl) |
---|
1334 | ENDDO |
---|
1335 | ENDDO |
---|
1336 | ENDDO |
---|
1337 | !------------------------- |
---|
1338 | END SUBROUTINE restget_r3d |
---|
1339 | !=== |
---|
1340 | SUBROUTINE restget_real & |
---|
1341 | (fid,vname_q,iim,jjm,llm,itau,def_beha,var) |
---|
1342 | !--------------------------------------------------------------------- |
---|
1343 | !- This subroutine is for getting a variable from the restart file. |
---|
1344 | !- A number of verifications will be made : |
---|
1345 | !- - Is this the first time we read this variable ? |
---|
1346 | !- - Are the dimensions correct ? |
---|
1347 | !- - Is the correct time step present in the file |
---|
1348 | !- - is a default behaviour possible. If not the model is stoped. |
---|
1349 | !- Default procedure is to write the content of val_exp on all values. |
---|
1350 | !- |
---|
1351 | !- INPUT |
---|
1352 | !- |
---|
1353 | !- fid : Identification of the file |
---|
1354 | !- vname_q : Name of the variable to be read |
---|
1355 | !- iim, jjm ,llm : Dimensions of the variable that should be read |
---|
1356 | !- itau : Time step at whcih we are when we want |
---|
1357 | !- to read the variable |
---|
1358 | !- def_beha : If the model can restart without this variable |
---|
1359 | !- then some strange value is given. |
---|
1360 | !- |
---|
1361 | !- OUTPUT |
---|
1362 | !- |
---|
1363 | !- var : Variable in which the data is put |
---|
1364 | !--------------------------------------------------------------------- |
---|
1365 | IMPLICIT NONE |
---|
1366 | !- |
---|
1367 | INTEGER :: fid |
---|
1368 | CHARACTER(LEN=*) :: vname_q |
---|
1369 | INTEGER :: iim,jjm,llm,itau |
---|
1370 | LOGICAL :: def_beha |
---|
1371 | REAL :: var(:) |
---|
1372 | !- |
---|
1373 | ! LOCAL |
---|
1374 | !- |
---|
1375 | INTEGER :: vid,vnb,ncfid |
---|
1376 | INTEGER :: iret,index,it,ndim,ia |
---|
1377 | CHARACTER(LEN=70) str,str2 |
---|
1378 | CHARACTER(LEN=80) attname |
---|
1379 | INTEGER,DIMENSION(4) :: corner,edge |
---|
1380 | !--------------------------------------------------------------------- |
---|
1381 | !- |
---|
1382 | ncfid = netcdf_id(fid,1) |
---|
1383 | !- |
---|
1384 | CALL find_str (varname_in(fid,1:nbvar_in(fid)),vname_q,vnb) |
---|
1385 | !- |
---|
1386 | ! 1.0 If the variable is not present then ERROR or filled up |
---|
1387 | ! by default values if allowed |
---|
1388 | !- |
---|
1389 | IF (vnb < 0) THEN |
---|
1390 | IF (def_beha) THEN |
---|
1391 | !----- |
---|
1392 | lock_valexp = .TRUE. |
---|
1393 | var(:) = val_exp |
---|
1394 | !---- |
---|
1395 | str = 'Variable '//TRIM(vname_q) & |
---|
1396 | //' is not present in the restart file' |
---|
1397 | CALL ipslerr (1,'restget', & |
---|
1398 | & str,'but default values are used to fill in',' ') |
---|
1399 | !---- |
---|
1400 | IF (nbvar_in(fid) >= max_var) THEN |
---|
1401 | CALL ipslerr (3,'restget', & |
---|
1402 | 'Too many variables for the restcom module', & |
---|
1403 | 'Please increase the value of max_var',' ') |
---|
1404 | ENDIF |
---|
1405 | nbvar_in(fid) = nbvar_in(fid)+1 |
---|
1406 | vnb = nbvar_in(fid) |
---|
1407 | varname_in(fid,vnb) = vname_q |
---|
1408 | touched_in(fid,vnb) = .TRUE. |
---|
1409 | !----- |
---|
1410 | CALL restdefv (fid,vname_q,iim,jjm,llm,.TRUE.) |
---|
1411 | !----- |
---|
1412 | ELSE |
---|
1413 | str = 'Variable '//TRIM(vname_q) & |
---|
1414 | //' is not present in the restart file' |
---|
1415 | CALL ipslerr (3,'restget', & |
---|
1416 | & str,'but it is need to restart the model',' ') |
---|
1417 | ENDIF |
---|
1418 | !--- |
---|
1419 | ELSE |
---|
1420 | !--- |
---|
1421 | !-- 2.0 Check if the variable has not yet been read |
---|
1422 | !-- and that the time is OK |
---|
1423 | !--- |
---|
1424 | vid = varid_in(fid,vnb) |
---|
1425 | !--- |
---|
1426 | nbvar_read(fid) = nbvar_read(fid)+1 |
---|
1427 | !--- |
---|
1428 | IF (touched_in(fid,vnb)) THEN |
---|
1429 | str = 'Variable '//TRIM(vname_q) & |
---|
1430 | //' has already been read from file' |
---|
1431 | CALL ipslerr (3,'restget',str,' ',' ') |
---|
1432 | ENDIF |
---|
1433 | !--- |
---|
1434 | !-- 3.0 get the time step of the restart file |
---|
1435 | !-- and check if it is correct |
---|
1436 | !--- |
---|
1437 | index = -1 |
---|
1438 | DO it=1,tax_size_in(fid) |
---|
1439 | IF (t_index(fid,it) == itau) index = it |
---|
1440 | ENDDO |
---|
1441 | !--- |
---|
1442 | IF (index < 0) THEN |
---|
1443 | str = 'The time step requested for variable '//TRIM(vname_q) |
---|
1444 | CALL ipslerr (3,'restget', & |
---|
1445 | & str,'is not available in the current file',' ') |
---|
1446 | ENDIF |
---|
1447 | !--- |
---|
1448 | !-- 4.0 Read the data. Note that the variables in the restart files |
---|
1449 | !-- have no time axis is and thus we write -1 |
---|
1450 | !--- |
---|
1451 | str='Incorrect dimension for '//TRIM(vname_q) |
---|
1452 | ndim = 0 |
---|
1453 | IF (iim > 0) THEN |
---|
1454 | ndim = ndim+1 |
---|
1455 | IF (vardims_in(fid,vnb,ndim) == iim) THEN |
---|
1456 | corner(ndim) = 1 |
---|
1457 | edge(ndim) = iim |
---|
1458 | ELSE |
---|
1459 | WRITE (str2,'("Incompatibility for iim : ",I6,I6)') & |
---|
1460 | iim,vardims_in(fid,vnb,ndim) |
---|
1461 | CALL ipslerr (3,'restget',str,str2,' ') |
---|
1462 | ENDIF |
---|
1463 | ENDIF |
---|
1464 | !--- |
---|
1465 | IF (jjm > 0) THEN |
---|
1466 | ndim = ndim+1 |
---|
1467 | IF (vardims_in(fid,vnb,ndim) == jjm) THEN |
---|
1468 | corner(ndim) = 1 |
---|
1469 | edge(ndim) = jjm |
---|
1470 | ELSE |
---|
1471 | WRITE (str2,'("Incompatibility for jjm : ",I6,I6)') & |
---|
1472 | jjm,vardims_in(fid,vnb,ndim) |
---|
1473 | CALL ipslerr (3,'restget',str,str2,' ') |
---|
1474 | ENDIF |
---|
1475 | ENDIF |
---|
1476 | !--- |
---|
1477 | IF (llm > 0) THEN |
---|
1478 | ndim = ndim+1 |
---|
1479 | IF (vardims_in(fid,vnb,ndim) == llm) THEN |
---|
1480 | corner(ndim) = 1 |
---|
1481 | edge(ndim) = llm |
---|
1482 | ELSE |
---|
1483 | WRITE (str2,'("Incompatibility for llm : ",I6,I6)') & |
---|
1484 | llm,vardims_in(fid,vnb,ndim) |
---|
1485 | CALL ipslerr (3,'restget',str,str2,' ') |
---|
1486 | ENDIF |
---|
1487 | ENDIF |
---|
1488 | !--- |
---|
1489 | !-- Time |
---|
1490 | !--- |
---|
1491 | ndim = ndim+1 |
---|
1492 | corner(ndim) = index |
---|
1493 | !!????? edge(ndim) = index |
---|
1494 | edge(ndim) = 1 |
---|
1495 | !--- |
---|
1496 | iret = NF90_GET_VAR(ncfid,vid,var, & |
---|
1497 | & start=corner(1:ndim),count=edge(1:ndim)) |
---|
1498 | !--- |
---|
1499 | !-- 5.0 The variable we have just read is created |
---|
1500 | !-- in the next restart file |
---|
1501 | !--- |
---|
1502 | IF ( (netcdf_id(fid,1) /= netcdf_id(fid,2)) & |
---|
1503 | & .AND.(netcdf_id(fid,2) > 0) ) THEN |
---|
1504 | !----- |
---|
1505 | CALL restdefv (fid,vname_q,iim,jjm,llm,.FALSE.) |
---|
1506 | !----- |
---|
1507 | DO ia = 1,varatt_in(fid,vnb) |
---|
1508 | iret = NF90_INQ_ATTNAME(ncfid,vid,ia,attname) |
---|
1509 | iret = NF90_COPY_ATT(ncfid,vid,attname, & |
---|
1510 | & netcdf_id(fid,2),varid_out(fid,nbvar_out(fid))) |
---|
1511 | ENDDO |
---|
1512 | !----- |
---|
1513 | IF (itau_out(fid) >= 0) THEN |
---|
1514 | iret = NF90_ENDDEF(netcdf_id(fid,2)) |
---|
1515 | ENDIF |
---|
1516 | ENDIF |
---|
1517 | !--- |
---|
1518 | ENDIF |
---|
1519 | !-------------------------- |
---|
1520 | END SUBROUTINE restget_real |
---|
1521 | !=== |
---|
1522 | SUBROUTINE restput_opp_r1d & |
---|
1523 | & (fid,vname_q,iim,jjm,llm,itau,var,MY_OPERATOR,nbindex,ijndex) |
---|
1524 | !--------------------------------------------------------------------- |
---|
1525 | !- This subroutine is the interface to restput_real which allows |
---|
1526 | !- to re-index data onto the original grid of the restart file. |
---|
1527 | !- The logic we use is still fuzzy in my mind but that is probably |
---|
1528 | !- only because I have not yet though through everything. |
---|
1529 | !- |
---|
1530 | !- In the case iim = nbindex it means that the user attempts |
---|
1531 | !- to project a vector back onto the original 2D or 3D field. |
---|
1532 | !- This requires that jjm and llm be equal to 1 or 0, |
---|
1533 | !- else I would not know what it means. |
---|
1534 | !--------------------------------------------------------------------- |
---|
1535 | IMPLICIT NONE |
---|
1536 | !- |
---|
1537 | INTEGER :: fid |
---|
1538 | CHARACTER(LEN=*) :: vname_q |
---|
1539 | INTEGER :: iim,jjm,llm,itau |
---|
1540 | REAL :: var(:) |
---|
1541 | CHARACTER(LEN=*) :: MY_OPERATOR |
---|
1542 | INTEGER :: nbindex,ijndex(nbindex) |
---|
1543 | !- |
---|
1544 | ! LOCAL |
---|
1545 | !- |
---|
1546 | INTEGER :: req_sz,siz1 |
---|
1547 | REAL :: scal |
---|
1548 | CHARACTER(LEN=7) :: topp |
---|
1549 | LOGICAL :: check = .FALSE. |
---|
1550 | !--------------------------------------------------------------------- |
---|
1551 | !- |
---|
1552 | ! 0.0 What size should be the data in the file |
---|
1553 | !- |
---|
1554 | req_sz = 1 |
---|
1555 | IF ( nbindex == iim .AND. jjm <= 1 .AND. llm <= 1) THEN |
---|
1556 | IF (xax_infs(fid,1,1) > 0) req_sz = req_sz*xax_infs(fid,1,1) |
---|
1557 | IF (yax_infs(fid,1,1) > 0) req_sz = req_sz*yax_infs(fid,1,1) |
---|
1558 | IF (zax_infs(fid,1,1) > 0) req_sz = req_sz*zax_infs(fid,1,1) |
---|
1559 | ELSE |
---|
1560 | CALL ipslerr (3,'restput_opp_r1d', & |
---|
1561 | 'Unable to performe an operation on this variable as it has', & |
---|
1562 | 'a second and third dimension',vname_q) |
---|
1563 | ENDIF |
---|
1564 | !- |
---|
1565 | ! 1.0 Allocate the temporary buffer we need |
---|
1566 | ! to put the variable in right dimension |
---|
1567 | !- |
---|
1568 | siz1 = SIZE(var) |
---|
1569 | CALL rest_alloc (1,siz1,check,'restput_opp_r1d') |
---|
1570 | CALL rest_alloc (2,req_sz,check,'restput_opp_r1d') |
---|
1571 | !- |
---|
1572 | ! 2.0 We do the operation needed. |
---|
1573 | ! It can only be a re-indexing operation. |
---|
1574 | ! You would not want to change the values in a restart file or ? |
---|
1575 | !- |
---|
1576 | topp = MY_OPERATOR(1:MIN(LEN_TRIM(MY_OPERATOR),7)) |
---|
1577 | !- |
---|
1578 | IF (INDEX(indchfun,topp(:LEN_TRIM(topp))) > 0) THEN |
---|
1579 | scal = missing_val |
---|
1580 | buff_tmp1(1:siz1) = var(:) |
---|
1581 | CALL mathop & |
---|
1582 | & (topp,siz1,buff_tmp1,missing_val,nbindex,ijndex, & |
---|
1583 | & scal,req_sz,buff_tmp2) |
---|
1584 | ELSE |
---|
1585 | CALL ipslerr (3,'restput_opp_r1d', & |
---|
1586 | & 'The operation you wish to do on the variable for the ', & |
---|
1587 | & 'restart file is not allowed.',topp) |
---|
1588 | ENDIF |
---|
1589 | !- |
---|
1590 | CALL restput_real & |
---|
1591 | & (fid,vname_q,xax_infs(fid,1,1),yax_infs(fid,1,1), & |
---|
1592 | & zax_infs(fid,1,1),itau,buff_tmp2) |
---|
1593 | !----------------------------- |
---|
1594 | END SUBROUTINE restput_opp_r1d |
---|
1595 | !=== |
---|
1596 | SUBROUTINE restput_opp_r2d & |
---|
1597 | & (fid,vname_q,iim,jjm,llm,itau,var,MY_OPERATOR,nbindex,ijndex) |
---|
1598 | !--------------------------------------------------------------------- |
---|
1599 | !- This subroutine is the interface to restput_real which allows |
---|
1600 | !- to re-index data onto the original grid of the restart file. |
---|
1601 | !- The logic we use is still fuzzy in my mind but that is probably |
---|
1602 | !- only because I have not yet though through everything. |
---|
1603 | !- |
---|
1604 | !- In the case iim = nbindex it means that the user attempts |
---|
1605 | !- to project the first dimension of the matrix back onto a 3D field |
---|
1606 | !- where jjm will be the third dimension. |
---|
1607 | !- Here we do not allow for 4D data, thus we will take the first |
---|
1608 | !- two dimensions in the file and require that llm = 1. |
---|
1609 | !- These are pretty heavy constraints but I do not know how |
---|
1610 | !- to make it more general. I need to think about it some more. |
---|
1611 | !--------------------------------------------------------------------- |
---|
1612 | IMPLICIT NONE |
---|
1613 | !- |
---|
1614 | INTEGER :: fid |
---|
1615 | CHARACTER(LEN=*) :: vname_q |
---|
1616 | INTEGER :: iim,jjm,llm,itau |
---|
1617 | REAL :: var(:,:) |
---|
1618 | CHARACTER(LEN=*) :: MY_OPERATOR |
---|
1619 | INTEGER :: nbindex,ijndex(nbindex) |
---|
1620 | !- |
---|
1621 | ! LOCAL |
---|
1622 | !- |
---|
1623 | INTEGER :: jj,req_sz,ist,siz1 |
---|
1624 | REAL :: scal |
---|
1625 | CHARACTER(LEN=7) :: topp |
---|
1626 | LOGICAL :: check = .FALSE. |
---|
1627 | !--------------------------------------------------------------------- |
---|
1628 | !- |
---|
1629 | ! 0.0 What size should be the data in the file |
---|
1630 | !- |
---|
1631 | req_sz = 1 |
---|
1632 | IF ( nbindex == iim .AND. llm <= 1) THEN |
---|
1633 | IF (xax_infs(fid,1,1) > 0) req_sz = req_sz*xax_infs(fid,1,1) |
---|
1634 | IF (yax_infs(fid,1,1) > 0) req_sz = req_sz*yax_infs(fid,1,1) |
---|
1635 | ELSE |
---|
1636 | CALL ipslerr (3,'restput_opp_r2d', & |
---|
1637 | 'Unable to performe an operation on this variable as it has', & |
---|
1638 | 'a second and third dimension',vname_q) |
---|
1639 | ENDIF |
---|
1640 | !- |
---|
1641 | IF (jjm < 1) THEN |
---|
1642 | CALL ipslerr (3,'restput_opp_r2d', & |
---|
1643 | 'Please specify a second dimension which is the', & |
---|
1644 | 'layer on which the operations are performed',vname_q) |
---|
1645 | ENDIF |
---|
1646 | !- |
---|
1647 | ! 1.0 Allocate the temporary buffer we need |
---|
1648 | ! to put the variable in right dimension |
---|
1649 | !- |
---|
1650 | siz1 = SIZE(var,1) |
---|
1651 | CALL rest_alloc (1,siz1,check,'restput_opp_r2d') |
---|
1652 | CALL rest_alloc (2,req_sz*jjm,check,'restput_opp_r2d') |
---|
1653 | !- |
---|
1654 | ! 2.0 We do the operation needed. |
---|
1655 | ! It can only be a re-indexing operation. |
---|
1656 | ! You would not want to change the values in a restart file or ? |
---|
1657 | !- |
---|
1658 | topp = MY_OPERATOR(1:MIN(LEN_TRIM(MY_OPERATOR),7)) |
---|
1659 | !- |
---|
1660 | IF (INDEX(indchfun,topp(:LEN_TRIM(topp))) > 0) THEN |
---|
1661 | scal = missing_val |
---|
1662 | DO jj = 1,jjm |
---|
1663 | buff_tmp1(1:siz1) = var(:,jj) |
---|
1664 | ist = (jj-1)*req_sz+1 |
---|
1665 | CALL mathop & |
---|
1666 | & (topp,siz1,buff_tmp1,missing_val,nbindex,ijndex, & |
---|
1667 | & scal,req_sz,buff_tmp2(ist:ist+req_sz-1)) |
---|
1668 | ENDDO |
---|
1669 | ELSE |
---|
1670 | CALL ipslerr (3,'restput_opp_r2d', & |
---|
1671 | & 'The operation you wish to do on the variable for the ', & |
---|
1672 | & 'restart file is not allowed.',topp) |
---|
1673 | ENDIF |
---|
1674 | !- |
---|
1675 | CALL restput_real & |
---|
1676 | & (fid,vname_q,xax_infs(fid,1,1),yax_infs(fid,1,1), & |
---|
1677 | & jjm,itau,buff_tmp2) |
---|
1678 | !----------------------------- |
---|
1679 | END SUBROUTINE restput_opp_r2d |
---|
1680 | !=== |
---|
1681 | SUBROUTINE restput_r1d (fid,vname_q,iim,jjm,llm,itau,var) |
---|
1682 | !--------------------------------------------------------------------- |
---|
1683 | !- This subroutine serves as an interface to restput_real |
---|
1684 | !--------------------------------------------------------------------- |
---|
1685 | IMPLICIT NONE |
---|
1686 | !- |
---|
1687 | INTEGER :: fid |
---|
1688 | CHARACTER(LEN=*) :: vname_q |
---|
1689 | INTEGER :: iim,jjm,llm,itau |
---|
1690 | REAL :: var(:) |
---|
1691 | !- |
---|
1692 | ! LOCAL |
---|
1693 | !- |
---|
1694 | INTEGER :: ji,jl,req_sz,var_sz,siz1 |
---|
1695 | CHARACTER(LEN=70) :: str,str2 |
---|
1696 | LOGICAL :: check = .FALSE. |
---|
1697 | !--------------------------------------------------------------------- |
---|
1698 | !- |
---|
1699 | ! 1.0 Allocate the temporary buffer we need |
---|
1700 | ! to put the variable in right dimension |
---|
1701 | !- |
---|
1702 | siz1 = SIZE(var) |
---|
1703 | var_sz = siz1 |
---|
1704 | CALL rest_alloc (1,var_sz,check,'restput_r1d') |
---|
1705 | !- |
---|
1706 | ! 2.0 Here we could check if the sizes specified agree |
---|
1707 | ! with the size of the variable provided |
---|
1708 | !- |
---|
1709 | req_sz = 1 |
---|
1710 | IF (iim > 0) req_sz = req_sz*iim |
---|
1711 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1712 | IF (llm > 0) req_sz = req_sz*llm |
---|
1713 | IF (req_sz > var_sz) THEN |
---|
1714 | WRITE(str, & |
---|
1715 | & '("Size of variable put to the file should be ",I6)') req_sz |
---|
1716 | WRITE(str2, & |
---|
1717 | & '("but the provided variable is of size ",I6)') var_sz |
---|
1718 | CALL ipslerr (3,'restput_r1d',str,str2,' ') |
---|
1719 | ENDIF |
---|
1720 | IF (req_sz < var_sz) THEN |
---|
1721 | WRITE(str,'("the size of variable put to the file is ",I6)') req_sz |
---|
1722 | WRITE(str2,'("but the provided variable is larger ",I6)') var_sz |
---|
1723 | CALL ipslerr (2,'restput_r1d', & |
---|
1724 | 'There could be a problem here :',str,str2) |
---|
1725 | ENDIF |
---|
1726 | !- |
---|
1727 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1728 | ! into the variable the model expects |
---|
1729 | !- |
---|
1730 | jl=0 |
---|
1731 | DO ji=1,siz1 |
---|
1732 | jl=jl+1 |
---|
1733 | buff_tmp1(jl) = var(ji) |
---|
1734 | ENDDO |
---|
1735 | !- |
---|
1736 | CALL restput_real (fid,vname_q,iim,jjm,llm,itau,buff_tmp1) |
---|
1737 | !------------------------- |
---|
1738 | END SUBROUTINE restput_r1d |
---|
1739 | !=== |
---|
1740 | SUBROUTINE restput_r2d (fid,vname_q,iim,jjm,llm,itau,var) |
---|
1741 | !--------------------------------------------------------------------- |
---|
1742 | !- This subroutine serves as an interface to restput_real |
---|
1743 | !--------------------------------------------------------------------- |
---|
1744 | IMPLICIT NONE |
---|
1745 | !- |
---|
1746 | INTEGER :: fid |
---|
1747 | CHARACTER(LEN=*) :: vname_q |
---|
1748 | INTEGER :: iim,jjm,llm,itau |
---|
1749 | REAL :: var(:,:) |
---|
1750 | !- |
---|
1751 | ! LOCAL |
---|
1752 | !- |
---|
1753 | INTEGER :: ji,jj,jl,req_sz,var_sz,siz1,siz2 |
---|
1754 | CHARACTER(LEN=70) :: str,str2 |
---|
1755 | LOGICAL :: check = .FALSE. |
---|
1756 | !--------------------------------------------------------------------- |
---|
1757 | !- |
---|
1758 | ! 1.0 Allocate the temporary buffer we need |
---|
1759 | ! to put the variable in right dimension |
---|
1760 | !- |
---|
1761 | siz1 = SIZE(var,1) |
---|
1762 | siz2 = SIZE(var,2) |
---|
1763 | var_sz = siz1*siz2 |
---|
1764 | CALL rest_alloc (1,var_sz,check,'restput_r2d') |
---|
1765 | !- |
---|
1766 | ! 2.0 Here we could check if the sizes specified agree |
---|
1767 | ! with the size of the variable provided |
---|
1768 | !- |
---|
1769 | req_sz = 1 |
---|
1770 | IF (iim > 0) req_sz = req_sz*iim |
---|
1771 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1772 | IF (llm > 0) req_sz = req_sz*llm |
---|
1773 | IF (req_sz > var_sz) THEN |
---|
1774 | WRITE(str, & |
---|
1775 | & '("Size of variable put to the file should be ",I6)') req_sz |
---|
1776 | WRITE(str2,'("but the provided variable is of size ",I6)') var_sz |
---|
1777 | CALL ipslerr (3,'restput_r2d',str,str2,' ') |
---|
1778 | ENDIF |
---|
1779 | IF (req_sz < var_sz) THEN |
---|
1780 | WRITE(str,'("the size of variable put to the file is ",I6)') req_sz |
---|
1781 | WRITE(str2,'("but the provided variable is larger ",I6)') var_sz |
---|
1782 | CALL ipslerr (2,'restput_r2d', & |
---|
1783 | 'There could be a problem here :',str,str2) |
---|
1784 | ENDIF |
---|
1785 | !- |
---|
1786 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1787 | ! into the variable the model expects |
---|
1788 | !- |
---|
1789 | jl=0 |
---|
1790 | DO jj=1,siz2 |
---|
1791 | DO ji=1,siz1 |
---|
1792 | jl=jl+1 |
---|
1793 | buff_tmp1(jl) = var(ji,jj) |
---|
1794 | ENDDO |
---|
1795 | ENDDO |
---|
1796 | !- |
---|
1797 | CALL restput_real(fid,vname_q,iim,jjm,llm,itau,buff_tmp1) |
---|
1798 | !------------------------- |
---|
1799 | END SUBROUTINE restput_r2d |
---|
1800 | !=== |
---|
1801 | SUBROUTINE restput_r3d (fid,vname_q,iim,jjm,llm,itau,var) |
---|
1802 | !--------------------------------------------------------------------- |
---|
1803 | !- This subroutine serves as an interface to restput_real |
---|
1804 | !--------------------------------------------------------------------- |
---|
1805 | IMPLICIT NONE |
---|
1806 | !- |
---|
1807 | INTEGER :: fid |
---|
1808 | CHARACTER(LEN=*) :: vname_q |
---|
1809 | INTEGER :: iim,jjm,llm,itau |
---|
1810 | REAL :: var(:,:,:) |
---|
1811 | !- |
---|
1812 | ! LOCAL |
---|
1813 | !- |
---|
1814 | INTEGER :: ji,jj,jk,jl,req_sz,var_sz,siz1,siz2,siz3 |
---|
1815 | CHARACTER(LEN=70) :: str,str2 |
---|
1816 | LOGICAL :: check = .FALSE. |
---|
1817 | !--------------------------------------------------------------------- |
---|
1818 | !- |
---|
1819 | ! 1.0 Allocate the temporary buffer we need |
---|
1820 | ! to put the variable in right dimension |
---|
1821 | !- |
---|
1822 | siz1 = SIZE(var,1) |
---|
1823 | siz2 = SIZE(var,2) |
---|
1824 | siz3 = SIZE(var,3) |
---|
1825 | var_sz = siz1*siz2*siz3 |
---|
1826 | CALL rest_alloc (1,var_sz,check,'restput_r3d') |
---|
1827 | !- |
---|
1828 | ! 2.0 Here we could check if the sizes specified agree |
---|
1829 | ! with the size of the variable provided |
---|
1830 | !- |
---|
1831 | req_sz = 1 |
---|
1832 | IF (iim > 0) req_sz = req_sz*iim |
---|
1833 | IF (jjm > 0) req_sz = req_sz*jjm |
---|
1834 | IF (llm > 0) req_sz = req_sz*llm |
---|
1835 | IF (req_sz > var_sz) THEN |
---|
1836 | WRITE(str, & |
---|
1837 | & '("Size of variable put to the file should be ",I6)') req_sz |
---|
1838 | WRITE(str2, & |
---|
1839 | & '("but the provided variable is of size ",I6)') var_sz |
---|
1840 | CALL ipslerr (3,'restput_r3d',str,str2,' ') |
---|
1841 | ENDIF |
---|
1842 | IF (req_sz < var_sz) THEN |
---|
1843 | WRITE(str,'("the size of variable put to the file is ",I6)') req_sz |
---|
1844 | WRITE(str2,'("but the provided variable is larger ",I6)') var_sz |
---|
1845 | CALL ipslerr (2,'restput_r3d', & |
---|
1846 | 'There could be a problem here :',str,str2) |
---|
1847 | ENDIF |
---|
1848 | !- |
---|
1849 | ! 4.0 Transfer the buffer obtained from the restart file |
---|
1850 | ! into the variable the model expects |
---|
1851 | !- |
---|
1852 | jl=0 |
---|
1853 | DO jk=1,siz3 |
---|
1854 | DO jj=1,siz2 |
---|
1855 | DO ji=1,siz1 |
---|
1856 | jl=jl+1 |
---|
1857 | buff_tmp1(jl) = var(ji,jj,jk) |
---|
1858 | ENDDO |
---|
1859 | ENDDO |
---|
1860 | ENDDO |
---|
1861 | !- |
---|
1862 | CALL restput_real(fid,vname_q,iim,jjm,llm,itau,buff_tmp1) |
---|
1863 | !------------------------- |
---|
1864 | END SUBROUTINE restput_r3d |
---|
1865 | !=== |
---|
1866 | SUBROUTINE restput_real (fid,vname_q,iim,jjm,llm,itau,var) |
---|
1867 | !--------------------------------------------------------------------- |
---|
1868 | !- This subroutine will put a variable into the restart file. |
---|
1869 | !- But it will do a lot of other things if needed : |
---|
1870 | !- - Open a file if non is opened for this time-step |
---|
1871 | !- and all variables were written. |
---|
1872 | !- - Add an axis if needed |
---|
1873 | !- - verify that the variable has the right time step for this file |
---|
1874 | !- - If it is time for a new file then it is opened |
---|
1875 | !- and the old one closed |
---|
1876 | !- This requires that variables read from the last restart file were all |
---|
1877 | !- written |
---|
1878 | !- |
---|
1879 | !- INPUT |
---|
1880 | !- |
---|
1881 | !- fid : Id of the file in which we will write the variable |
---|
1882 | !- vname_q : Name of the variable to be written |
---|
1883 | !- iim,jjm,llm : Size in 3D of the variable |
---|
1884 | !- itau : Time step at which the variable is written |
---|
1885 | !- var : Variable |
---|
1886 | !- |
---|
1887 | !- OUTPUT |
---|
1888 | !- |
---|
1889 | !- NONE |
---|
1890 | !--------------------------------------------------------------------- |
---|
1891 | IMPLICIT NONE |
---|
1892 | !- |
---|
1893 | CHARACTER(LEN=*) :: vname_q |
---|
1894 | INTEGER :: fid,iim,jjm,llm,itau |
---|
1895 | REAL :: var(:) |
---|
1896 | !- |
---|
1897 | ! LOCAL |
---|
1898 | !- |
---|
1899 | INTEGER :: iret,vid,ncid,iv,vnb |
---|
1900 | INTEGER :: ierr |
---|
1901 | REAL :: secsince,one_day,one_year |
---|
1902 | INTEGER :: ndims |
---|
1903 | INTEGER,DIMENSION(4) :: corner,edge |
---|
1904 | !- |
---|
1905 | LOGICAL :: check = .FALSE. |
---|
1906 | !--------------------------------------------------------------------- |
---|
1907 | !- |
---|
1908 | ! 0.0 Get some variables |
---|
1909 | !- |
---|
1910 | ncid = netcdf_id(fid,2) |
---|
1911 | IF (netcdf_id(fid,2) < 0) THEN |
---|
1912 | CALL ipslerr (3,'restput', & |
---|
1913 | & 'The output restart file is undefined.',' ',' ') |
---|
1914 | ENDIF |
---|
1915 | CALL ioget_calendar (one_year,one_day) |
---|
1916 | !- |
---|
1917 | ! 1.0 Check if the variable is already present |
---|
1918 | !- |
---|
1919 | IF (check) WRITE(*,*) 'RESTPUT 1.0 : ',TRIM(vname_q) |
---|
1920 | !- |
---|
1921 | CALL find_str (varname_out(fid,1:nbvar_out(fid)),vname_q,vnb) |
---|
1922 | !- |
---|
1923 | IF (check) THEN |
---|
1924 | WRITE(*,*) 'RESTPUT 1.1 : ',varname_out(fid,1:nbvar_out(fid)),vnb |
---|
1925 | ENDIF |
---|
1926 | !- |
---|
1927 | ! 2.0 If variable is not present then declare it |
---|
1928 | ! and add extra dimensions if needed. |
---|
1929 | !- |
---|
1930 | IF (vnb <= 0) THEN |
---|
1931 | CALL restdefv (fid,vname_q,iim,jjm,llm,.TRUE.) |
---|
1932 | vnb = nbvar_out(fid) |
---|
1933 | ENDIF |
---|
1934 | vid = varid_out(fid,vnb) |
---|
1935 | !- |
---|
1936 | IF (check) WRITE(*,*) 'RESTPUT 2.0 : ',vnb,vid |
---|
1937 | !- |
---|
1938 | ! 2.1 Is this file already in write mode ? |
---|
1939 | ! If itau_out is still negative then we have |
---|
1940 | ! never written to it and we need to go into write mode. |
---|
1941 | !- |
---|
1942 | IF (itau_out(fid) < 0) THEN |
---|
1943 | iret = NF90_ENDDEF(ncid) |
---|
1944 | ENDIF |
---|
1945 | !- |
---|
1946 | ! 3.0 Is this itau already on the axis ? |
---|
1947 | ! If not then check that all variables of previous time is OK. |
---|
1948 | !- |
---|
1949 | IF (check) WRITE(*,*) 'RESTPUT 3.0 : ',itau,itau_out(fid) |
---|
1950 | !- |
---|
1951 | IF (itau /= itau_out(fid)) THEN |
---|
1952 | !--- |
---|
1953 | !-- If it is the first time step written on the restart |
---|
1954 | !-- then we only check the number |
---|
1955 | !-- Else we see if every variable was written |
---|
1956 | !--- |
---|
1957 | IF (tstp_out(fid) == 0) THEN |
---|
1958 | IF (nbvar_out(fid) < nbvar_read(fid)) THEN |
---|
1959 | WRITE(*,*) "ERROR :",tstp_out(fid), & |
---|
1960 | nbvar_out(fid),nbvar_read(fid) |
---|
1961 | CALL ipslerr (1,'restput', & |
---|
1962 | & 'There are fewer variables read from the output file', & |
---|
1963 | & 'than written onto the input file.', & |
---|
1964 | & 'We trust you know what you are doing') |
---|
1965 | ENDIF |
---|
1966 | ELSE |
---|
1967 | ierr = 0 |
---|
1968 | DO iv=1,nbvar_out(fid) |
---|
1969 | IF (.NOT.touched_out(fid,iv)) ierr = ierr+1 |
---|
1970 | ENDDO |
---|
1971 | IF (ierr > 0) THEN |
---|
1972 | WRITE(*,*) "ERROR :",nbvar_out(fid) |
---|
1973 | CALL ipslerr (1,'restput', & |
---|
1974 | & 'There are fewer variables in the output file for this', & |
---|
1975 | & 'time step than for the previous one',' ') |
---|
1976 | ELSE |
---|
1977 | touched_out(fid,:) = .FALSE. |
---|
1978 | ENDIF |
---|
1979 | ENDIF |
---|
1980 | !--- |
---|
1981 | secsince = itau*deltat(fid) |
---|
1982 | corner(1) = tstp_out(fid)+1 |
---|
1983 | edge(1) = 1 |
---|
1984 | !--- |
---|
1985 | !-- 3.1 Here we add the values to the time axes |
---|
1986 | !--- |
---|
1987 | IF (check) THEN |
---|
1988 | WRITE(*,*) 'RESTPUT 3.1 : ',itau,secsince,corner(1),edge(1) |
---|
1989 | ENDIF |
---|
1990 | !--- |
---|
1991 | iret = NF90_PUT_VAR(ncid,tind_varid_out(fid),itau, & |
---|
1992 | & start=corner(1:1)) |
---|
1993 | iret = NF90_PUT_VAR(ncid,tax_varid_out(fid),secsince, & |
---|
1994 | & start=corner(1:1)) |
---|
1995 | !--- |
---|
1996 | tstp_out(fid) = tstp_out(fid)+1 |
---|
1997 | itau_out(fid) = itau |
---|
1998 | ENDIF |
---|
1999 | !- |
---|
2000 | ! 4.0 Variable and time step should be present |
---|
2001 | ! now so we can dump variable |
---|
2002 | !- |
---|
2003 | ndims = 0 |
---|
2004 | IF (iim > 0) THEN |
---|
2005 | ndims = ndims+1 |
---|
2006 | corner(ndims) = 1 |
---|
2007 | edge(ndims) = iim |
---|
2008 | ENDIF |
---|
2009 | IF (jjm > 0) THEN |
---|
2010 | ndims = ndims+1 |
---|
2011 | corner(ndims) = 1 |
---|
2012 | edge(ndims) = jjm |
---|
2013 | ENDIF |
---|
2014 | IF (llm > 0) THEN |
---|
2015 | ndims = ndims+1 |
---|
2016 | corner(ndims) = 1 |
---|
2017 | edge(ndims) = llm |
---|
2018 | ENDIF |
---|
2019 | ndims = ndims+1 |
---|
2020 | corner(ndims) = tstp_out(fid) |
---|
2021 | edge(ndims) = 1 |
---|
2022 | !- |
---|
2023 | iret = NF90_PUT_VAR(ncid,vid,var, & |
---|
2024 | & start=corner(1:ndims),count=edge(1:ndims)) |
---|
2025 | !- |
---|
2026 | IF (iret /= NF90_NOERR) THEN |
---|
2027 | CALL ipslerr (2,'restput_real',NF90_STRERROR(iret), & |
---|
2028 | & 'Bug in restput.',& |
---|
2029 | & 'Please, verify compatibility between get and put commands.') |
---|
2030 | ENDIF |
---|
2031 | !- |
---|
2032 | ! 5.0 Note that the variables was treated |
---|
2033 | !- |
---|
2034 | touched_out(fid,vnb) = .TRUE. |
---|
2035 | !--------------------------- |
---|
2036 | END SUBROUTINE restput_real |
---|
2037 | !=== |
---|
2038 | SUBROUTINE restdefv (fid,varname,iim,jjm,llm,write_att) |
---|
2039 | !--------------------------------------------------------------------- |
---|
2040 | ! This subroutine adds a variable to the output file. |
---|
2041 | ! The attributes are either taken from. |
---|
2042 | !--------------------------------------------------------------------- |
---|
2043 | IMPLICIT NONE |
---|
2044 | !- |
---|
2045 | INTEGER ::fid |
---|
2046 | CHARACTER(LEN=*) :: varname |
---|
2047 | INTEGER :: iim,jjm,llm |
---|
2048 | LOGICAL :: write_att |
---|
2049 | !- |
---|
2050 | ! Local |
---|
2051 | !- |
---|
2052 | INTEGER :: dims(4),ic,xloc,ndim,ncfid |
---|
2053 | INTEGER :: iret,ax_id |
---|
2054 | CHARACTER(LEN=3) :: str |
---|
2055 | !- |
---|
2056 | LOGICAL :: check = .FALSE. |
---|
2057 | !--------------------------------------------------------------------- |
---|
2058 | ncfid = netcdf_id(fid,2) |
---|
2059 | IF (nbvar_out(fid) >= max_var) THEN |
---|
2060 | CALL ipslerr (3,'restdefv', & |
---|
2061 | 'Too many variables for the restcom module', & |
---|
2062 | 'Please increase the value of max_var',' ') |
---|
2063 | ENDIF |
---|
2064 | nbvar_out(fid) = nbvar_out(fid)+1 |
---|
2065 | varname_out(fid,nbvar_out(fid)) = varname |
---|
2066 | !- |
---|
2067 | ! 0.0 Put the file in define mode if needed |
---|
2068 | !- |
---|
2069 | IF (itau_out(fid) >= 0) THEN |
---|
2070 | iret = NF90_REDEF(ncfid) |
---|
2071 | ENDIF |
---|
2072 | !- |
---|
2073 | ! 1.0 Do we have all dimensions and can we go ahead |
---|
2074 | !- |
---|
2075 | IF (check) THEN |
---|
2076 | WRITE(*,*) 'restdefv 1.0 :',TRIM(varname),nbvar_out(fid) |
---|
2077 | ENDIF |
---|
2078 | !- |
---|
2079 | ndim = 0 |
---|
2080 | !- |
---|
2081 | ! 1.1 Work on x |
---|
2082 | !- |
---|
2083 | IF (iim > 0) THEN |
---|
2084 | ndim = ndim+1 |
---|
2085 | xloc = 0 |
---|
2086 | DO ic=1,xax_nb(fid) |
---|
2087 | IF (xax_infs(fid,ic,1) == iim) xloc = ic |
---|
2088 | ENDDO |
---|
2089 | !--- |
---|
2090 | IF (xloc > 0) THEN |
---|
2091 | dims(ndim) = xax_infs(fid,xloc,2) |
---|
2092 | ELSE |
---|
2093 | str='x_'//CHAR(96+xax_nb(fid)) |
---|
2094 | iret = NF90_DEF_DIM(ncfid,str,iim,ax_id) |
---|
2095 | xax_nb(fid) = xax_nb(fid)+1 |
---|
2096 | xax_infs(fid,xax_nb(fid),1) = iim |
---|
2097 | xax_infs(fid,xax_nb(fid),2) = ax_id |
---|
2098 | dims(ndim) = ax_id |
---|
2099 | ENDIF |
---|
2100 | ENDIF |
---|
2101 | !- |
---|
2102 | ! 1.2 Work on y |
---|
2103 | !- |
---|
2104 | IF (jjm > 0) THEN |
---|
2105 | ndim = ndim+1 |
---|
2106 | xloc = 0 |
---|
2107 | DO ic=1,yax_nb(fid) |
---|
2108 | IF (yax_infs(fid,ic,1) == jjm) xloc = ic |
---|
2109 | ENDDO |
---|
2110 | !--- |
---|
2111 | IF (xloc > 0) THEN |
---|
2112 | dims(ndim) = yax_infs(fid,xloc,2) |
---|
2113 | ELSE |
---|
2114 | str='y_'//CHAR(96+yax_nb(fid)) |
---|
2115 | iret = NF90_DEF_DIM(ncfid,str,jjm,ax_id) |
---|
2116 | yax_nb(fid) = yax_nb(fid)+1 |
---|
2117 | yax_infs(fid,yax_nb(fid),1) = jjm |
---|
2118 | yax_infs(fid,yax_nb(fid),2) = ax_id |
---|
2119 | dims(ndim) = ax_id |
---|
2120 | ENDIF |
---|
2121 | ENDIF |
---|
2122 | !- |
---|
2123 | ! 1.3 Work on z |
---|
2124 | !- |
---|
2125 | IF (llm > 0) THEN |
---|
2126 | ndim = ndim+1 |
---|
2127 | xloc = 0 |
---|
2128 | DO ic=1,zax_nb(fid) |
---|
2129 | IF (zax_infs(fid,ic,1) == llm) xloc = ic |
---|
2130 | ENDDO |
---|
2131 | !--- |
---|
2132 | IF (xloc > 0) THEN |
---|
2133 | dims(ndim) = zax_infs(fid,xloc,2) |
---|
2134 | ELSE |
---|
2135 | str='z_'//CHAR(96+zax_nb(fid)) |
---|
2136 | iret = NF90_DEF_DIM(ncfid,str,llm,ax_id) |
---|
2137 | zax_nb(fid) = zax_nb(fid)+1 |
---|
2138 | zax_infs(fid,zax_nb(fid),1) = llm |
---|
2139 | zax_infs(fid,zax_nb(fid),2) = ax_id |
---|
2140 | dims(ndim) = ax_id |
---|
2141 | ENDIF |
---|
2142 | ENDIF |
---|
2143 | !- |
---|
2144 | ! 1.4 Time needs to be added |
---|
2145 | !- |
---|
2146 | ndim = ndim+1 |
---|
2147 | dims(ndim) = tdimid_out(fid) |
---|
2148 | !- |
---|
2149 | ! 2.0 Declare the variable |
---|
2150 | !- |
---|
2151 | IF (check) THEN |
---|
2152 | WRITE(*,*) 'restdefv 2.0 :',ndim,' :: ',dims(1:ndim),tdimid_out(fid) |
---|
2153 | ENDIF |
---|
2154 | !- |
---|
2155 | iret = NF90_DEF_VAR(ncfid,varname,NF90_DOUBLE,dims(1:ndim), & |
---|
2156 | & varid_out(fid,nbvar_out(fid))) |
---|
2157 | IF (iret /= NF90_NOERR) THEN |
---|
2158 | CALL ipslerr (3,'restdefv', & |
---|
2159 | 'Could not define new variable in file', & |
---|
2160 | NF90_STRERROR(iret),varname) |
---|
2161 | ENDIF |
---|
2162 | !- |
---|
2163 | ! 3.0 Add the attributes if requested |
---|
2164 | !- |
---|
2165 | IF (write_att) THEN |
---|
2166 | IF (rest_units /= 'XXXXX') THEN |
---|
2167 | iret = NF90_PUT_ATT(ncfid,varid_out(fid,nbvar_out(fid)), & |
---|
2168 | & 'units',TRIM(rest_units)) |
---|
2169 | rest_units = 'XXXXX' |
---|
2170 | ENDIF |
---|
2171 | !--- |
---|
2172 | IF (rest_lname /= 'XXXXX') THEN |
---|
2173 | iret = NF90_PUT_ATT(ncfid,varid_out(fid,nbvar_out(fid)), & |
---|
2174 | & 'long_name',TRIM(rest_lname)) |
---|
2175 | rest_lname = 'XXXXX' |
---|
2176 | ENDIF |
---|
2177 | !--- |
---|
2178 | iret = NF90_PUT_ATT(ncfid,varid_out(fid,nbvar_out(fid)), & |
---|
2179 | & 'missing_value',REAL(missing_val,KIND=4)) |
---|
2180 | !--- |
---|
2181 | IF (itau_out(fid) >= 0) THEN |
---|
2182 | iret = NF90_ENDDEF(ncfid) |
---|
2183 | ENDIF |
---|
2184 | ENDIF |
---|
2185 | !- |
---|
2186 | IF (check) THEN |
---|
2187 | WRITE(*,*) & |
---|
2188 | & 'restdefv 3.0 : LIST OF VARS ',varname_out(fid,1:nbvar_out(fid)) |
---|
2189 | ENDIF |
---|
2190 | !---------------------- |
---|
2191 | END SUBROUTINE restdefv |
---|
2192 | !=== |
---|
2193 | SUBROUTINE rest_atim (l_msg,c_p) |
---|
2194 | !--------------------------------------------------------------------- |
---|
2195 | ! Called by "c_p", [re]allocate the time axes |
---|
2196 | !--------------------------------------------------------------------- |
---|
2197 | IMPLICIT NONE |
---|
2198 | !- |
---|
2199 | LOGICAL,INTENT(IN) :: l_msg |
---|
2200 | CHARACTER(LEN=*),INTENT(IN) :: c_p |
---|
2201 | !- |
---|
2202 | INTEGER :: i_err,tszij |
---|
2203 | INTEGER,ALLOCATABLE :: tmp_index(:,:) |
---|
2204 | REAL,ALLOCATABLE :: tmp_julian(:,:) |
---|
2205 | !--------------------------------------------------------------------- |
---|
2206 | !- |
---|
2207 | ! Allocate the space we need for the time axes |
---|
2208 | !- |
---|
2209 | IF (.NOT.ALLOCATED(t_index) .AND. .NOT.ALLOCATED(t_julian)) THEN |
---|
2210 | IF (l_msg) THEN |
---|
2211 | WRITE(*,*) TRIM(c_p)//' : Allocate times axes at :', & |
---|
2212 | & max_file,tax_size_in(nb_fi) |
---|
2213 | ENDIF |
---|
2214 | !--- |
---|
2215 | ALLOCATE(t_index(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2216 | IF (i_err/=0) THEN |
---|
2217 | WRITE(*,*) "ERROR IN ALLOCATION of t_index : ",i_err |
---|
2218 | CALL ipslerr (3,TRIM(c_p), & |
---|
2219 | & 'Problem in allocation of t_index','', & |
---|
2220 | & '(you must increase memory)') |
---|
2221 | ENDIF |
---|
2222 | t_index (:,:) = 0 |
---|
2223 | !--- |
---|
2224 | ALLOCATE(t_julian(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2225 | IF (i_err/=0) THEN |
---|
2226 | WRITE(*,*) "ERROR IN ALLOCATION of t_julian : ",i_err |
---|
2227 | CALL ipslerr (3,TRIM(c_p), & |
---|
2228 | & 'Problem in allocation of max_file,tax_size_in','', & |
---|
2229 | & '(you must increase memory)') |
---|
2230 | ENDIF |
---|
2231 | t_julian (:,:) = 0.0 |
---|
2232 | ELSE IF ( (SIZE(t_index,DIM=2) < tax_size_in(nb_fi)) & |
---|
2233 | & .OR.(SIZE(t_julian,DIM=2) < tax_size_in(nb_fi)) ) THEN |
---|
2234 | IF (l_msg) THEN |
---|
2235 | WRITE(*,*) TRIM(c_p)//' : Reallocate times axes at :', & |
---|
2236 | & max_file,tax_size_in(nb_fi) |
---|
2237 | ENDIF |
---|
2238 | !--- |
---|
2239 | ALLOCATE (tmp_index(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2240 | IF (i_err/=0) THEN |
---|
2241 | WRITE(*,*) "ERROR IN ALLOCATION of tmp_index : ",i_err |
---|
2242 | CALL ipslerr (3,TRIM(c_p), & |
---|
2243 | & 'Problem in allocation of tmp_index','', & |
---|
2244 | & '(you must increase memory)') |
---|
2245 | ENDIF |
---|
2246 | tszij = SIZE(t_index,DIM=2) |
---|
2247 | tmp_index(:,1:tszij) = t_index(:,1:tszij) |
---|
2248 | DEALLOCATE(t_index) |
---|
2249 | ALLOCATE (t_index(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2250 | IF (i_err/=0) THEN |
---|
2251 | WRITE(*,*) "ERROR IN ALLOCATION of t_index : ",i_err |
---|
2252 | CALL ipslerr (3,TRIM(c_p), & |
---|
2253 | & 'Problem in reallocation of t_index','', & |
---|
2254 | & '(you must increase memory)') |
---|
2255 | ENDIF |
---|
2256 | t_index(:,1:tszij) = tmp_index(:,1:tszij) |
---|
2257 | !--- |
---|
2258 | ALLOCATE (tmp_julian(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2259 | IF (i_err/=0) THEN |
---|
2260 | WRITE(*,*) "ERROR IN ALLOCATION of tmp_julian : ",i_err |
---|
2261 | CALL ipslerr (3,TRIM(c_p), & |
---|
2262 | & 'Problem in allocation of tmp_julian','', & |
---|
2263 | & '(you must increase memory)') |
---|
2264 | ENDIF |
---|
2265 | tszij = SIZE(t_julian,DIM=2) |
---|
2266 | tmp_julian(:,1:tszij) = t_julian(:,1:tszij) |
---|
2267 | DEALLOCATE(t_julian) |
---|
2268 | ALLOCATE (t_julian(max_file,tax_size_in(nb_fi)),STAT=i_err) |
---|
2269 | IF (i_err/=0) THEN |
---|
2270 | WRITE(*,*) "ERROR IN ALLOCATION of t_julian : ",i_err |
---|
2271 | CALL ipslerr (3,TRIM(c_p), & |
---|
2272 | & 'Problem in reallocation of t_julian','', & |
---|
2273 | & '(you must increase memory)') |
---|
2274 | ENDIF |
---|
2275 | t_julian(:,1:tszij) = tmp_julian(:,1:tszij) |
---|
2276 | ENDIF |
---|
2277 | !----------------------- |
---|
2278 | END SUBROUTINE rest_atim |
---|
2279 | !=== |
---|
2280 | SUBROUTINE rest_alloc (i_buff,i_qsz,l_msg,c_p) |
---|
2281 | !--------------------------------------------------------------------- |
---|
2282 | ! Called by "c_p", allocate a temporary buffer |
---|
2283 | ! (buff_tmp[1/2] depending on "i_buff" value) to the size "i_qsz". |
---|
2284 | !--------------------------------------------------------------------- |
---|
2285 | IMPLICIT NONE |
---|
2286 | !- |
---|
2287 | INTEGER,INTENT(IN) :: i_buff,i_qsz |
---|
2288 | LOGICAL,INTENT(IN) :: l_msg |
---|
2289 | CHARACTER(LEN=*),INTENT(IN) :: c_p |
---|
2290 | !- |
---|
2291 | INTEGER :: i_bsz,i_err |
---|
2292 | LOGICAL :: l_alloc1,l_alloc2 |
---|
2293 | CHARACTER(LEN=9) :: cbn |
---|
2294 | CHARACTER(LEN=5) :: c_err |
---|
2295 | !--------------------------------------------------------------------- |
---|
2296 | IF (i_buff == 1) THEN |
---|
2297 | IF (ALLOCATED(buff_tmp1)) THEN |
---|
2298 | i_bsz = SIZE(buff_tmp1) |
---|
2299 | ELSE |
---|
2300 | i_bsz = 0 |
---|
2301 | ENDIF |
---|
2302 | l_alloc1 = (.NOT.ALLOCATED(buff_tmp1)) & |
---|
2303 | & .OR.((ALLOCATED(buff_tmp1)).AND.(i_qsz > i_bsz)) |
---|
2304 | l_alloc2 = .FALSE. |
---|
2305 | cbn = 'buff_tmp1' |
---|
2306 | ELSE IF (i_buff == 2) THEN |
---|
2307 | IF (ALLOCATED(buff_tmp2)) THEN |
---|
2308 | i_bsz = SIZE(buff_tmp2) |
---|
2309 | ELSE |
---|
2310 | i_bsz = 0 |
---|
2311 | ENDIF |
---|
2312 | l_alloc1 = .FALSE. |
---|
2313 | l_alloc2 = (.NOT.ALLOCATED(buff_tmp2)) & |
---|
2314 | & .OR.((ALLOCATED(buff_tmp2)).AND.(i_qsz > i_bsz)) |
---|
2315 | cbn = 'buff_tmp2' |
---|
2316 | ELSE |
---|
2317 | CALL ipslerr (3,'rest_alloc', & |
---|
2318 | & 'Called by '//TRIM(c_p),'with a wrong value of i_buff','') |
---|
2319 | ENDIF |
---|
2320 | !- |
---|
2321 | !- |
---|
2322 | IF (l_alloc1.OR.l_alloc2) THEN |
---|
2323 | IF (l_msg) THEN |
---|
2324 | IF ( (l_alloc1.AND.ALLOCATED(buff_tmp1)) & |
---|
2325 | & .OR.(l_alloc2.AND.ALLOCATED(buff_tmp2)) ) THEN |
---|
2326 | WRITE(*,*) TRIM(c_p)//" : re_allocate "//TRIM(cbn)//"=",i_qsz |
---|
2327 | ELSE |
---|
2328 | WRITE(*,*) TRIM(c_p)//" : allocate "//TRIM(cbn)//"=",i_qsz |
---|
2329 | ENDIF |
---|
2330 | ENDIF |
---|
2331 | IF (l_alloc1) THEN |
---|
2332 | IF (ALLOCATED(buff_tmp1)) THEN |
---|
2333 | DEALLOCATE(buff_tmp1) |
---|
2334 | ENDIF |
---|
2335 | ALLOCATE (buff_tmp1(i_qsz),STAT=i_err) |
---|
2336 | ELSE |
---|
2337 | IF (ALLOCATED(buff_tmp2)) THEN |
---|
2338 | DEALLOCATE(buff_tmp2) |
---|
2339 | ENDIF |
---|
2340 | ALLOCATE (buff_tmp2(i_qsz),STAT=i_err) |
---|
2341 | ENDIF |
---|
2342 | IF (i_err /= 0) THEN |
---|
2343 | WRITE (UNIT=c_err,FMT='(I5)') i_err |
---|
2344 | CALL ipslerr (3,TRIM(c_p), & |
---|
2345 | & 'Problem in allocation of',TRIM(cbn), & |
---|
2346 | & 'Error : '//TRIM(c_err)//' (you must increase memory)') |
---|
2347 | ENDIF |
---|
2348 | ENDIF |
---|
2349 | !------------------------ |
---|
2350 | END SUBROUTINE rest_alloc |
---|
2351 | !=== |
---|
2352 | SUBROUTINE ioconf_setatt (attname,value) |
---|
2353 | !--------------------------------------------------------------------- |
---|
2354 | IMPLICIT NONE |
---|
2355 | !- |
---|
2356 | CHARACTER(LEN=*) :: attname,value |
---|
2357 | !- |
---|
2358 | ! LOCAL |
---|
2359 | !- |
---|
2360 | CHARACTER(LEN=LEN_TRIM(attname)) :: tmp_str |
---|
2361 | !--------------------------------------------------------------------- |
---|
2362 | tmp_str = attname |
---|
2363 | CALL strlowercase (tmp_str) |
---|
2364 | !- |
---|
2365 | SELECT CASE(tmp_str) |
---|
2366 | CASE('units') |
---|
2367 | rest_units = value |
---|
2368 | CASE('long_name') |
---|
2369 | rest_lname = value |
---|
2370 | CASE DEFAULT |
---|
2371 | CALL ipslerr (2,'ioconf_restatt', & |
---|
2372 | 'The attribute name provided is unknown',attname,' ') |
---|
2373 | END SELECT |
---|
2374 | !--------------------------- |
---|
2375 | END SUBROUTINE ioconf_setatt |
---|
2376 | !=== |
---|
2377 | SUBROUTINE ioget_vdim (fid,vname_q,varnbdim_max,varnbdim,vardims) |
---|
2378 | !--------------------------------------------------------------------- |
---|
2379 | !- This routine allows the user to get the dimensions |
---|
2380 | !- of a field in the restart file. |
---|
2381 | !- This is the file which is read. |
---|
2382 | !--------------------------------------------------------------------- |
---|
2383 | IMPLICIT NONE |
---|
2384 | !- |
---|
2385 | INTEGER,INTENT(IN) :: fid |
---|
2386 | CHARACTER(LEN=*) :: vname_q |
---|
2387 | INTEGER,INTENT(IN) :: varnbdim_max |
---|
2388 | INTEGER,INTENT(OUT) :: varnbdim |
---|
2389 | INTEGER,DIMENSION(varnbdim_max),INTENT(OUT) :: vardims |
---|
2390 | !- |
---|
2391 | ! LOCAL |
---|
2392 | !- |
---|
2393 | INTEGER :: vnb |
---|
2394 | !--------------------------------------------------------------------- |
---|
2395 | ! Find the index of the variable |
---|
2396 | CALL find_str (varname_in(fid,1:nbvar_in(fid)),vname_q,vnb) |
---|
2397 | !- |
---|
2398 | IF (vnb > 0) THEN |
---|
2399 | varnbdim = varnbdim_in(fid,vnb) |
---|
2400 | IF (varnbdim_max < varnbdim) THEN |
---|
2401 | CALL ipslerr (3,'ioget_vdim', & |
---|
2402 | 'The provided array for the variable dimensions is too small', & |
---|
2403 | '','') |
---|
2404 | ELSE |
---|
2405 | vardims(1:varnbdim) = vardims_in(fid,vnb,1:varnbdim) |
---|
2406 | ENDIF |
---|
2407 | ELSE |
---|
2408 | varnbdim = 0 |
---|
2409 | CALL ipslerr (2,'ioget_vdim', & |
---|
2410 | 'Variable '//TRIM(vname_q)//' not found','','') |
---|
2411 | ENDIF |
---|
2412 | !------------------------ |
---|
2413 | END SUBROUTINE ioget_vdim |
---|
2414 | !=== |
---|
2415 | SUBROUTINE ioget_vname (fid,nbvar,varnames) |
---|
2416 | !--------------------------------------------------------------------- |
---|
2417 | !- This routine allows the user to extract the list |
---|
2418 | !- of variables in an opened restart file. |
---|
2419 | !- This is the file which is read |
---|
2420 | !--------------------------------------------------------------------- |
---|
2421 | IMPLICIT NONE |
---|
2422 | !- |
---|
2423 | INTEGER,INTENT(IN) :: fid |
---|
2424 | INTEGER,INTENT(OUT) :: nbvar |
---|
2425 | CHARACTER(LEN=*),INTENT(OUT) :: varnames(:) |
---|
2426 | !--------------------------------------------------------------------- |
---|
2427 | nbvar = nbvar_in(fid) |
---|
2428 | !- |
---|
2429 | IF (SIZE(varnames) < nbvar) THEN |
---|
2430 | CALL ipslerr (3,'ioget_vname', & |
---|
2431 | 'The provided array for the variable names is too small','','') |
---|
2432 | ELSE |
---|
2433 | varnames(1:nbvar) = varname_in(fid,1:nbvar) |
---|
2434 | ENDIF |
---|
2435 | !------------------------- |
---|
2436 | END SUBROUTINE ioget_vname |
---|
2437 | !=== |
---|
2438 | SUBROUTINE ioconf_expval (new_exp_val) |
---|
2439 | !--------------------------------------------------------------------- |
---|
2440 | !- The default value written into the variables which are not |
---|
2441 | !- in the restart file can only be changed once. |
---|
2442 | !- This avoids further complications. |
---|
2443 | !--------------------------------------------------------------------- |
---|
2444 | IMPLICIT NONE |
---|
2445 | !- |
---|
2446 | REAL :: new_exp_val |
---|
2447 | !--------------------------------------------------------------------- |
---|
2448 | IF (.NOT.lock_valexp) THEN |
---|
2449 | lock_valexp = .TRUE. |
---|
2450 | val_exp = new_exp_val |
---|
2451 | ELSE |
---|
2452 | CALL ipslerr (2,'ioconf_expval', & |
---|
2453 | 'The default value for variable' & |
---|
2454 | //'not available in the restart file ', & |
---|
2455 | 'has already been locked and can not be changed at this point', & |
---|
2456 | ' ') |
---|
2457 | ENDIF |
---|
2458 | !--------------------------- |
---|
2459 | END SUBROUTINE ioconf_expval |
---|
2460 | !=== |
---|
2461 | SUBROUTINE ioget_expval (get_exp_val) |
---|
2462 | !--------------------------------------------------------------------- |
---|
2463 | !- Once the user has extracted the default value, |
---|
2464 | !- we lock it so that it can not be changed anymore. |
---|
2465 | !--------------------------------------------------------------------- |
---|
2466 | IMPLICIT NONE |
---|
2467 | !- |
---|
2468 | REAL :: get_exp_val |
---|
2469 | !--------------------------------------------------------------------- |
---|
2470 | get_exp_val = val_exp |
---|
2471 | lock_valexp = .TRUE. |
---|
2472 | !-------------------------- |
---|
2473 | END SUBROUTINE ioget_expval |
---|
2474 | !=== |
---|
2475 | SUBROUTINE restclo (fid) |
---|
2476 | !--------------------------------------------------------------------- |
---|
2477 | !- This subroutine closes one or any opened restart file. |
---|
2478 | !- |
---|
2479 | !- INPUT |
---|
2480 | !- |
---|
2481 | !- fid : File ID in the restcom system (not the netCDF ID)(optional) |
---|
2482 | !- |
---|
2483 | !- OUTPUT |
---|
2484 | !- |
---|
2485 | !- NONE |
---|
2486 | !--------------------------------------------------------------------- |
---|
2487 | IMPLICIT NONE |
---|
2488 | !- |
---|
2489 | INTEGER,INTENT(in),OPTIONAL :: fid |
---|
2490 | !- |
---|
2491 | !- LOCAL |
---|
2492 | !- |
---|
2493 | INTEGER :: iret,ifnc |
---|
2494 | CHARACTER(LEN=6) :: n_e |
---|
2495 | CHARACTER(LEN=3) :: n_f |
---|
2496 | LOGICAL :: check = .FALSE. |
---|
2497 | !--------------------------------------------------------------------- |
---|
2498 | IF (PRESENT(fid)) THEN |
---|
2499 | !--- |
---|
2500 | IF (check) THEN |
---|
2501 | WRITE(*,*) & |
---|
2502 | 'restclo : Closing specified restart file number :', & |
---|
2503 | fid,netcdf_id(fid,1:2) |
---|
2504 | ENDIF |
---|
2505 | !--- |
---|
2506 | IF (netcdf_id(fid,1) > 0) THEN |
---|
2507 | iret = NF90_CLOSE(netcdf_id(fid,1)) |
---|
2508 | IF (iret /= NF90_NOERR) THEN |
---|
2509 | WRITE (n_e,'(I6)') iret |
---|
2510 | WRITE (n_f,'(I3)') netcdf_id(fid,1) |
---|
2511 | CALL ipslerr (2,'restclo', & |
---|
2512 | "Error "//n_e//" in closing file : "//n_f,'',' ') |
---|
2513 | ENDIF |
---|
2514 | IF (netcdf_id(fid,1) == netcdf_id(fid,2)) THEN |
---|
2515 | netcdf_id(fid,2) = -1 |
---|
2516 | ENDIF |
---|
2517 | netcdf_id(fid,1) = -1 |
---|
2518 | ENDIF |
---|
2519 | !--- |
---|
2520 | IF (netcdf_id(fid,2) > 0) THEN |
---|
2521 | iret = NF90_CLOSE(netcdf_id(fid,2)) |
---|
2522 | IF (iret /= NF90_NOERR) THEN |
---|
2523 | WRITE (n_e,'(I6)') iret |
---|
2524 | WRITE (n_f,'(I3)') netcdf_id(fid,2) |
---|
2525 | CALL ipslerr (2,'restclo', & |
---|
2526 | "Error "//n_e//" in closing file : "//n_f,'',' ') |
---|
2527 | ENDIF |
---|
2528 | netcdf_id(fid,2) = -1 |
---|
2529 | ENDIF |
---|
2530 | !--- |
---|
2531 | ELSE |
---|
2532 | !--- |
---|
2533 | IF (check) WRITE(*,*) 'restclo : Closing all files' |
---|
2534 | !--- |
---|
2535 | DO ifnc=1,nb_fi |
---|
2536 | IF (netcdf_id(ifnc,1) > 0) THEN |
---|
2537 | iret = NF90_CLOSE(netcdf_id(ifnc,1)) |
---|
2538 | IF (iret /= NF90_NOERR) THEN |
---|
2539 | WRITE (n_e,'(I6)') iret |
---|
2540 | WRITE (n_f,'(I3)') netcdf_id(ifnc,1) |
---|
2541 | CALL ipslerr (2,'restclo', & |
---|
2542 | "Error "//n_e//" in closing file : "//n_f,'',' ') |
---|
2543 | ENDIF |
---|
2544 | IF (netcdf_id(ifnc,1) == netcdf_id(ifnc,2)) THEN |
---|
2545 | netcdf_id(ifnc,2) = -1 |
---|
2546 | ENDIF |
---|
2547 | netcdf_id(ifnc,1) = -1 |
---|
2548 | ENDIF |
---|
2549 | !----- |
---|
2550 | IF (netcdf_id(ifnc,2) > 0) THEN |
---|
2551 | iret = NF90_CLOSE(netcdf_id(ifnc,2)) |
---|
2552 | IF (iret /= NF90_NOERR) THEN |
---|
2553 | WRITE (n_e,'(I6)') iret |
---|
2554 | WRITE (n_f,'(I3)') netcdf_id(ifnc,2) |
---|
2555 | CALL ipslerr (2,'restclo', & |
---|
2556 | "Error "//n_e//" in closing file : "//n_f,'',' ') |
---|
2557 | END IF |
---|
2558 | netcdf_id(ifnc,2) = -1 |
---|
2559 | ENDIF |
---|
2560 | ENDDO |
---|
2561 | ENDIF |
---|
2562 | !--------------------- |
---|
2563 | END SUBROUTINE restclo |
---|
2564 | !=== |
---|
2565 | !----------------- |
---|
2566 | END MODULE restcom |
---|