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

Last change on this file since 3562 was 3562, checked in by mmaurice, 2 weeks ago

Generic PCM

1D restart operational: a restart1D.nc file is created that contains
psurf, tracers, winds and temperature profiles. A retartfi.nc file is
also created. Move those to and start1D.nc and startfi.nc and set
"restart" flag to .true. in rcm1d.def to restart from the files (also
make sure that day0 corresponds to the value in startfi.nc).

MM

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