source: trunk/LMDZ.GENERIC/libf/phystd/iostart.F90 @ 3809

Last change on this file since 3809 was 3747, checked in by afalco, 2 months ago

Generic/Pluto?: create restartfi if it does not exist when opening it.
AF

File size: 42.7 KB
Line 
1MODULE iostart
2
3    IMPLICIT NONE
4PRIVATE
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!$OMP THREADPRIVATE(nid_start)
8
9    ! restartfi.nc file dimension identifiers: (see open_restartphy())
10    INTEGER,SAVE :: idim1 ! "index" dimension
11    INTEGER,SAVE :: idim2 ! "physical_points" dimension
12    INTEGER,SAVE :: idim3 ! "subsurface_layers" dimension
13    INTEGER,SAVE :: idim4 ! "nlayer_plus_1" dimension
14    INTEGER,SAVE :: idim5 ! "number_of_advected_fields" dimension
15    INTEGER,SAVE :: idim6 ! "nlayer" dimension
16    INTEGER,SAVE :: idim7 ! "Time" dimension
17    INTEGER,SAVE :: idim8 ! "ocean_layers" dimension
18!$OMP THREADPRIVATE(idim1,idim2,idim3,idim4,idim5,idim6,idim7,idim8)
19    INTEGER,SAVE :: idim10 ! "descriptor" dimension
20    INTEGER,SAVE :: idim11 ! "description_size" dimension
21    INTEGER,SAVE :: timeindex ! current time index (for time-dependent fields)
22!$OMP THREADPRIVATE(idim10,idim11,timeindex)
23    INTEGER,PARAMETER :: length=100 ! size of tab_cntrl array
24    INTEGER,PARAMETER :: ldscrpt = 35 ! size of dscrpt_tab_cntrl array
25    INTEGER,PARAMETER :: ndscrpt = 50 ! size of characters in dscrpt_tab_cntrl array
26
27    INTERFACE get_field
28      MODULE PROCEDURE Get_field_r1,Get_field_r2,Get_field_r3
29    END INTERFACE get_field
30
31    INTERFACE get_var
32      MODULE PROCEDURE get_var_r0,Get_var_r1,Get_var_r2,Get_var_r3
33    END INTERFACE get_var
34
35    INTERFACE put_field
36      MODULE PROCEDURE put_field_r1,put_field_r2,put_field_r3
37    END INTERFACE put_field
38
39    INTERFACE put_var
40      MODULE PROCEDURE put_var_r0,put_var_r1,put_var_r2,put_var_r3,put_var_c1
41    END INTERFACE put_var
42
43    PUBLIC nid_start, length, ldscrpt, ndscrpt
44    PUBLIC get_field,get_var,put_field,put_var
45    PUBLIC inquire_dimension, inquire_dimension_length
46    PUBLIC inquire_field, inquire_field_ndims
47    PUBLIC open_startphy,close_startphy,create_restartphy,open_restartphy,close_restartphy
48
49CONTAINS
50
51  SUBROUTINE open_startphy(filename,nid_start)
52  USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_NOWRITE, nf90_strerror
53  USE mod_phys_lmdz_para, only: is_master, bcast
54  IMPLICIT NONE
55    CHARACTER(LEN=*)      :: filename
56    INTEGER,INTENT(INOUT) :: nid_start
57    INTEGER               :: ierr
58
59    IF (is_master) THEN
60      ierr = NF90_OPEN (filename, NF90_NOWRITE, nid_start)
61      IF (ierr.NE.NF90_NOERR) THEN
62        write(*,*)'open_startphy: problem opening file '//trim(filename)
63        write(*,*)trim(nf90_strerror(ierr))
64        CALL abort_physic("open_startphy","Cannot open file",1)
65      ENDIF
66    ENDIF
67
68    CALL bcast(nid_start) ! tell all procs about nid_start
69
70  END SUBROUTINE open_startphy
71
72  SUBROUTINE close_startphy(nid_start)
73  USE netcdf, only: NF90_CLOSE
74  USE mod_phys_lmdz_para, only: is_master
75  IMPLICIT NONE
76    INTEGER,INTENT(IN) :: nid_start
77    INTEGER            :: ierr
78
79    IF (is_master) THEN
80        ierr = NF90_CLOSE (nid_start)
81    ENDIF
82
83  END SUBROUTINE close_startphy
84
85
86  FUNCTION inquire_field(nid_start,Field_name)
87  ! check if a given field is present in the input file
88  USE netcdf, only: NF90_INQ_VARID, NF90_NOERR
89  USE mod_phys_lmdz_para, only: is_master, bcast
90  IMPLICIT NONE
91    INTEGER,INTENT(IN)          :: nid_start
92    CHARACTER(LEN=*),INTENT(IN) :: Field_name
93    LOGICAL :: inquire_field
94    INTEGER :: varid
95    INTEGER :: ierr
96
97    IF (is_master) THEN
98      ierr=NF90_INQ_VARID(nid_start,Field_name,varid)
99      IF (ierr==NF90_NOERR) THEN
100        Inquire_field=.TRUE.
101      ELSE
102        Inquire_field=.FALSE.
103      ENDIF
104    ENDIF
105
106    CALL bcast(inquire_field)
107
108  END FUNCTION inquire_field
109
110
111  FUNCTION inquire_field_ndims(nid_start,Field_name)
112  ! give the number of dimensions of "Field_name" stored in the input file
113  USE netcdf, only: nf90_inq_varid, nf90_inquire_variable, &
114                    NF90_NOERR, nf90_strerror
115  USE mod_phys_lmdz_para, only: is_master, bcast
116  IMPLICIT NONE
117    INTEGER,INTENT(IN)          :: nid_start
118    CHARACTER(LEN=*),INTENT(IN) :: Field_name
119    INTEGER :: inquire_field_ndims
120    INTEGER :: varid
121    INTEGER :: ierr
122
123    IF (is_master) THEN
124      ierr=nf90_inq_varid(nid_start,Field_name,varid)
125      ierr=nf90_inquire_variable(nid_start,varid,&
126                                  ndims=inquire_field_ndims)
127      IF (ierr.NE.NF90_NOERR) THEN
128        write(*,*)'inquire_field_ndims: problem obtaining ndims of '&
129                  //trim(field_name)
130        write(*,*)trim(nf90_strerror(ierr))
131        CALL abort_physic("inquire_field_ndims","Failed to get ndims",1)
132      ENDIF
133    ENDIF
134
135    CALL bcast(inquire_field_ndims)
136
137  END FUNCTION inquire_field_ndims
138
139
140  FUNCTION inquire_dimension(nid_start,Field_name)
141  ! check if a given dimension is present in the input file
142  USE netcdf, only: nf90_inq_dimid, NF90_NOERR
143  USE mod_phys_lmdz_para, only: is_master, bcast
144  IMPLICIT NONE
145    INTEGER,INTENT(IN)          :: nid_start
146    CHARACTER(LEN=*),INTENT(IN) :: Field_name
147    LOGICAL :: inquire_dimension
148    INTEGER :: varid
149    INTEGER :: ierr
150
151    IF (is_master) THEN
152      ierr=NF90_INQ_DIMID(nid_start,Field_name,varid)
153      IF (ierr==NF90_NOERR) THEN
154        Inquire_dimension=.TRUE.
155      ELSE
156        Inquire_dimension=.FALSE.
157      ENDIF
158    ENDIF
159
160    CALL bcast(inquire_dimension)
161
162  END FUNCTION inquire_dimension
163
164  FUNCTION inquire_dimension_length(nid_start,Field_name)
165  ! give the length of the "Field_name" dimension stored in the input file
166  USE netcdf, only: nf90_inquire_dimension, nf90_inq_dimid, &
167                    NF90_NOERR, nf90_strerror
168  USE mod_phys_lmdz_para, only: is_master, bcast
169  IMPLICIT NONE
170    INTEGER,INTENT(IN)          :: nid_start
171    CHARACTER(LEN=*),INTENT(IN) :: Field_name
172    INTEGER :: inquire_dimension_length
173    INTEGER :: varid
174    INTEGER :: ierr
175
176    IF (is_master) THEN
177      ierr=nf90_inq_dimid(nid_start,Field_name,varid)
178      ierr=nf90_inquire_dimension(nid_start,varid,&
179                                  len=inquire_dimension_length)
180      IF (ierr.NE.NF90_NOERR) THEN
181        write(*,*)'inquire_field_length: problem obtaining length of '&
182                  //trim(field_name)
183        write(*,*)trim(nf90_strerror(ierr))
184        CALL abort_physic("inquire_field_ndims","Failed to get length",1)
185      ENDIF
186    ENDIF
187
188    CALL bcast(inquire_dimension_length)
189
190  END FUNCTION inquire_dimension_length
191
192
193
194  SUBROUTINE Get_Field_r1(nid_start,field_name,field,found,timeindex)
195  ! For a surface field
196  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
197  IMPLICIT NONE
198    INTEGER,INTENT(IN)             :: nid_start
199    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
200    REAL,INTENT(INOUT)               :: Field(:)
201    LOGICAL,INTENT(OUT),OPTIONAL   :: found
202    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
203
204    integer :: edges(4), corners(4)
205
206    edges(1)=klon_glo
207    edges(2:4)=1
208    corners(1:4)=1
209    if (PRESENT(timeindex)) then
210      corners(2)=timeindex
211    endif
212
213    IF (PRESENT(found)) THEN
214      CALL Get_field_rgen(nid_start,field_name,field,1,corners,edges,found)
215    ELSE
216      CALL Get_field_rgen(nid_start,field_name,field,1,corners,edges)
217    ENDIF
218
219  END SUBROUTINE Get_Field_r1
220
221  SUBROUTINE Get_Field_r2(nid_start,field_name,field,found,timeindex)
222  ! For a "3D" horizontal-vertical field
223  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
224  IMPLICIT NONE
225    INTEGER,INTENT(IN)             :: nid_start
226    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
227    REAL,INTENT(INOUT)               :: Field(:,:)
228    LOGICAL,INTENT(OUT),OPTIONAL   :: found
229    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
230
231    integer :: edges(4), corners(4)
232
233    edges(1)=klon_glo
234    edges(2)=size(field,2)
235    edges(3:4)=1
236    corners(1:4)=1
237    if (PRESENT(timeindex)) then
238      corners(3)=timeindex
239    endif
240
241    IF (PRESENT(found)) THEN
242      CALL Get_field_rgen(nid_start,field_name,field,size(field,2),&
243                          corners,edges,found)
244    ELSE
245      CALL Get_field_rgen(nid_start,field_name,field,size(field,2),&
246                          corners,edges)
247    ENDIF
248
249
250  END SUBROUTINE Get_Field_r2
251
252  SUBROUTINE Get_Field_r3(nid_start,field_name,field,found,timeindex)
253  ! for a "4D" field surf/alt/??
254  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
255  IMPLICIT NONE
256    INTEGER,INTENT(IN)             :: nid_start
257    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
258    REAL,INTENT(INOUT)               :: Field(:,:,:)
259    LOGICAL,INTENT(OUT),OPTIONAL   :: found
260    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
261
262    integer :: edges(4), corners(4)
263
264    edges(1)=klon_glo
265    edges(2)=size(field,2)
266    edges(3)=size(field,3)
267    edges(4)=1
268    corners(1:4)=1
269    if (PRESENT(timeindex)) then
270      corners(4)=timeindex
271    endif
272
273    IF (PRESENT(found)) THEN
274      CALL Get_field_rgen(nid_start,field_name,field,size(field,2)*size(field,3),&
275                          corners,edges,found)
276    ELSE
277      CALL Get_field_rgen(nid_start,field_name,field,size(field,2)*size(field,3),&
278                          corners,edges)
279    ENDIF
280
281  END SUBROUTINE Get_Field_r3
282
283  SUBROUTINE Get_field_rgen(nid_start,field_name,field,field_size, &
284                            corners,edges,found)
285  USE netcdf, ONLY: NF90_INQ_VARID, NF90_GET_VAR, NF90_NOERR
286  USE dimphy, ONLY: klon ! number of columns on local grid
287  USE geometry_mod, ONLY: ind_cell_glo
288  USE mod_grid_phy_lmdz, ONLY: klon_glo ! number of columns on global grid
289  USE mod_phys_lmdz_para, ONLY: is_master, bcast, scatter, gather
290  IMPLICIT NONE
291    INTEGER,INTENT(IN) :: nid_start
292    CHARACTER(LEN=*),INTENT(IN) :: Field_name
293    INTEGER,INTENT(IN) :: field_size
294    REAL,INTENT(OUT) :: field(klon,field_size)
295    INTEGER,INTENT(IN) :: corners(4)
296    INTEGER,INTENT(IN) :: edges(4)
297    LOGICAL,OPTIONAL,INTENT(OUT) :: found
298
299    REAL    :: field_glo(klon_glo,field_size) ! field on global grid
300    REAL :: field_glo_tmp(klon_glo,field_size)
301    INTEGER :: ind_cell_glo_glo(klon_glo) ! cell indexes on global grid
302    LOGICAL :: tmp_found
303    INTEGER :: varid
304    INTEGER :: ierr, i
305
306    ! gather columns indexes on global grid
307    CALL gather(ind_cell_glo,ind_cell_glo_glo)
308
309    IF (is_master) THEN
310
311      ierr=NF90_INQ_VARID(nid_start,Field_name,varid)
312
313      IF (ierr==NF90_NOERR) THEN
314        CALL body(field_glo_tmp,nid_start)
315        tmp_found=.TRUE.
316      ELSE
317        tmp_found=.FALSE.
318      ENDIF
319
320    ENDIF ! of IF (is_master)
321
322    CALL bcast(tmp_found)
323
324    IF (tmp_found) THEN
325      IF (is_master) THEN
326        ! reorder columns according to ind_cell_glo(:) indexes
327        DO i=1,klon_glo
328          field_glo(i,:)=field_glo_tmp(ind_cell_glo_glo(i),:)
329        ENDDO
330      ENDIF
331      CALL scatter(field_glo,field)
332    ENDIF
333
334    IF (PRESENT(found)) THEN
335      found=tmp_found
336    ELSE
337      IF (.NOT. tmp_found) THEN
338        PRINT*, 'get_field_rgen: Field <'//field_name//'> not found'
339        CALL abort
340      ENDIF
341    ENDIF
342
343
344    CONTAINS
345
346     SUBROUTINE body(field_glo,nid_start)
347       REAL :: field_glo(klon_glo*field_size)
348       INTEGER,INTENT(IN) :: nid_start
349         ierr=NF90_GET_VAR(nid_start,varid,field_glo,corners,edges)
350         IF (ierr/=NF90_NOERR) THEN
351           ! La variable exist dans le fichier mais la lecture a echouee.
352           PRINT*, 'get_field_rgen: Failed reading <'//field_name//'>'
353
354           CALL abort_physic("get_field_rgen","Failed to read field",1)
355         ENDIF
356
357     END SUBROUTINE body
358
359  END SUBROUTINE Get_field_rgen
360
361
362  SUBROUTINE get_var_r0(nid_start,var_name,var,found)
363  ! Get a scalar from input file
364  IMPLICIT NONE
365    INTEGER,INTENT(IN)           :: nid_start
366    CHARACTER(LEN=*),INTENT(IN)  :: var_name
367    REAL,INTENT(INOUT)             :: var
368    LOGICAL,OPTIONAL,INTENT(OUT) :: found
369
370    REAL                         :: varout(1)
371
372    IF (PRESENT(found)) THEN
373      CALL Get_var_rgen(nid_start,var_name,varout,size(varout),found)
374    ELSE
375      CALL Get_var_rgen(nid_start,var_name,varout,size(varout))
376    ENDIF
377    var=varout(1)
378
379  END SUBROUTINE get_var_r0
380
381  SUBROUTINE get_var_r1(nid_start,var_name,var,found)
382  ! Get a vector from input file
383  IMPLICIT NONE
384    CHARACTER(LEN=*),INTENT(IN)  :: var_name
385    REAL,INTENT(INOUT)             :: var(:)
386    LOGICAL,OPTIONAL,INTENT(OUT) :: found
387    INTEGER,INTENT(IN)           :: nid_start
388
389    IF (PRESENT(found)) THEN
390      CALL Get_var_rgen(nid_start,var_name,var,size(var),found)
391    ELSE
392      CALL Get_var_rgen(nid_start,var_name,var,size(var))
393    ENDIF
394
395  END SUBROUTINE get_var_r1
396
397  SUBROUTINE get_var_r2(nid_start,var_name,var,found)
398  ! Get a 2D field from input file
399  IMPLICIT NONE
400    CHARACTER(LEN=*),INTENT(IN)  :: var_name
401    REAL,INTENT(OUT)             :: var(:,:)
402    LOGICAL,OPTIONAL,INTENT(OUT) :: found
403    INTEGER,INTENT(IN)           :: nid_start
404
405    IF (PRESENT(found)) THEN
406      CALL Get_var_rgen(nid_start,var_name,var,size(var),found)
407    ELSE
408      CALL Get_var_rgen(nid_start,var_name,var,size(var))
409    ENDIF
410
411  END SUBROUTINE get_var_r2
412
413  SUBROUTINE get_var_r3(nid_start,var_name,var,found)
414  ! Get a 3D field frominput file
415  IMPLICIT NONE
416    CHARACTER(LEN=*),INTENT(IN)  :: var_name
417    REAL,INTENT(INOUT)             :: var(:,:,:)
418    LOGICAL,OPTIONAL,INTENT(OUT) :: found
419    INTEGER,INTENT(IN)           :: nid_start
420
421    IF (PRESENT(found)) THEN
422      CALL Get_var_rgen(nid_start,var_name,var,size(var),found)
423    ELSE
424      CALL Get_var_rgen(nid_start,var_name,var,size(var))
425    ENDIF
426
427  END SUBROUTINE get_var_r3
428
429  SUBROUTINE Get_var_rgen(nid_start,var_name,var,var_size,found)
430  USE netcdf, ONLY: NF90_INQ_VARID, NF90_GET_VAR, NF90_NOERR
431  USE mod_phys_lmdz_para, ONLY: is_master, bcast
432  IMPLICIT NONE
433    INTEGER,INTENT(IN) :: nid_start
434    CHARACTER(LEN=*),INTENT(IN) :: var_name
435    INTEGER,INTENT(IN) :: var_size
436    REAL,INTENT(OUT) :: var(var_size)
437    LOGICAL,OPTIONAL,INTENT(OUT) :: found
438
439    LOGICAL :: tmp_found
440    INTEGER :: varid
441    INTEGER :: ierr
442
443    IF (is_master) THEN
444
445      ierr=NF90_INQ_VARID(nid_start,var_name,varid)
446
447      IF (ierr==NF90_NOERR) THEN
448        ierr=NF90_GET_VAR(nid_start,varid,var)
449        IF (ierr/=NF90_NOERR) THEN
450          PRINT*, 'phyetat0: Failed loading <'//trim(var_name)//'>'
451          CALL abort_physic("get_var_rgen","Failed to read variable",1)
452        ENDIF
453        tmp_found=.TRUE.
454      ELSE
455        tmp_found=.FALSE.
456      ENDIF
457
458    ENDIF
459
460    CALL bcast(tmp_found)
461
462    IF (tmp_found) THEN
463      CALL bcast(var)
464    ENDIF
465
466    IF (PRESENT(found)) THEN
467      found=tmp_found
468    ELSE
469      IF (.NOT. tmp_found) THEN
470        PRINT*, 'phyetat0: Variable <'//trim(var_name)//'> not found'
471        CALL abort_physic("get_var_rgen","Failed to read variable",1)
472      ENDIF
473    ENDIF
474
475  END SUBROUTINE Get_var_rgen
476
477!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
478
479  SUBROUTINE create_restartphy(filename,nid_restart)
480  USE netcdf, only: NF90_CREATE, NF90_CLOBBER, NF90_64BIT_OFFSET, &
481                    NF90_NOERR, nf90_strerror, &
482                    NF90_PUT_ATT, NF90_GLOBAL, NF90_DEF_DIM, &
483                    NF90_UNLIMITED, NF90_ENDDEF
484  USE mod_phys_lmdz_para, only: is_master
485  USE mod_grid_phy_lmdz, only: klon_glo
486  USE dimphy, only: klev, klevp1
487  USE tracer_h, only: nqtot
488  USE comsoil_h, only: nsoilmx
489!  USE slab_ice_h, only: noceanmx
490  USE ocean_slab_mod, ONLY: nslay
491
492  IMPLICIT NONE
493    CHARACTER(LEN=*),INTENT(IN) :: filename
494    INTEGER,INTENT(INOUT)       :: nid_restart
495    INTEGER                     :: ierr
496
497    IF (is_master) THEN
498
499        ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), &
500                          nid_restart)
501        IF (ierr/=NF90_NOERR) THEN
502          write(*,*)'create_restartphy: problem creating file '//trim(filename)
503          write(*,*)trim(nf90_strerror(ierr))
504          CALL abort_physic("create_restartphy","Failed creating file",1)
505        ENDIF
506
507      ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",&
508                        "Physics start file")
509      IF (ierr/=NF90_NOERR) THEN
510        write(*,*)'create_restartphy: problem writing title '
511        write(*,*)trim(nf90_strerror(ierr))
512      ENDIF
513
514      ierr=NF90_DEF_DIM(nid_restart,"index",length,idim1)
515      IF (ierr/=NF90_NOERR) THEN
516        write(*,*)'create_restartphy: problem defining index dimension '
517        write(*,*)trim(nf90_strerror(ierr))
518        CALL abort_physic("create_restartphy","Failed defining index",1)
519      ENDIF
520
521      ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2)
522      IF (ierr/=NF90_NOERR) THEN
523        write(*,*)'create_restartphy: problem defining physical_points dimension '
524        write(*,*)trim(nf90_strerror(ierr))
525        CALL abort_physic("create_restartphy","Failed defining physical_points",1)
526      ENDIF
527
528      ierr=NF90_DEF_DIM(nid_restart,"subsurface_layers",nsoilmx,idim3)
529      IF (ierr/=NF90_NOERR) THEN
530        write(*,*)'create_restartphy: problem defining subsurface_layers dimension '
531        write(*,*)trim(nf90_strerror(ierr))
532        CALL abort_physic("create_restartphy","Failed defining subsurface_layers",1)
533      ENDIF
534
535      ierr=NF90_DEF_DIM(nid_restart,"nlayer_plus_1",klevp1,idim4)
536      IF (ierr/=NF90_NOERR) THEN
537        write(*,*)'create_restartphy: problem defining nlayer_plus_1 dimension '
538        write(*,*)trim(nf90_strerror(ierr))
539        CALL abort_physic("create_restartphy","Failed defining nlayer_plus_1",1)
540      ENDIF
541
542      if (nqtot>0) then
543        ! only define a tracer dimension if there are tracers
544        ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",nqtot,idim5)
545        IF (ierr/=NF90_NOERR) THEN
546          write(*,*)'create_restartphy: problem defining number_of_advected_fields dimension '
547          write(*,*)trim(nf90_strerror(ierr))
548          CALL abort_physic("create_restartphy","Failed defining number_of_advected_fields",1)
549        ENDIF
550      endif
551
552      ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
553      IF (ierr/=NF90_NOERR) THEN
554        write(*,*)'create_restartphy: problem defining nlayer dimension '
555        write(*,*)trim(nf90_strerror(ierr))
556        CALL abort_physic("create_restartphy","Failed defining nlayer",1)
557      ENDIF
558
559      ierr=NF90_DEF_DIM(nid_restart,"Time",NF90_UNLIMITED,idim7)
560      IF (ierr/=NF90_NOERR) THEN
561        write(*,*)'create_restartphy: problem defining Time dimension '
562        write(*,*)trim(nf90_strerror(ierr))
563        CALL abort_physic("create_restartphy","Failed defining Time",1)
564      ENDIF
565
566      ierr=NF90_DEF_DIM(nid_restart,"ocean_layers",nslay,idim8)
567      IF (ierr/=NF90_NOERR) THEN
568        write(*,*)'create_restartphy: problem defining oceanic layer dimension '
569        write(*,*)trim(nf90_strerror(ierr))
570        CALL abort_physic("create_restartphy","Failed defining nslay",1)
571      ENDIF
572
573      ierr=NF90_DEF_DIM(nid_restart,"descriptor",ldscrpt,idim10)
574      IF (ierr/=NF90_NOERR) THEN
575        write(*,*)'create_restartphy: problem defining descriptor dimension '
576        write(*,*)trim(nf90_strerror(ierr))
577        CALL abort_physic("create_restartphy","Failed defining descriptor",1)
578      ENDIF
579
580      ierr=NF90_DEF_DIM(nid_restart,"description_size",ndscrpt,idim11)
581      IF (ierr/=NF90_NOERR) THEN
582        write(*,*)'create_restartphy: problem defining description_size dimension '
583        write(*,*)trim(nf90_strerror(ierr))
584        CALL abort_physic("create_restartphy","Failed defining description_size",1)
585      ENDIF
586
587      ierr=NF90_ENDDEF(nid_restart)
588      IF (ierr/=NF90_NOERR) THEN
589        write(*,*)'create_restartphy: problem ending definition mode '
590        write(*,*)trim(nf90_strerror(ierr))
591        CALL abort_physic("create_restartphy","Failed ending definition mode",1)
592      ENDIF
593    ENDIF
594
595  END SUBROUTINE create_restartphy
596
597  SUBROUTINE open_restartphy(filename,nid_restart)
598  USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_WRITE, nf90_strerror
599  USE mod_phys_lmdz_para, only: is_master
600
601  IMPLICIT NONE
602    CHARACTER(LEN=*),INTENT(IN) :: filename
603    INTEGER,INTENT(INOUT)       :: nid_restart
604    INTEGER                     :: ierr
605
606    IF (is_master) THEN
607        ierr=NF90_OPEN(filename,NF90_WRITE,nid_restart)
608        IF (ierr/=NF90_NOERR) THEN
609          if (ierr.eq.2) then ! file does not exist
610            call create_restartphy(filename,nid_restart)
611            return
612          endif
613          write(*,*)'open_restartphy: problem opening file '//trim(filename)
614          write(*,*)trim(nf90_strerror(ierr))
615          CALL abort_physic("open_restartphy","Failed opening file",1)
616        ENDIF
617    ENDIF
618
619  END SUBROUTINE open_restartphy
620
621  SUBROUTINE close_restartphy(nid_restart)
622  USE netcdf, only: NF90_CLOSE
623  USE mod_phys_lmdz_para, only: is_master
624  IMPLICIT NONE
625    INTEGER,INTENT(IN) :: nid_restart
626    INTEGER            :: ierr
627
628    IF (is_master) THEN
629      ierr = NF90_CLOSE (nid_restart)
630    ENDIF
631
632  END SUBROUTINE close_restartphy
633
634  SUBROUTINE put_field_r1(nid_restart,field_name,title,field,time)
635  ! For a surface field
636  IMPLICIT NONE
637  INTEGER,INTENT(IN)             :: nid_restart
638  CHARACTER(LEN=*),INTENT(IN)    :: field_name
639  CHARACTER(LEN=*),INTENT(IN)    :: title
640  REAL,INTENT(IN)                :: field(:)
641  REAL,OPTIONAL,INTENT(IN)       :: time
642
643  IF (present(time)) THEN
644    ! if timeindex is present, it is a time-dependent variable
645    CALL put_field_rgen(nid_restart,field_name,title,field,1,time)
646  ELSE
647    CALL put_field_rgen(nid_restart,field_name,title,field,1)
648  ENDIF
649
650  END SUBROUTINE put_field_r1
651
652  SUBROUTINE put_field_r2(nid_restart,field_name,title,field,time)
653  ! For a "3D" horizontal-vertical field
654  IMPLICIT NONE
655  INTEGER,INTENT(IN)             :: nid_restart
656  CHARACTER(LEN=*),INTENT(IN)    :: field_name
657  CHARACTER(LEN=*),INTENT(IN)    :: title
658  REAL,INTENT(IN)                :: field(:,:)
659  REAL,OPTIONAL,INTENT(IN)       :: time
660
661  IF (present(time)) THEN
662    ! if timeindex is present, it is a time-dependent variable
663    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2),time)
664  ELSE
665    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2))
666  ENDIF
667
668  END SUBROUTINE put_field_r2
669
670  SUBROUTINE put_field_r3(nid_restart,field_name,title,field,time)
671  ! For a "4D" field surf/alt/??
672  IMPLICIT NONE
673  INTEGER,INTENT(IN)             :: nid_restart
674  CHARACTER(LEN=*),INTENT(IN)    :: field_name
675  CHARACTER(LEN=*),INTENT(IN)    :: title
676  REAL,INTENT(IN)                :: field(:,:,:)
677  REAL,OPTIONAL,INTENT(IN)       :: time
678
679  IF (present(time)) THEN
680    ! if timeindex is present, it is a time-dependent variable
681    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2)*size(field,3),&
682                        time)
683  ELSE
684    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2)*size(field,3))
685  ENDIF
686
687  END SUBROUTINE put_field_r3
688
689  SUBROUTINE put_field_rgen(nid_restart,field_name,title,field,field_size,time)
690  USE netcdf, ONLY: NF90_REDEF, NF90_ENDDEF, NF90_DEF_VAR, NF90_PUT_ATT, &
691                    NF90_INQ_VARID, NF90_PUT_VAR, NF90_STRERROR, &
692                    NF90_NOERR, NF90_FLOAT, NF90_DOUBLE
693  USE dimphy, ONLY: klon, klev, klevp1
694  USE comsoil_h, only: nsoilmx
695  USE mod_grid_phy_lmdz, ONLY: klon_glo
696  USE mod_phys_lmdz_para, ONLY: is_master, gather
697  USE geometry_mod, ONLY: ind_cell_glo
698!  USE slab_ice_h, only: noceanmx
699  USE ocean_slab_mod, ONLY: nslay
700
701  IMPLICIT NONE
702  INTEGER,INTENT(IN)             :: nid_restart
703  CHARACTER(LEN=*),INTENT(IN)    :: field_name
704  CHARACTER(LEN=*),INTENT(IN)    :: title
705  INTEGER,INTENT(IN)             :: field_size
706  REAL,INTENT(IN)                :: field(klon,field_size)
707  REAL,OPTIONAL,INTENT(IN)       :: time
708
709  REAL :: field_glo(klon_glo,field_size)
710  REAL :: field_glo_tmp(klon_glo,field_size)
711  INTEGER :: ind_cell_glo_glo(klon_glo) ! cell indexes on global grid
712
713  INTEGER                        :: ierr
714  INTEGER                        :: nvarid
715  INTEGER                        :: idim
716  INTEGER :: i
717
718    ! gather indexes on global grid
719    CALL gather(ind_cell_glo,ind_cell_glo_glo)
720    ! gather field on master
721    CALL gather(field,field_glo_tmp)
722
723    IF (is_master) THEN
724      ! reorder columns
725      DO i=1,klon_glo
726        field_glo(ind_cell_glo_glo(i),:)=field_glo_tmp(i,:)
727      ENDDO
728
729      IF (field_size==1) THEN
730        ! input is a 1D "surface field" array
731        if (.not.present(time)) then ! for a time-independent field
732          ierr=NF90_REDEF(nid_restart)
733#ifdef NC_DOUBLE
734          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
735                            (/idim2/),nvarid)
736#else
737          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
738                            (/idim2/),nvarid)
739#endif
740          if (ierr.ne.NF90_NOERR) then
741            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
742            write(*,*)trim(nf90_strerror(ierr))
743          endif
744          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
745          ierr=NF90_ENDDEF(nid_restart)
746          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
747        else
748          ! check if the variable has already been defined:
749          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
750          if (ierr/=NF90_NOERR) then ! variable not found, define it
751            ierr=NF90_REDEF(nid_restart)
752#ifdef NC_DOUBLE
753            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
754                              (/idim2,idim7/),nvarid)
755#else
756            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
757                              (/idim2,idim7/),nvarid)
758#endif
759            if (ierr.ne.NF90_NOERR) then
760              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
761              write(*,*)trim(nf90_strerror(ierr))
762            endif
763            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
764            ierr=NF90_ENDDEF(nid_restart)
765          endif
766          ! Write the variable
767          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
768                            start=(/1,timeindex/))
769        endif ! of if (.not.present(timeindex))
770
771      ELSE IF (field_size==klev) THEN
772        ! input is a 2D "atmospheric field" array
773        if (.not.present(time)) then ! for a time-independent field
774          ierr=NF90_REDEF(nid_restart)
775#ifdef NC_DOUBLE
776          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
777                            (/idim2,idim6/),nvarid)
778#else
779          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
780                            (/idim2,idim6/),nvarid)
781#endif
782          if (ierr.ne.NF90_NOERR) then
783            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
784            write(*,*)trim(nf90_strerror(ierr))
785          endif
786          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
787          ierr=NF90_ENDDEF(nid_restart)
788          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
789        else
790          ! check if the variable has already been defined:
791          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
792          if (ierr/=NF90_NOERR) then ! variable not found, define it
793            ierr=NF90_REDEF(nid_restart)
794#ifdef NC_DOUBLE
795            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
796                              (/idim2,idim6,idim7/),nvarid)
797#else
798            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
799                              (/idim2,idim6,idim7/),nvarid)
800#endif
801            if (ierr.ne.NF90_NOERR) then
802              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
803              write(*,*)trim(nf90_strerror(ierr))
804            endif
805            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
806            ierr=NF90_ENDDEF(nid_restart)
807          endif
808          ! Write the variable
809          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
810                            start=(/1,1,timeindex/))
811        endif ! of if (.not.present(time))
812
813      ELSE IF (field_size==klevp1) THEN
814        ! input is a 2D "interlayer atmospheric field" array
815        if (.not.present(time)) then ! for a time-independent field
816          ierr=NF90_REDEF(nid_restart)
817#ifdef NC_DOUBLE
818          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
819                            (/idim2,idim4/),nvarid)
820#else
821          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
822                            (/idim2,idim4/),nvarid)
823#endif
824          if (ierr.ne.NF90_NOERR) then
825            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
826            write(*,*)trim(nf90_strerror(ierr))
827          endif
828          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
829          ierr = NF90_ENDDEF(nid_restart)
830          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
831        else
832          ! check if the variable has already been defined:
833          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
834          if (ierr/=NF90_NOERR) then ! variable not found, define it
835            ierr=NF90_REDEF(nid_restart)
836#ifdef NC_DOUBLE
837            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
838                              (/idim2,idim4,idim7/),nvarid)
839#else
840            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
841                              (/idim2,idim4,idim7/),nvarid)
842#endif
843            if (ierr.ne.NF90_NOERR) then
844              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
845              write(*,*)trim(nf90_strerror(ierr))
846            endif
847            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
848            ierr=NF90_ENDDEF(nid_restart)
849          endif
850          ! Write the variable
851          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
852                            start=(/1,1,timeindex/))
853        endif ! of if (.not.present(timeindex))
854
855      ELSE IF (field_size==nsoilmx) THEN
856        ! input is a 2D "subsurface field" array
857        if (.not.present(time)) then ! for a time-independent field
858          ierr = NF90_REDEF(nid_restart)
859#ifdef NC_DOUBLE
860          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
861                            (/idim2,idim3/),nvarid)
862#else
863          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
864                            (/idim2,idim3/),nvarid)
865#endif
866          if (ierr.ne.NF90_NOERR) then
867            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
868            write(*,*)trim(nf90_strerror(ierr))
869          endif
870          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
871          ierr = NF90_ENDDEF(nid_restart)
872          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
873        else
874          ! check if the variable has already been defined:
875          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
876          if (ierr/=NF90_NOERR) then ! variable not found, define it
877            ierr=NF90_REDEF(nid_restart)
878#ifdef NC_DOUBLE
879            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
880                              (/idim2,idim3,idim7/),nvarid)
881#else
882            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
883                              (/idim2,idim3,idim7/),nvarid)
884#endif
885           if (ierr.ne.NF90_NOERR) then
886              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
887              write(*,*)trim(nf90_strerror(ierr))
888            endif
889            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
890            ierr=NF90_ENDDEF(nid_restart)
891          endif
892          ! Write the variable
893          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
894                            start=(/1,1,timeindex/))
895
896        endif ! of if (.not.present(time))
897
898      ELSE IF (field_size==nslay) THEN
899        ! input is a 2D "oceanic field" array
900        if (.not.present(time)) then ! for a time-independent field
901          ierr = NF90_REDEF(nid_restart)
902#ifdef NC_DOUBLE
903          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
904                            (/idim2,idim8/),nvarid)
905#else
906          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
907                            (/idim2,idim8/),nvarid)
908#endif
909          if (ierr.ne.NF90_NOERR) then
910            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
911            write(*,*)trim(nf90_strerror(ierr))
912          endif
913          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
914          ierr = NF90_ENDDEF(nid_restart)
915          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
916        else
917          ! check if the variable has already been defined:
918          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
919          if (ierr/=NF90_NOERR) then ! variable not found, define it
920            ierr=NF90_REDEF(nid_restart)
921#ifdef NC_DOUBLE
922            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
923                              (/idim2,idim8,idim7/),nvarid)
924#else
925            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
926                              (/idim2,idim8,idim7/),nvarid)
927#endif
928           if (ierr.ne.NF90_NOERR) then
929              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
930              write(*,*)trim(nf90_strerror(ierr))
931            endif
932            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
933            ierr=NF90_ENDDEF(nid_restart)
934          endif
935          ! Write the variable
936          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
937                            start=(/1,1,timeindex/))
938
939        endif ! of if (.not.present(time))
940
941
942      ELSE
943        PRINT *, "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name)
944        write(*,*) "  field_size =",field_size
945        CALL abort_physic("put_field_rgen","wrong field dimension",1)
946      ENDIF
947
948      ! Check the writting of field to file went OK
949      if (ierr.ne.NF90_NOERR) then
950        write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name)
951        write(*,*)trim(nf90_strerror(ierr))
952        CALL abort_physic("put_field_rgen","Failed writing field",1)
953      endif
954
955    ENDIF ! of IF (is_master)
956
957  END SUBROUTINE put_field_rgen
958
959  SUBROUTINE put_var_r0(nid_restart,var_name,title,var)
960  ! Put a scalar in file
961   IMPLICIT NONE
962     INTEGER,INTENT(IN)          :: nid_restart
963     CHARACTER(LEN=*),INTENT(IN) :: var_name
964     CHARACTER(LEN=*),INTENT(IN) :: title
965     REAL,INTENT(IN)             :: var
966     REAL                        :: varin(1)
967
968     varin(1)=var
969
970     CALL put_var_rgen(nid_restart,var_name,title,varin,size(varin))
971
972  END SUBROUTINE put_var_r0
973
974
975  SUBROUTINE put_var_r1(nid_restart,var_name,title,var)
976  ! Put a vector in file
977   IMPLICIT NONE
978     INTEGER,INTENT(IN)          :: nid_restart
979     CHARACTER(LEN=*),INTENT(IN) :: var_name
980     CHARACTER(LEN=*),INTENT(IN) :: title
981     REAL,INTENT(IN)             :: var(:)
982
983     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
984
985  END SUBROUTINE put_var_r1
986
987  SUBROUTINE put_var_r2(nid_restart,var_name,title,var)
988  ! Put a 2D field in file
989   IMPLICIT NONE
990     INTEGER,INTENT(IN)          :: nid_restart
991     CHARACTER(LEN=*),INTENT(IN) :: var_name
992     CHARACTER(LEN=*),INTENT(IN) :: title
993     REAL,INTENT(IN)             :: var(:,:)
994
995     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
996
997  END SUBROUTINE put_var_r2
998
999  SUBROUTINE put_var_r3(nid_restart,var_name,title,var)
1000  ! Put a 3D field in file
1001   IMPLICIT NONE
1002     INTEGER,INTENT(IN)          :: nid_restart
1003     CHARACTER(LEN=*),INTENT(IN) :: var_name
1004     CHARACTER(LEN=*),INTENT(IN) :: title
1005     REAL,INTENT(IN)             :: var(:,:,:)
1006
1007     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
1008
1009  END SUBROUTINE put_var_r3
1010
1011  SUBROUTINE put_var_rgen(nid_restart,var_name,title,var,var_size)
1012  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
1013                    NF90_FLOAT, NF90_DOUBLE, &
1014                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
1015                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
1016  USE comsoil_h, only: nsoilmx
1017  USE mod_phys_lmdz_para, only: is_master
1018!  USE slab_ice_h, only: noceanmx
1019  USE ocean_slab_mod, ONLY: nslay
1020  IMPLICIT NONE
1021    INTEGER,INTENT(IN)          :: nid_restart
1022     CHARACTER(LEN=*),INTENT(IN) :: var_name
1023     CHARACTER(LEN=*),INTENT(IN) :: title
1024     INTEGER,INTENT(IN)          :: var_size
1025     REAL,INTENT(IN)             :: var(var_size)
1026
1027     INTEGER :: ierr
1028     INTEGER :: nvarid
1029     INTEGER :: idim1d
1030     logical,save :: firsttime=.true.
1031!$OMP THREADPRIVATE(firsttime)
1032
1033    IF (is_master) THEN
1034
1035      IF (var_name=="Time") THEN
1036        ! Very specific case of "Time" variable
1037        if (firsttime) then
1038          ! Create the "Time variable"
1039          ierr=NF90_REDEF(nid_restart)
1040#ifdef NC_DOUBLE
1041          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,&
1042                            (/idim7/),nvarid)
1043#else
1044          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,&
1045                            (/idim7/),nvarid)
1046#endif
1047          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1048          ierr=NF90_ENDDEF(nid_restart)
1049
1050          firsttime=.false.
1051        endif
1052        ! Append "time" value
1053        ! get current length of "Time" dimension
1054        ierr=nf90_inq_dimid(nid_restart,var_name,idim1d)
1055        ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex)
1056        timeindex=timeindex+1
1057        ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid)
1058        ierr=NF90_PUT_VAR(nid_restart,nvarid,var,&
1059                           start=(/timeindex/))
1060        IF (ierr/=NF90_NOERR) THEN
1061          write(*,*)'put_var_rgen: problem writing Time'
1062          write(*,*)trim(nf90_strerror(ierr))
1063          CALL abort_physic("put_var_rgen","Failed to write Time",1)
1064        ENDIF
1065        return ! nothing left to do
1066      ELSEIF (var_size==length) THEN
1067        ! We know it is a "controle" kind of 1D array
1068        idim1d=idim1
1069      ELSEIF (var_size==nsoilmx) THEN
1070        ! We know it is an  "mlayer" kind of 1D array
1071        idim1d=idim3
1072      ELSEIF (var_size==nslay) THEN
1073        ! We know it is an  "mlayer" kind of 1D array
1074        idim1d=idim8
1075      ELSE
1076        PRINT *, "put_var_rgen error : wrong dimension"
1077        write(*,*) "  var_size =",var_size
1078        CALL abort_physic("put_var_rgen","Wrong variable dimension",1)
1079
1080      ENDIF ! of IF (var_size==length) THEN
1081
1082      ! Swich to NetCDF define mode
1083      ierr=NF90_REDEF (nid_restart)
1084      ! Define the variable
1085#ifdef NC_DOUBLE
1086      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid)
1087#else
1088      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid)
1089#endif
1090      ! Add a "title" attribute
1091      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1092      ! Swich out of define mode
1093      ierr=NF90_ENDDEF(nid_restart)
1094      ! Write variable to file
1095      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
1096      IF (ierr/=NF90_NOERR) THEN
1097        write(*,*)'put_var_rgen: problem writing '//trim(var_name)
1098        write(*,*)trim(nf90_strerror(ierr))
1099        CALL abort_physic("put_var_rgen","Failed writing variable",1)
1100      ENDIF
1101    ENDIF ! of IF (is_master)
1102
1103  END SUBROUTINE put_var_rgen
1104
1105  SUBROUTINE put_var_c1(nid_restart,var_name,title,var)
1106  ! Put a vector of characters in file
1107
1108  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
1109                    NF90_CHAR, &
1110                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
1111                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
1112  USE mod_phys_lmdz_para, only: is_master
1113
1114   IMPLICIT NONE
1115     INTEGER,INTENT(IN)          :: nid_restart
1116     CHARACTER(LEN=*),INTENT(IN) :: var_name
1117     CHARACTER(LEN=*),INTENT(IN) :: title
1118     CHARACTER(LEN=*),INTENT(IN) :: var(:)
1119
1120     INTEGER :: ierr
1121     INTEGER :: nvarid
1122     INTEGER :: idim1d_1, idim1d_2
1123     INTEGER :: var_size
1124
1125    IF (is_master) THEN
1126
1127      var_size = size(var)
1128      IF (var_size==ldscrpt) THEN
1129        ! We know it is a "controle descriptor" kind of 1D array
1130        idim1d_1=idim11
1131        idim1d_2=idim10
1132      ELSE
1133        PRINT *, "put_var_cgen error : wrong dimension"
1134        write(*,*) "  var_size =",var_size
1135        CALL abort_physic("put_var_cgen","Wrong variable dimension",1)
1136
1137      ENDIF ! of IF (var_size==length) THEN
1138
1139      ! Swich to NetCDF define mode
1140      ierr=NF90_REDEF (nid_restart)
1141      ! Define the variable
1142      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_CHAR,(/idim1d_1,idim1d_2/),nvarid)
1143      ! Add a "title" attribute
1144      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1145      ! Swich out of define mode
1146      ierr=NF90_ENDDEF(nid_restart)
1147      ! Write variable to file
1148      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
1149      IF (ierr/=NF90_NOERR) THEN
1150        write(*,*)'put_var_cgen: problem writing '//trim(var_name)
1151        write(*,*)trim(nf90_strerror(ierr))
1152        CALL abort_physic("put_var_cgen","Failed writing variable",1)
1153      ENDIF
1154    ENDIF ! of IF (is_master)
1155
1156  END SUBROUTINE put_var_c1
1157
1158!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1159
1160  SUBROUTINE create_restart1D(filename,nid_restart)
1161    USE netcdf, only: NF90_CREATE, NF90_CLOBBER, NF90_64BIT_OFFSET, &
1162                      NF90_NOERR, nf90_strerror, &
1163                      NF90_PUT_ATT, NF90_GLOBAL, NF90_DEF_DIM, &
1164                      NF90_UNLIMITED, NF90_ENDDEF
1165    USE mod_phys_lmdz_para, only: is_master
1166    USE mod_grid_phy_lmdz, only: klon_glo
1167    USE dimphy, only: klev, klevp1
1168    USE tracer_h, only: nqtot
1169    USE comsoil_h, only: nsoilmx
1170  !  USE slab_ice_h, only: noceanmx
1171    USE ocean_slab_mod, ONLY: nslay
1172
1173    IMPLICIT NONE
1174      CHARACTER(LEN=*),INTENT(IN) :: filename
1175      INTEGER,INTENT(INOUT)       :: nid_restart
1176      INTEGER                     :: ierr
1177
1178      IF (is_master) THEN
1179
1180          ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), &
1181                            nid_restart)
1182          IF (ierr/=NF90_NOERR) THEN
1183            write(*,*)'create_restart1D: problem creating file '//trim(filename)
1184            write(*,*)trim(nf90_strerror(ierr))
1185            CALL abort_physic("create_restart1D","Failed creating file",1)
1186          ENDIF
1187
1188        ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",&
1189                          "Physics start file")
1190        IF (ierr/=NF90_NOERR) THEN
1191          write(*,*)'create_restart1D: problem writing title '
1192          write(*,*)trim(nf90_strerror(ierr))
1193        ENDIF
1194
1195        ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2)
1196        IF (ierr/=NF90_NOERR) THEN
1197          write(*,*)'create_restart1D: problem defining physical_points dimension '
1198          write(*,*)trim(nf90_strerror(ierr))
1199          CALL abort_physic("create_restart1D","Failed defining physical_points",1)
1200        ENDIF
1201
1202        ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
1203        IF (ierr/=NF90_NOERR) THEN
1204          write(*,*)'create_restart1D: problem defining nlayer dimension '
1205          write(*,*)trim(nf90_strerror(ierr))
1206          CALL abort_physic("create_restart1D","Failed defining nlayer",1)
1207        ENDIF
1208
1209        ierr=NF90_ENDDEF(nid_restart)
1210        IF (ierr/=NF90_NOERR) THEN
1211          write(*,*)'create_restart1D: problem ending definition mode '
1212          write(*,*)trim(nf90_strerror(ierr))
1213          CALL abort_physic("create_restart1D","Failed ending definition mode",1)
1214        ENDIF
1215      ENDIF
1216
1217    END SUBROUTINE create_restart1D
1218
1219END MODULE iostart
Note: See TracBrowser for help on using the repository browser.