1 | MODULE iostart |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | PRIVATE |
---|
5 | INTEGER,SAVE :: nid_start ! NetCDF file identifier for startfi.nc file |
---|
6 | INTEGER,SAVE :: nid_restart ! NetCDF file identifier for restartfi.nc file |
---|
7 | |
---|
8 | ! restartfi.nc file dimension identifiers: (see open_restartphy()) |
---|
9 | INTEGER,SAVE :: idim1 ! "index" dimension |
---|
10 | INTEGER,SAVE :: idim2 ! "physical_points" dimension |
---|
11 | INTEGER,SAVE :: idim3 ! "subsurface_layers" dimension |
---|
12 | INTEGER,SAVE :: idim4 ! "nlayer_plus_1" dimension |
---|
13 | INTEGER,SAVE :: idim5 ! "number_of_advected_fields" dimension |
---|
14 | INTEGER,SAVE :: idim6 ! "nlayer" dimension |
---|
15 | INTEGER,SAVE :: idim7 ! "Time" dimension |
---|
16 | INTEGER,SAVE :: timeindex ! current time index (for time-dependent fields) |
---|
17 | INTEGER,PARAMETER :: length=100 ! size of tab_cntrl array |
---|
18 | |
---|
19 | INTERFACE get_field |
---|
20 | MODULE PROCEDURE Get_field_r1,Get_field_r2,Get_field_r3 |
---|
21 | END INTERFACE get_field |
---|
22 | |
---|
23 | INTERFACE get_var |
---|
24 | MODULE PROCEDURE get_var_r0,Get_var_r1,Get_var_r2,Get_var_r3 |
---|
25 | END INTERFACE get_var |
---|
26 | |
---|
27 | INTERFACE put_field |
---|
28 | MODULE PROCEDURE put_field_r1,put_field_r2,put_field_r3 |
---|
29 | END INTERFACE put_field |
---|
30 | |
---|
31 | INTERFACE put_var |
---|
32 | MODULE PROCEDURE put_var_r0,put_var_r1,put_var_r2,put_var_r3 |
---|
33 | END INTERFACE put_var |
---|
34 | |
---|
35 | PUBLIC nid_start, length |
---|
36 | PUBLIC get_field,get_var,put_field,put_var |
---|
37 | PUBLIC inquire_dimension, inquire_dimension_length |
---|
38 | PUBLIC inquire_field, inquire_field_ndims |
---|
39 | PUBLIC open_startphy,close_startphy,open_restartphy,close_restartphy |
---|
40 | |
---|
41 | CONTAINS |
---|
42 | |
---|
43 | SUBROUTINE open_startphy(filename) |
---|
44 | USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_NOWRITE, nf90_strerror |
---|
45 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
46 | IMPLICIT NONE |
---|
47 | CHARACTER(LEN=*) :: filename |
---|
48 | INTEGER :: ierr |
---|
49 | |
---|
50 | IF (is_master) THEN |
---|
51 | ierr = NF90_OPEN (filename, NF90_NOWRITE, nid_start) |
---|
52 | IF (ierr.NE.NF90_NOERR) THEN |
---|
53 | write(*,*)'open_startphy: problem opening file '//trim(filename) |
---|
54 | write(*,*)trim(nf90_strerror(ierr)) |
---|
55 | CALL abort_physic("open_startphy","Cannot open file",1) |
---|
56 | ENDIF |
---|
57 | ENDIF |
---|
58 | |
---|
59 | CALL bcast(nid_start) ! tell all procs about nid_start |
---|
60 | |
---|
61 | END SUBROUTINE open_startphy |
---|
62 | |
---|
63 | SUBROUTINE close_startphy |
---|
64 | USE netcdf, only: NF90_CLOSE |
---|
65 | USE mod_phys_lmdz_para, only: is_master |
---|
66 | IMPLICIT NONE |
---|
67 | INTEGER :: ierr |
---|
68 | |
---|
69 | IF (is_master) THEN |
---|
70 | ierr = NF90_CLOSE (nid_start) |
---|
71 | ENDIF |
---|
72 | |
---|
73 | END SUBROUTINE close_startphy |
---|
74 | |
---|
75 | |
---|
76 | FUNCTION inquire_field(Field_name) |
---|
77 | ! check if a given field is present in the input file |
---|
78 | USE netcdf, only: NF90_INQ_VARID, NF90_NOERR |
---|
79 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
80 | IMPLICIT NONE |
---|
81 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
82 | LOGICAL :: inquire_field |
---|
83 | INTEGER :: varid |
---|
84 | INTEGER :: ierr |
---|
85 | |
---|
86 | IF (is_master) THEN |
---|
87 | ierr=NF90_INQ_VARID(nid_start,Field_name,varid) |
---|
88 | IF (ierr==NF90_NOERR) THEN |
---|
89 | Inquire_field=.TRUE. |
---|
90 | ELSE |
---|
91 | Inquire_field=.FALSE. |
---|
92 | ENDIF |
---|
93 | ENDIF |
---|
94 | |
---|
95 | CALL bcast(inquire_field) |
---|
96 | |
---|
97 | END FUNCTION inquire_field |
---|
98 | |
---|
99 | |
---|
100 | FUNCTION inquire_field_ndims(Field_name) |
---|
101 | ! give the number of dimensions of "Field_name" stored in the input file |
---|
102 | USE netcdf, only: nf90_inq_varid, nf90_inquire_variable, & |
---|
103 | NF90_NOERR, nf90_strerror |
---|
104 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
105 | IMPLICIT NONE |
---|
106 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
107 | INTEGER :: inquire_field_ndims |
---|
108 | INTEGER :: varid |
---|
109 | INTEGER :: ierr |
---|
110 | |
---|
111 | IF (is_master) THEN |
---|
112 | ierr=nf90_inq_varid(nid_start,Field_name,varid) |
---|
113 | ierr=nf90_inquire_variable(nid_start,varid,& |
---|
114 | ndims=inquire_field_ndims) |
---|
115 | IF (ierr.NE.NF90_NOERR) THEN |
---|
116 | write(*,*)'inquire_field_ndims: problem obtaining ndims of '& |
---|
117 | //trim(field_name) |
---|
118 | write(*,*)trim(nf90_strerror(ierr)) |
---|
119 | CALL abort_physic("inquire_field_ndims","Failed to get ndims",1) |
---|
120 | ENDIF |
---|
121 | ENDIF |
---|
122 | |
---|
123 | CALL bcast(inquire_field_ndims) |
---|
124 | |
---|
125 | END FUNCTION inquire_field_ndims |
---|
126 | |
---|
127 | |
---|
128 | FUNCTION inquire_dimension(Field_name) |
---|
129 | ! check if a given dimension is present in the input file |
---|
130 | USE netcdf, only: nf90_inq_dimid, NF90_NOERR |
---|
131 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
132 | IMPLICIT NONE |
---|
133 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
134 | LOGICAL :: inquire_dimension |
---|
135 | INTEGER :: varid |
---|
136 | INTEGER :: ierr |
---|
137 | |
---|
138 | IF (is_master) THEN |
---|
139 | ierr=NF90_INQ_DIMID(nid_start,Field_name,varid) |
---|
140 | IF (ierr==NF90_NOERR) THEN |
---|
141 | Inquire_dimension=.TRUE. |
---|
142 | ELSE |
---|
143 | Inquire_dimension=.FALSE. |
---|
144 | ENDIF |
---|
145 | ENDIF |
---|
146 | |
---|
147 | CALL bcast(inquire_dimension) |
---|
148 | |
---|
149 | END FUNCTION inquire_dimension |
---|
150 | |
---|
151 | FUNCTION inquire_dimension_length(Field_name) |
---|
152 | ! give the length of the "Field_name" dimension stored in the input file |
---|
153 | USE netcdf, only: nf90_inquire_dimension, nf90_inq_dimid, & |
---|
154 | NF90_NOERR, nf90_strerror |
---|
155 | USE mod_phys_lmdz_para, only: is_master, bcast |
---|
156 | IMPLICIT NONE |
---|
157 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
158 | INTEGER :: inquire_dimension_length |
---|
159 | INTEGER :: varid |
---|
160 | INTEGER :: ierr |
---|
161 | |
---|
162 | IF (is_master) THEN |
---|
163 | ierr=nf90_inq_dimid(nid_start,Field_name,varid) |
---|
164 | ierr=nf90_inquire_dimension(nid_start,varid,& |
---|
165 | len=inquire_dimension_length) |
---|
166 | IF (ierr.NE.NF90_NOERR) THEN |
---|
167 | write(*,*)'inquire_field_length: problem obtaining length of '& |
---|
168 | //trim(field_name) |
---|
169 | write(*,*)trim(nf90_strerror(ierr)) |
---|
170 | CALL abort_physic("inquire_field_ndims","Failed to get length",1) |
---|
171 | ENDIF |
---|
172 | ENDIF |
---|
173 | |
---|
174 | CALL bcast(inquire_dimension_length) |
---|
175 | |
---|
176 | END FUNCTION inquire_dimension_length |
---|
177 | |
---|
178 | |
---|
179 | |
---|
180 | SUBROUTINE Get_Field_r1(field_name,field,found,timeindex) |
---|
181 | ! For a surface field |
---|
182 | use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid) |
---|
183 | IMPLICIT NONE |
---|
184 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
185 | REAL,INTENT(INOUT) :: Field(:) |
---|
186 | LOGICAL,INTENT(OUT),OPTIONAL :: found |
---|
187 | INTEGER,INTENT(IN),OPTIONAL :: timeindex ! time index of sought data |
---|
188 | |
---|
189 | integer :: edges(4), corners(4) |
---|
190 | |
---|
191 | edges(1)=klon_glo |
---|
192 | edges(2:4)=1 |
---|
193 | corners(1:4)=1 |
---|
194 | if (PRESENT(timeindex)) then |
---|
195 | corners(2)=timeindex |
---|
196 | endif |
---|
197 | |
---|
198 | IF (PRESENT(found)) THEN |
---|
199 | CALL Get_field_rgen(field_name,field,1,corners,edges,found) |
---|
200 | ELSE |
---|
201 | CALL Get_field_rgen(field_name,field,1,corners,edges) |
---|
202 | ENDIF |
---|
203 | |
---|
204 | END SUBROUTINE Get_Field_r1 |
---|
205 | |
---|
206 | SUBROUTINE Get_Field_r2(field_name,field,found,timeindex) |
---|
207 | ! For a "3D" horizontal-vertical field |
---|
208 | use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid) |
---|
209 | IMPLICIT NONE |
---|
210 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
211 | REAL,INTENT(INOUT) :: Field(:,:) |
---|
212 | LOGICAL,INTENT(OUT),OPTIONAL :: found |
---|
213 | INTEGER,INTENT(IN),OPTIONAL :: timeindex ! time index of sought data |
---|
214 | |
---|
215 | integer :: edges(4), corners(4) |
---|
216 | |
---|
217 | edges(1)=klon_glo |
---|
218 | edges(2)=size(field,2) |
---|
219 | edges(3:4)=1 |
---|
220 | corners(1:4)=1 |
---|
221 | if (PRESENT(timeindex)) then |
---|
222 | corners(3)=timeindex |
---|
223 | endif |
---|
224 | |
---|
225 | IF (PRESENT(found)) THEN |
---|
226 | CALL Get_field_rgen(field_name,field,size(field,2),& |
---|
227 | corners,edges,found) |
---|
228 | ELSE |
---|
229 | CALL Get_field_rgen(field_name,field,size(field,2),& |
---|
230 | corners,edges) |
---|
231 | ENDIF |
---|
232 | |
---|
233 | |
---|
234 | END SUBROUTINE Get_Field_r2 |
---|
235 | |
---|
236 | SUBROUTINE Get_Field_r3(field_name,field,found,timeindex) |
---|
237 | ! for a "4D" field surf/alt/?? |
---|
238 | use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid) |
---|
239 | IMPLICIT NONE |
---|
240 | CHARACTER(LEN=*),INTENT(IN) :: Field_name |
---|
241 | REAL,INTENT(INOUT) :: Field(:,:,:) |
---|
242 | LOGICAL,INTENT(OUT),OPTIONAL :: found |
---|
243 | INTEGER,INTENT(IN),OPTIONAL :: timeindex ! time index of sought data |
---|
244 | |
---|
245 | integer :: edges(4), corners(4) |
---|
246 | |
---|
247 | edges(1)=klon_glo |
---|
248 | edges(2)=size(field,2) |
---|
249 | edges(3)=size(field,3) |
---|
250 | edges(4)=1 |
---|
251 | corners(1:4)=1 |
---|
252 | if (PRESENT(timeindex)) then |
---|
253 | corners(4)=timeindex |
---|
254 | endif |
---|
255 | |
---|
256 | IF (PRESENT(found)) THEN |
---|
257 | CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),& |
---|
258 | corners,edges,found) |
---|
259 | ELSE |
---|
260 | CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),& |
---|
261 | corners,edges) |
---|
262 | ENDIF |
---|
263 | |
---|
264 | END SUBROUTINE Get_Field_r3 |
---|
265 | |
---|
266 | SUBROUTINE Get_field_rgen(field_name,field,field_size, & |
---|
267 | corners,edges,found) |
---|
268 | USE netcdf |
---|
269 | USE dimphy |
---|
270 | USE mod_grid_phy_lmdz |
---|
271 | USE mod_phys_lmdz_para |
---|
272 | IMPLICIT NONE |
---|
273 | CHARACTER(LEN=*) :: Field_name |
---|
274 | INTEGER :: field_size |
---|
275 | REAL :: field(klon,field_size) |
---|
276 | INTEGER,INTENT(IN) :: corners(4) |
---|
277 | INTEGER,INTENT(IN) :: edges(4) |
---|
278 | LOGICAL,OPTIONAL :: found |
---|
279 | |
---|
280 | REAL :: field_glo(klon_glo,field_size) |
---|
281 | LOGICAL :: tmp_found |
---|
282 | INTEGER :: varid |
---|
283 | INTEGER :: ierr |
---|
284 | |
---|
285 | IF (is_master) THEN |
---|
286 | |
---|
287 | ierr=NF90_INQ_VARID(nid_start,Field_name,varid) |
---|
288 | |
---|
289 | IF (ierr==NF90_NOERR) THEN |
---|
290 | CALL body(field_glo) |
---|
291 | tmp_found=.TRUE. |
---|
292 | ELSE |
---|
293 | tmp_found=.FALSE. |
---|
294 | ENDIF |
---|
295 | |
---|
296 | ENDIF |
---|
297 | |
---|
298 | CALL bcast(tmp_found) |
---|
299 | |
---|
300 | IF (tmp_found) THEN |
---|
301 | CALL scatter(field_glo,field) |
---|
302 | ENDIF |
---|
303 | |
---|
304 | IF (PRESENT(found)) THEN |
---|
305 | found=tmp_found |
---|
306 | ELSE |
---|
307 | IF (.NOT. tmp_found) THEN |
---|
308 | PRINT*, 'get_field_rgen: Field <'//field_name//'> not found' |
---|
309 | CALL abort_physic("get_field_rgen","Failed to get field",1) |
---|
310 | ENDIF |
---|
311 | ENDIF |
---|
312 | |
---|
313 | |
---|
314 | CONTAINS |
---|
315 | |
---|
316 | SUBROUTINE body(field_glo) |
---|
317 | REAL :: field_glo(klon_glo*field_size) |
---|
318 | ierr=NF90_GET_VAR(nid_start,varid,field_glo,corners,edges) |
---|
319 | IF (ierr/=NF90_NOERR) THEN |
---|
320 | ! La variable exist dans le fichier mais la lecture a echouee. |
---|
321 | PRINT*, 'get_field_rgen: Failed reading <'//field_name//'>' |
---|
322 | |
---|
323 | ! IF (field_name=='CLWCON' .OR. field_name=='RNEBCON' .OR. field_name=='RATQS') THEN |
---|
324 | ! ! Essaye de lire le variable sur surface uniqument, comme fait avant |
---|
325 | ! field_glo(:)=0. |
---|
326 | ! ierr=NF90_GET_VAR(nid_start,varid,field_glo(1:klon_glo)) |
---|
327 | ! IF (ierr/=NF90_NOERR) THEN |
---|
328 | ! PRINT*, 'phyetat0: Lecture echouee aussi en 2D pour <'//field_name//'>' |
---|
329 | ! CALL abort |
---|
330 | ! ELSE |
---|
331 | ! PRINT*, 'phyetat0: La variable <'//field_name//'> lu sur surface seulement'!, selon ancien format, le reste mis a zero' |
---|
332 | ! END IF |
---|
333 | ! ELSE |
---|
334 | CALL abort_physic("get_field_rgen","Failed to read field",1) |
---|
335 | ! ENDIF |
---|
336 | ENDIF |
---|
337 | |
---|
338 | END SUBROUTINE body |
---|
339 | |
---|
340 | END SUBROUTINE Get_field_rgen |
---|
341 | |
---|
342 | |
---|
343 | SUBROUTINE get_var_r0(var_name,var,found) |
---|
344 | ! Get a scalar from input file |
---|
345 | IMPLICIT NONE |
---|
346 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
347 | REAL,INTENT(INOUT) :: var |
---|
348 | LOGICAL,OPTIONAL,INTENT(OUT) :: found |
---|
349 | |
---|
350 | REAL :: varout(1) |
---|
351 | |
---|
352 | IF (PRESENT(found)) THEN |
---|
353 | CALL Get_var_rgen(var_name,varout,size(varout),found) |
---|
354 | ELSE |
---|
355 | CALL Get_var_rgen(var_name,varout,size(varout)) |
---|
356 | ENDIF |
---|
357 | var=varout(1) |
---|
358 | |
---|
359 | END SUBROUTINE get_var_r0 |
---|
360 | |
---|
361 | SUBROUTINE get_var_r1(var_name,var,found) |
---|
362 | ! Get a vector from input file |
---|
363 | IMPLICIT NONE |
---|
364 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
365 | REAL,INTENT(INOUT) :: var(:) |
---|
366 | LOGICAL,OPTIONAL,INTENT(OUT) :: found |
---|
367 | |
---|
368 | IF (PRESENT(found)) THEN |
---|
369 | CALL Get_var_rgen(var_name,var,size(var),found) |
---|
370 | ELSE |
---|
371 | CALL Get_var_rgen(var_name,var,size(var)) |
---|
372 | ENDIF |
---|
373 | |
---|
374 | END SUBROUTINE get_var_r1 |
---|
375 | |
---|
376 | SUBROUTINE get_var_r2(var_name,var,found) |
---|
377 | ! Get a 2D field from input file |
---|
378 | IMPLICIT NONE |
---|
379 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
380 | REAL,INTENT(OUT) :: var(:,:) |
---|
381 | LOGICAL,OPTIONAL,INTENT(OUT) :: found |
---|
382 | |
---|
383 | IF (PRESENT(found)) THEN |
---|
384 | CALL Get_var_rgen(var_name,var,size(var),found) |
---|
385 | ELSE |
---|
386 | CALL Get_var_rgen(var_name,var,size(var)) |
---|
387 | ENDIF |
---|
388 | |
---|
389 | END SUBROUTINE get_var_r2 |
---|
390 | |
---|
391 | SUBROUTINE get_var_r3(var_name,var,found) |
---|
392 | ! Get a 3D field frominput file |
---|
393 | IMPLICIT NONE |
---|
394 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
395 | REAL,INTENT(INOUT) :: var(:,:,:) |
---|
396 | LOGICAL,OPTIONAL,INTENT(OUT) :: found |
---|
397 | |
---|
398 | IF (PRESENT(found)) THEN |
---|
399 | CALL Get_var_rgen(var_name,var,size(var),found) |
---|
400 | ELSE |
---|
401 | CALL Get_var_rgen(var_name,var,size(var)) |
---|
402 | ENDIF |
---|
403 | |
---|
404 | END SUBROUTINE get_var_r3 |
---|
405 | |
---|
406 | SUBROUTINE Get_var_rgen(var_name,var,var_size,found) |
---|
407 | USE netcdf |
---|
408 | USE dimphy |
---|
409 | USE mod_grid_phy_lmdz |
---|
410 | USE mod_phys_lmdz_para |
---|
411 | IMPLICIT NONE |
---|
412 | CHARACTER(LEN=*) :: var_name |
---|
413 | INTEGER :: var_size |
---|
414 | REAL :: var(var_size) |
---|
415 | LOGICAL,OPTIONAL :: found |
---|
416 | |
---|
417 | LOGICAL :: tmp_found |
---|
418 | INTEGER :: varid |
---|
419 | INTEGER :: ierr |
---|
420 | |
---|
421 | IF (is_mpi_root .AND. is_omp_root) THEN |
---|
422 | |
---|
423 | ierr=NF90_INQ_VARID(nid_start,var_name,varid) |
---|
424 | |
---|
425 | IF (ierr==NF90_NOERR) THEN |
---|
426 | ierr=NF90_GET_VAR(nid_start,varid,var) |
---|
427 | IF (ierr/=NF90_NOERR) THEN |
---|
428 | PRINT*, 'phyetat0: Failed loading <'//trim(var_name)//'>' |
---|
429 | CALL abort_physic("get_var_rgen","Failed to read variable",1) |
---|
430 | ENDIF |
---|
431 | tmp_found=.TRUE. |
---|
432 | ELSE |
---|
433 | tmp_found=.FALSE. |
---|
434 | ENDIF |
---|
435 | |
---|
436 | ENDIF |
---|
437 | |
---|
438 | CALL bcast(tmp_found) |
---|
439 | |
---|
440 | IF (tmp_found) THEN |
---|
441 | CALL bcast(var) |
---|
442 | ENDIF |
---|
443 | |
---|
444 | IF (PRESENT(found)) THEN |
---|
445 | found=tmp_found |
---|
446 | ELSE |
---|
447 | IF (.NOT. tmp_found) THEN |
---|
448 | PRINT*, 'phyetat0: Variable <'//trim(var_name)//'> not found' |
---|
449 | CALL abort_physic("get_var_rgen","Failed to read variable",1) |
---|
450 | ENDIF |
---|
451 | ENDIF |
---|
452 | |
---|
453 | END SUBROUTINE Get_var_rgen |
---|
454 | |
---|
455 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
456 | |
---|
457 | SUBROUTINE open_restartphy(filename) |
---|
458 | USE netcdf, only: NF90_CREATE, NF90_CLOBBER, NF90_64BIT_OFFSET, & |
---|
459 | NF90_NOERR, nf90_strerror, & |
---|
460 | NF90_PUT_ATT, NF90_GLOBAL, NF90_DEF_DIM, & |
---|
461 | NF90_UNLIMITED, NF90_ENDDEF, & |
---|
462 | NF90_WRITE, NF90_OPEN |
---|
463 | USE mod_phys_lmdz_para, only: is_master |
---|
464 | USE mod_grid_phy_lmdz, only: klon_glo |
---|
465 | USE dimphy, only: klev, klevp1 |
---|
466 | USE tracer_mod, only: nqmx |
---|
467 | USE comsoil_h, only: nsoilmx |
---|
468 | IMPLICIT NONE |
---|
469 | CHARACTER(LEN=*),INTENT(IN) :: filename |
---|
470 | INTEGER :: ierr |
---|
471 | LOGICAL,SAVE :: already_created=.false. |
---|
472 | |
---|
473 | IF (is_master) THEN |
---|
474 | if (.not.already_created) then |
---|
475 | ! At the very first call, create the file |
---|
476 | ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), & |
---|
477 | nid_restart) |
---|
478 | IF (ierr/=NF90_NOERR) THEN |
---|
479 | write(*,*)'open_restartphy: problem creating file '//trim(filename) |
---|
480 | write(*,*)trim(nf90_strerror(ierr)) |
---|
481 | CALL abort_physic("open_restartphy","Failed creating file",1) |
---|
482 | ENDIF |
---|
483 | already_created=.true. |
---|
484 | else |
---|
485 | ! Just open the file |
---|
486 | ierr=NF90_OPEN(filename,NF90_WRITE,nid_restart) |
---|
487 | IF (ierr/=NF90_NOERR) THEN |
---|
488 | write(*,*)'open_restartphy: problem opening file '//trim(filename) |
---|
489 | write(*,*)trim(nf90_strerror(ierr)) |
---|
490 | CALL abort_physic("open_restartphy","Failed opening file",1) |
---|
491 | ENDIF |
---|
492 | return |
---|
493 | endif ! of if (.not.already_created) |
---|
494 | |
---|
495 | ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",& |
---|
496 | "Physics start file") |
---|
497 | IF (ierr/=NF90_NOERR) THEN |
---|
498 | write(*,*)'open_restartphy: problem writing title ' |
---|
499 | write(*,*)trim(nf90_strerror(ierr)) |
---|
500 | ENDIF |
---|
501 | |
---|
502 | ierr=NF90_DEF_DIM(nid_restart,"index",length,idim1) |
---|
503 | IF (ierr/=NF90_NOERR) THEN |
---|
504 | write(*,*)'open_restartphy: problem defining index dimension ' |
---|
505 | write(*,*)trim(nf90_strerror(ierr)) |
---|
506 | CALL abort_physic("open_restartphy","Failed defining index",1) |
---|
507 | ENDIF |
---|
508 | |
---|
509 | ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2) |
---|
510 | IF (ierr/=NF90_NOERR) THEN |
---|
511 | write(*,*)'open_restartphy: problem defining physical_points dimension ' |
---|
512 | write(*,*)trim(nf90_strerror(ierr)) |
---|
513 | CALL abort_physic("open_restartphy","Failed defining physical_points",1) |
---|
514 | ENDIF |
---|
515 | |
---|
516 | ierr=NF90_DEF_DIM(nid_restart,"subsurface_layers",nsoilmx,idim3) |
---|
517 | IF (ierr/=NF90_NOERR) THEN |
---|
518 | write(*,*)'open_restartphy: problem defining subsurface_layers dimension ' |
---|
519 | write(*,*)trim(nf90_strerror(ierr)) |
---|
520 | CALL abort_physic("open_restartphy","Failed defining subsurface_layers",1) |
---|
521 | ENDIF |
---|
522 | |
---|
523 | ierr=NF90_DEF_DIM(nid_restart,"nlayer_plus_1",klevp1,idim4) |
---|
524 | IF (ierr/=NF90_NOERR) THEN |
---|
525 | write(*,*)'open_restartphy: problem defining nlayer_plus_1 dimension ' |
---|
526 | write(*,*)trim(nf90_strerror(ierr)) |
---|
527 | CALL abort_physic("open_restartphy","Failed defining nlayer_plus_1",1) |
---|
528 | ENDIF |
---|
529 | |
---|
530 | if (nqmx.ne.0) then |
---|
531 | ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",nqmx,idim5) |
---|
532 | else |
---|
533 | ! pretend nqmx=1 because 0 size implies unlimited dimension for NetCDF |
---|
534 | ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",1,idim5) |
---|
535 | endif |
---|
536 | IF (ierr/=NF90_NOERR) THEN |
---|
537 | write(*,*)'open_restartphy: problem defining number_of_advected_fields dimension ' |
---|
538 | write(*,*)trim(nf90_strerror(ierr)) |
---|
539 | CALL abort_physic("open_restartphy","Failed defining number_of_advected_fields",1) |
---|
540 | ENDIF |
---|
541 | |
---|
542 | ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6) |
---|
543 | IF (ierr/=NF90_NOERR) THEN |
---|
544 | write(*,*)'open_restartphy: problem defining nlayer dimension ' |
---|
545 | write(*,*)trim(nf90_strerror(ierr)) |
---|
546 | CALL abort_physic("open_restartphy","Failed defining nlayer",1) |
---|
547 | ENDIF |
---|
548 | |
---|
549 | ierr=NF90_DEF_DIM(nid_restart,"Time",NF90_UNLIMITED,idim7) |
---|
550 | IF (ierr/=NF90_NOERR) THEN |
---|
551 | write(*,*)'open_restartphy: problem defining Time dimension ' |
---|
552 | write(*,*)trim(nf90_strerror(ierr)) |
---|
553 | CALL abort_physic("open_restartphy","Failed defining Time",1) |
---|
554 | ENDIF |
---|
555 | |
---|
556 | ierr=NF90_ENDDEF(nid_restart) |
---|
557 | IF (ierr/=NF90_NOERR) THEN |
---|
558 | write(*,*)'open_restartphy: problem ending definition mode ' |
---|
559 | write(*,*)trim(nf90_strerror(ierr)) |
---|
560 | CALL abort_physic("open_restartphy","Failed ending definition mode",1) |
---|
561 | ENDIF |
---|
562 | ENDIF |
---|
563 | |
---|
564 | END SUBROUTINE open_restartphy |
---|
565 | |
---|
566 | SUBROUTINE close_restartphy |
---|
567 | USE netcdf, only: NF90_CLOSE |
---|
568 | USE mod_phys_lmdz_para, only: is_master |
---|
569 | IMPLICIT NONE |
---|
570 | INTEGER :: ierr |
---|
571 | |
---|
572 | IF (is_master) THEN |
---|
573 | ierr = NF90_CLOSE (nid_restart) |
---|
574 | ENDIF |
---|
575 | |
---|
576 | END SUBROUTINE close_restartphy |
---|
577 | |
---|
578 | SUBROUTINE put_field_r1(field_name,title,field,time) |
---|
579 | ! For a surface field |
---|
580 | IMPLICIT NONE |
---|
581 | CHARACTER(LEN=*),INTENT(IN) :: field_name |
---|
582 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
583 | REAL,INTENT(IN) :: field(:) |
---|
584 | REAL,OPTIONAL,INTENT(IN) :: time |
---|
585 | |
---|
586 | IF (present(time)) THEN |
---|
587 | ! if timeindex is present, it is a time-dependent variable |
---|
588 | CALL put_field_rgen(field_name,title,field,1,time) |
---|
589 | ELSE |
---|
590 | CALL put_field_rgen(field_name,title,field,1) |
---|
591 | ENDIF |
---|
592 | |
---|
593 | END SUBROUTINE put_field_r1 |
---|
594 | |
---|
595 | SUBROUTINE put_field_r2(field_name,title,field,time) |
---|
596 | ! For a "3D" horizontal-vertical field |
---|
597 | IMPLICIT NONE |
---|
598 | CHARACTER(LEN=*),INTENT(IN) :: field_name |
---|
599 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
600 | REAL,INTENT(IN) :: field(:,:) |
---|
601 | REAL,OPTIONAL,INTENT(IN) :: time |
---|
602 | |
---|
603 | IF (present(time)) THEN |
---|
604 | ! if timeindex is present, it is a time-dependent variable |
---|
605 | CALL put_field_rgen(field_name,title,field,size(field,2),time) |
---|
606 | ELSE |
---|
607 | CALL put_field_rgen(field_name,title,field,size(field,2)) |
---|
608 | ENDIF |
---|
609 | |
---|
610 | END SUBROUTINE put_field_r2 |
---|
611 | |
---|
612 | SUBROUTINE put_field_r3(field_name,title,field,time) |
---|
613 | ! For a "4D" field surf/alt/?? |
---|
614 | IMPLICIT NONE |
---|
615 | CHARACTER(LEN=*),INTENT(IN) :: field_name |
---|
616 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
617 | REAL,INTENT(IN) :: field(:,:,:) |
---|
618 | REAL,OPTIONAL,INTENT(IN) :: time |
---|
619 | |
---|
620 | IF (present(time)) THEN |
---|
621 | ! if timeindex is present, it is a time-dependent variable |
---|
622 | CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3),& |
---|
623 | time) |
---|
624 | ELSE |
---|
625 | CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3)) |
---|
626 | ENDIF |
---|
627 | |
---|
628 | END SUBROUTINE put_field_r3 |
---|
629 | |
---|
630 | SUBROUTINE put_field_rgen(field_name,title,field,field_size,time) |
---|
631 | USE netcdf |
---|
632 | USE dimphy |
---|
633 | USE comsoil_h, only: nsoilmx |
---|
634 | USE mod_grid_phy_lmdz |
---|
635 | USE mod_phys_lmdz_para |
---|
636 | IMPLICIT NONE |
---|
637 | CHARACTER(LEN=*),INTENT(IN) :: field_name |
---|
638 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
639 | INTEGER,INTENT(IN) :: field_size |
---|
640 | REAL,INTENT(IN) :: field(klon,field_size) |
---|
641 | REAL,OPTIONAL,INTENT(IN) :: time |
---|
642 | |
---|
643 | REAL :: field_glo(klon_glo,field_size) |
---|
644 | INTEGER :: ierr |
---|
645 | INTEGER :: nvarid |
---|
646 | INTEGER :: idim |
---|
647 | |
---|
648 | CALL gather(field,field_glo) |
---|
649 | |
---|
650 | IF (is_master) THEN |
---|
651 | |
---|
652 | IF (field_size==1) THEN |
---|
653 | ! input is a 1D "surface field" array |
---|
654 | if (.not.present(time)) then ! for a time-independent field |
---|
655 | ierr=NF90_REDEF(nid_restart) |
---|
656 | #ifdef NC_DOUBLE |
---|
657 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
658 | (/idim2/),nvarid) |
---|
659 | #else |
---|
660 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
661 | (/idim2/),nvarid) |
---|
662 | #endif |
---|
663 | if (ierr.ne.NF90_NOERR) then |
---|
664 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
665 | write(*,*)trim(nf90_strerror(ierr)) |
---|
666 | endif |
---|
667 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
668 | ierr=NF90_ENDDEF(nid_restart) |
---|
669 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo) |
---|
670 | else |
---|
671 | ! check if the variable has already been defined: |
---|
672 | ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid) |
---|
673 | if (ierr/=NF90_NOERR) then ! variable not found, define it |
---|
674 | ierr=NF90_REDEF(nid_restart) |
---|
675 | #ifdef NC_DOUBLE |
---|
676 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
677 | (/idim2,idim7/),nvarid) |
---|
678 | #else |
---|
679 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
680 | (/idim2,idim7/),nvarid) |
---|
681 | #endif |
---|
682 | if (ierr.ne.NF90_NOERR) then |
---|
683 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
684 | write(*,*)trim(nf90_strerror(ierr)) |
---|
685 | endif |
---|
686 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
687 | ierr=NF90_ENDDEF(nid_restart) |
---|
688 | endif |
---|
689 | ! Write the variable |
---|
690 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,& |
---|
691 | start=(/1,timeindex/)) |
---|
692 | endif ! of if (.not.present(timeindex)) |
---|
693 | |
---|
694 | ELSE IF (field_size==klev) THEN |
---|
695 | ! input is a 2D "atmospheric field" array |
---|
696 | if (.not.present(time)) then ! for a time-independent field |
---|
697 | ierr=NF90_REDEF(nid_restart) |
---|
698 | #ifdef NC_DOUBLE |
---|
699 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
700 | (/idim2,idim6/),nvarid) |
---|
701 | #else |
---|
702 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
703 | (/idim2,idim6/),nvarid) |
---|
704 | #endif |
---|
705 | if (ierr.ne.NF90_NOERR) then |
---|
706 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
707 | write(*,*)trim(nf90_strerror(ierr)) |
---|
708 | endif |
---|
709 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
710 | ierr=NF90_ENDDEF(nid_restart) |
---|
711 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo) |
---|
712 | else |
---|
713 | ! check if the variable has already been defined: |
---|
714 | ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid) |
---|
715 | if (ierr/=NF90_NOERR) then ! variable not found, define it |
---|
716 | ierr=NF90_REDEF(nid_restart) |
---|
717 | #ifdef NC_DOUBLE |
---|
718 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
719 | (/idim2,idim6,idim7/),nvarid) |
---|
720 | #else |
---|
721 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
722 | (/idim2,idim6,idim7/),nvarid) |
---|
723 | #endif |
---|
724 | if (ierr.ne.NF90_NOERR) then |
---|
725 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
726 | write(*,*)trim(nf90_strerror(ierr)) |
---|
727 | endif |
---|
728 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
729 | ierr=NF90_ENDDEF(nid_restart) |
---|
730 | endif |
---|
731 | ! Write the variable |
---|
732 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,& |
---|
733 | start=(/1,1,timeindex/)) |
---|
734 | endif ! of if (.not.present(time)) |
---|
735 | |
---|
736 | ELSE IF (field_size==klevp1) THEN |
---|
737 | ! input is a 2D "interlayer atmospheric field" array |
---|
738 | if (.not.present(time)) then ! for a time-independent field |
---|
739 | ierr=NF90_REDEF(nid_restart) |
---|
740 | #ifdef NC_DOUBLE |
---|
741 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
742 | (/idim2,idim4/),nvarid) |
---|
743 | #else |
---|
744 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
745 | (/idim2,idim4/),nvarid) |
---|
746 | #endif |
---|
747 | if (ierr.ne.NF90_NOERR) then |
---|
748 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
749 | write(*,*)trim(nf90_strerror(ierr)) |
---|
750 | endif |
---|
751 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
752 | ierr = NF90_ENDDEF(nid_restart) |
---|
753 | ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo) |
---|
754 | else |
---|
755 | ! check if the variable has already been defined: |
---|
756 | ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid) |
---|
757 | if (ierr/=NF90_NOERR) then ! variable not found, define it |
---|
758 | ierr=NF90_REDEF(nid_restart) |
---|
759 | #ifdef NC_DOUBLE |
---|
760 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
761 | (/idim2,idim4,idim7/),nvarid) |
---|
762 | #else |
---|
763 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
764 | (/idim2,idim4,idim7/),nvarid) |
---|
765 | #endif |
---|
766 | if (ierr.ne.NF90_NOERR) then |
---|
767 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
768 | write(*,*)trim(nf90_strerror(ierr)) |
---|
769 | endif |
---|
770 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
771 | ierr=NF90_ENDDEF(nid_restart) |
---|
772 | endif |
---|
773 | ! Write the variable |
---|
774 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,& |
---|
775 | start=(/1,1,timeindex/)) |
---|
776 | endif ! of if (.not.present(timeindex)) |
---|
777 | |
---|
778 | ELSE IF (field_size==nsoilmx) THEN |
---|
779 | ! input is a 2D "subsurface field" array |
---|
780 | if (.not.present(time)) then ! for a time-independent field |
---|
781 | ierr = NF90_REDEF(nid_restart) |
---|
782 | #ifdef NC_DOUBLE |
---|
783 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
784 | (/idim2,idim3/),nvarid) |
---|
785 | #else |
---|
786 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
787 | (/idim2,idim3/),nvarid) |
---|
788 | #endif |
---|
789 | if (ierr.ne.NF90_NOERR) then |
---|
790 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
791 | write(*,*)trim(nf90_strerror(ierr)) |
---|
792 | endif |
---|
793 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
794 | ierr = NF90_ENDDEF(nid_restart) |
---|
795 | ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo) |
---|
796 | else |
---|
797 | ! check if the variable has already been defined: |
---|
798 | ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid) |
---|
799 | if (ierr/=NF90_NOERR) then ! variable not found, define it |
---|
800 | ierr=NF90_REDEF(nid_restart) |
---|
801 | #ifdef NC_DOUBLE |
---|
802 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,& |
---|
803 | (/idim2,idim3,idim7/),nvarid) |
---|
804 | #else |
---|
805 | ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,& |
---|
806 | (/idim2,idim3,idim7/),nvarid) |
---|
807 | #endif |
---|
808 | if (ierr.ne.NF90_NOERR) then |
---|
809 | write(*,*)"put_field_rgen error: failed to define "//trim(field_name) |
---|
810 | write(*,*)trim(nf90_strerror(ierr)) |
---|
811 | endif |
---|
812 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
813 | ierr=NF90_ENDDEF(nid_restart) |
---|
814 | endif |
---|
815 | ! Write the variable |
---|
816 | ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,& |
---|
817 | start=(/1,1,timeindex/)) |
---|
818 | |
---|
819 | endif ! of if (.not.present(time)) |
---|
820 | |
---|
821 | ELSE |
---|
822 | PRINT *, "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name) |
---|
823 | write(*,*) " field_size =",field_size |
---|
824 | CALL abort_physic("put_field_rgen","wrong field dimension",1) |
---|
825 | ENDIF |
---|
826 | |
---|
827 | ! Check the writting of field to file went OK |
---|
828 | if (ierr.ne.NF90_NOERR) then |
---|
829 | write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name) |
---|
830 | write(*,*)trim(nf90_strerror(ierr)) |
---|
831 | CALL abort_physic("put_field_rgen","Failed writing field",1) |
---|
832 | endif |
---|
833 | |
---|
834 | ENDIF ! of IF (is_master) |
---|
835 | |
---|
836 | END SUBROUTINE put_field_rgen |
---|
837 | |
---|
838 | SUBROUTINE put_var_r0(var_name,title,var) |
---|
839 | ! Put a scalar in file |
---|
840 | IMPLICIT NONE |
---|
841 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
842 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
843 | REAL,INTENT(IN) :: var |
---|
844 | REAL :: varin(1) |
---|
845 | |
---|
846 | varin(1)=var |
---|
847 | |
---|
848 | CALL put_var_rgen(var_name,title,varin,size(varin)) |
---|
849 | |
---|
850 | END SUBROUTINE put_var_r0 |
---|
851 | |
---|
852 | |
---|
853 | SUBROUTINE put_var_r1(var_name,title,var) |
---|
854 | ! Put a vector in file |
---|
855 | IMPLICIT NONE |
---|
856 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
857 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
858 | REAL,INTENT(IN) :: var(:) |
---|
859 | |
---|
860 | CALL put_var_rgen(var_name,title,var,size(var)) |
---|
861 | |
---|
862 | END SUBROUTINE put_var_r1 |
---|
863 | |
---|
864 | SUBROUTINE put_var_r2(var_name,title,var) |
---|
865 | ! Put a 2D field in file |
---|
866 | IMPLICIT NONE |
---|
867 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
868 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
869 | REAL,INTENT(IN) :: var(:,:) |
---|
870 | |
---|
871 | CALL put_var_rgen(var_name,title,var,size(var)) |
---|
872 | |
---|
873 | END SUBROUTINE put_var_r2 |
---|
874 | |
---|
875 | SUBROUTINE put_var_r3(var_name,title,var) |
---|
876 | ! Put a 3D field in file |
---|
877 | IMPLICIT NONE |
---|
878 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
879 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
880 | REAL,INTENT(IN) :: var(:,:,:) |
---|
881 | |
---|
882 | CALL put_var_rgen(var_name,title,var,size(var)) |
---|
883 | |
---|
884 | END SUBROUTINE put_var_r3 |
---|
885 | |
---|
886 | SUBROUTINE put_var_rgen(var_name,title,var,var_size) |
---|
887 | USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, & |
---|
888 | NF90_FLOAT, NF90_DOUBLE, & |
---|
889 | NF90_PUT_ATT, NF90_NOERR, nf90_strerror, & |
---|
890 | nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID |
---|
891 | USE comsoil_h, only: nsoilmx |
---|
892 | USE mod_phys_lmdz_para, only: is_master |
---|
893 | IMPLICIT NONE |
---|
894 | CHARACTER(LEN=*),INTENT(IN) :: var_name |
---|
895 | CHARACTER(LEN=*),INTENT(IN) :: title |
---|
896 | INTEGER,INTENT(IN) :: var_size |
---|
897 | REAL,INTENT(IN) :: var(var_size) |
---|
898 | |
---|
899 | INTEGER :: ierr |
---|
900 | INTEGER :: nvarid |
---|
901 | INTEGER :: idim1d |
---|
902 | logical,save :: firsttime=.true. |
---|
903 | |
---|
904 | IF (is_master) THEN |
---|
905 | |
---|
906 | IF (var_name=="Time") THEN |
---|
907 | ! Very specific case of "Time" variable |
---|
908 | if (firsttime) then |
---|
909 | ! Create the "Time variable" |
---|
910 | ierr=NF90_REDEF(nid_restart) |
---|
911 | #ifdef NC_DOUBLE |
---|
912 | ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,& |
---|
913 | (/idim7/),nvarid) |
---|
914 | #else |
---|
915 | ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,& |
---|
916 | (/idim7/),nvarid) |
---|
917 | #endif |
---|
918 | IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
919 | ierr=NF90_ENDDEF(nid_restart) |
---|
920 | |
---|
921 | firsttime=.false. |
---|
922 | endif |
---|
923 | ! Append "time" value |
---|
924 | ! get current length of "Time" dimension |
---|
925 | ierr=nf90_inq_dimid(nid_restart,var_name,idim1d) |
---|
926 | ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex) |
---|
927 | timeindex=timeindex+1 |
---|
928 | ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid) |
---|
929 | ierr=NF90_PUT_VAR(nid_restart,nvarid,var,& |
---|
930 | start=(/timeindex/)) |
---|
931 | IF (ierr/=NF90_NOERR) THEN |
---|
932 | write(*,*)'put_var_rgen: problem writing Time' |
---|
933 | write(*,*)trim(nf90_strerror(ierr)) |
---|
934 | CALL abort_physic("get_var_rgen","Failed to write Time",1) |
---|
935 | ENDIF |
---|
936 | return ! nothing left to do |
---|
937 | ELSEIF (var_size==length) THEN |
---|
938 | ! We know it is a "controle" kind of 1D array |
---|
939 | idim1d=idim1 |
---|
940 | ELSEIF (var_size==nsoilmx) THEN |
---|
941 | ! We know it is an "mlayer" kind of 1D array |
---|
942 | idim1d=idim3 |
---|
943 | ELSE |
---|
944 | PRINT *, "put_var_rgen error : wrong dimension" |
---|
945 | write(*,*) " var_size =",var_size |
---|
946 | CALL abort_physic("get_var_rgen","Wrong variable dimension",1) |
---|
947 | |
---|
948 | ENDIF ! of IF (var_size==length) THEN |
---|
949 | |
---|
950 | ! Swich to NetCDF define mode |
---|
951 | ierr=NF90_REDEF (nid_restart) |
---|
952 | ! Define the variable |
---|
953 | #ifdef NC_DOUBLE |
---|
954 | ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid) |
---|
955 | #else |
---|
956 | ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid) |
---|
957 | #endif |
---|
958 | ! Add a "title" attribute |
---|
959 | IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title) |
---|
960 | ! Swich out of define mode |
---|
961 | ierr=NF90_ENDDEF(nid_restart) |
---|
962 | ! Write variable to file |
---|
963 | ierr=NF90_PUT_VAR(nid_restart,nvarid,var) |
---|
964 | IF (ierr/=NF90_NOERR) THEN |
---|
965 | write(*,*)'put_var_rgen: problem writing '//trim(var_name) |
---|
966 | write(*,*)trim(nf90_strerror(ierr)) |
---|
967 | CALL abort_physic("get_var_rgen","Failed writing variable",1) |
---|
968 | ENDIF |
---|
969 | ENDIF ! of IF (is_master) |
---|
970 | |
---|
971 | END SUBROUTINE put_var_rgen |
---|
972 | |
---|
973 | END MODULE iostart |
---|