source: trunk/LMDZ.PLUTO/libf/phypluto/iostart.F90

Last change on this file was 3750, checked in by afalco, 2 months ago

Pluto: write restartfi during the run, not just at the end. It is written every diagfi_output_rate.
AF

File size: 40.6 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
491  IMPLICIT NONE
492    CHARACTER(LEN=*),INTENT(IN) :: filename
493    INTEGER,INTENT(INOUT)       :: nid_restart
494    INTEGER                     :: ierr
495
496    IF (is_master) THEN
497
498        ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), &
499                          nid_restart)
500        IF (ierr/=NF90_NOERR) THEN
501          write(*,*)'create_restartphy: problem creating file '//trim(filename)
502          write(*,*)trim(nf90_strerror(ierr))
503          CALL abort_physic("create_restartphy","Failed creating file",1)
504        ENDIF
505
506      ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",&
507                        "Physics start file")
508      IF (ierr/=NF90_NOERR) THEN
509        write(*,*)'create_restartphy: problem writing title '
510        write(*,*)trim(nf90_strerror(ierr))
511      ENDIF
512
513      ierr=NF90_DEF_DIM(nid_restart,"index",length,idim1)
514      IF (ierr/=NF90_NOERR) THEN
515        write(*,*)'create_restartphy: problem defining index dimension '
516        write(*,*)trim(nf90_strerror(ierr))
517        CALL abort_physic("create_restartphy","Failed defining index",1)
518      ENDIF
519
520      ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2)
521      IF (ierr/=NF90_NOERR) THEN
522        write(*,*)'create_restartphy: problem defining physical_points dimension '
523        write(*,*)trim(nf90_strerror(ierr))
524        CALL abort_physic("create_restartphy","Failed defining physical_points",1)
525      ENDIF
526
527      ierr=NF90_DEF_DIM(nid_restart,"subsurface_layers",nsoilmx,idim3)
528      IF (ierr/=NF90_NOERR) THEN
529        write(*,*)'create_restartphy: problem defining subsurface_layers dimension '
530        write(*,*)trim(nf90_strerror(ierr))
531        CALL abort_physic("create_restartphy","Failed defining subsurface_layers",1)
532      ENDIF
533
534      ierr=NF90_DEF_DIM(nid_restart,"nlayer_plus_1",klevp1,idim4)
535      IF (ierr/=NF90_NOERR) THEN
536        write(*,*)'create_restartphy: problem defining nlayer_plus_1 dimension '
537        write(*,*)trim(nf90_strerror(ierr))
538        CALL abort_physic("create_restartphy","Failed defining nlayer_plus_1",1)
539      ENDIF
540
541      if (nqtot>0) then
542        ! only define a tracer dimension if there are tracers
543        ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",nqtot,idim5)
544        IF (ierr/=NF90_NOERR) THEN
545          write(*,*)'create_restartphy: problem defining number_of_advected_fields dimension '
546          write(*,*)trim(nf90_strerror(ierr))
547          CALL abort_physic("create_restartphy","Failed defining number_of_advected_fields",1)
548        ENDIF
549      endif
550
551      ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
552      IF (ierr/=NF90_NOERR) THEN
553        write(*,*)'create_restartphy: problem defining nlayer dimension '
554        write(*,*)trim(nf90_strerror(ierr))
555        CALL abort_physic("create_restartphy","Failed defining nlayer",1)
556      ENDIF
557
558      ierr=NF90_DEF_DIM(nid_restart,"Time",NF90_UNLIMITED,idim7)
559      IF (ierr/=NF90_NOERR) THEN
560        write(*,*)'create_restartphy: problem defining Time dimension '
561        write(*,*)trim(nf90_strerror(ierr))
562        CALL abort_physic("create_restartphy","Failed defining Time",1)
563      ENDIF
564
565      ierr=NF90_DEF_DIM(nid_restart,"descriptor",ldscrpt,idim10)
566      IF (ierr/=NF90_NOERR) THEN
567        write(*,*)'create_restartphy: problem defining descriptor dimension '
568        write(*,*)trim(nf90_strerror(ierr))
569        CALL abort_physic("create_restartphy","Failed defining descriptor",1)
570      ENDIF
571
572      ierr=NF90_DEF_DIM(nid_restart,"description_size",ndscrpt,idim11)
573      IF (ierr/=NF90_NOERR) THEN
574        write(*,*)'create_restartphy: problem defining description_size dimension '
575        write(*,*)trim(nf90_strerror(ierr))
576        CALL abort_physic("create_restartphy","Failed defining description_size",1)
577      ENDIF
578
579      ierr=NF90_ENDDEF(nid_restart)
580      IF (ierr/=NF90_NOERR) THEN
581        write(*,*)'create_restartphy: problem ending definition mode '
582        write(*,*)trim(nf90_strerror(ierr))
583        CALL abort_physic("create_restartphy","Failed ending definition mode",1)
584      ENDIF
585    ENDIF
586
587  END SUBROUTINE create_restartphy
588
589  SUBROUTINE open_restartphy(filename,nid_restart)
590  USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_WRITE, nf90_strerror
591  USE mod_phys_lmdz_para, only: is_master
592
593  IMPLICIT NONE
594    CHARACTER(LEN=*),INTENT(IN) :: filename
595    INTEGER,INTENT(INOUT)       :: nid_restart
596    INTEGER                     :: ierr
597
598    IF (is_master) THEN
599        ierr=NF90_OPEN(filename,NF90_WRITE,nid_restart)
600        IF (ierr/=NF90_NOERR) THEN
601          if (ierr.eq.2) then ! file does not exist
602            call create_restartphy(filename,nid_restart)
603            return
604          endif
605          write(*,*)'open_restartphy: problem opening file '//trim(filename)
606          write(*,*)trim(nf90_strerror(ierr))
607          CALL abort_physic("open_restartphy","Failed opening file",1)
608        ENDIF
609    ENDIF
610
611  END SUBROUTINE open_restartphy
612
613  SUBROUTINE close_restartphy(nid_restart)
614  USE netcdf, only: NF90_CLOSE
615  USE mod_phys_lmdz_para, only: is_master
616  IMPLICIT NONE
617    INTEGER,INTENT(IN) :: nid_restart
618    INTEGER            :: ierr
619
620    IF (is_master) THEN
621      ierr = NF90_CLOSE (nid_restart)
622    ENDIF
623
624  END SUBROUTINE close_restartphy
625
626  SUBROUTINE put_field_r1(nid_restart,field_name,title,field,time)
627  ! For a surface field
628  IMPLICIT NONE
629  INTEGER,INTENT(IN)             :: nid_restart
630  CHARACTER(LEN=*),INTENT(IN)    :: field_name
631  CHARACTER(LEN=*),INTENT(IN)    :: title
632  REAL,INTENT(IN)                :: field(:)
633  REAL,OPTIONAL,INTENT(IN)       :: time
634
635  IF (present(time)) THEN
636    ! if timeindex is present, it is a time-dependent variable
637    CALL put_field_rgen(nid_restart,field_name,title,field,1,time)
638  ELSE
639    CALL put_field_rgen(nid_restart,field_name,title,field,1)
640  ENDIF
641
642  END SUBROUTINE put_field_r1
643
644  SUBROUTINE put_field_r2(nid_restart,field_name,title,field,time)
645  ! For a "3D" horizontal-vertical field
646  IMPLICIT NONE
647  INTEGER,INTENT(IN)             :: nid_restart
648  CHARACTER(LEN=*),INTENT(IN)    :: field_name
649  CHARACTER(LEN=*),INTENT(IN)    :: title
650  REAL,INTENT(IN)                :: field(:,:)
651  REAL,OPTIONAL,INTENT(IN)       :: time
652
653  IF (present(time)) THEN
654    ! if timeindex is present, it is a time-dependent variable
655    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2),time)
656  ELSE
657    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2))
658  ENDIF
659
660  END SUBROUTINE put_field_r2
661
662  SUBROUTINE put_field_r3(nid_restart,field_name,title,field,time)
663  ! For a "4D" field surf/alt/??
664  IMPLICIT NONE
665  INTEGER,INTENT(IN)             :: nid_restart
666  CHARACTER(LEN=*),INTENT(IN)    :: field_name
667  CHARACTER(LEN=*),INTENT(IN)    :: title
668  REAL,INTENT(IN)                :: field(:,:,:)
669  REAL,OPTIONAL,INTENT(IN)       :: time
670
671  IF (present(time)) THEN
672    ! if timeindex is present, it is a time-dependent variable
673    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2)*size(field,3),&
674                        time)
675  ELSE
676    CALL put_field_rgen(nid_restart,field_name,title,field,size(field,2)*size(field,3))
677  ENDIF
678
679  END SUBROUTINE put_field_r3
680
681  SUBROUTINE put_field_rgen(nid_restart,field_name,title,field,field_size,time)
682  USE netcdf, ONLY: NF90_REDEF, NF90_ENDDEF, NF90_DEF_VAR, NF90_PUT_ATT, &
683                    NF90_INQ_VARID, NF90_PUT_VAR, NF90_STRERROR, &
684                    NF90_NOERR, NF90_FLOAT, NF90_DOUBLE
685  USE dimphy, ONLY: klon, klev, klevp1
686  USE comsoil_h, only: nsoilmx
687  USE mod_grid_phy_lmdz, ONLY: klon_glo
688  USE mod_phys_lmdz_para, ONLY: is_master, gather
689  USE geometry_mod, ONLY: ind_cell_glo
690!  USE slab_ice_h, only: noceanmx
691
692  IMPLICIT NONE
693  INTEGER,INTENT(IN)             :: nid_restart
694  CHARACTER(LEN=*),INTENT(IN)    :: field_name
695  CHARACTER(LEN=*),INTENT(IN)    :: title
696  INTEGER,INTENT(IN)             :: field_size
697  REAL,INTENT(IN)                :: field(klon,field_size)
698  REAL,OPTIONAL,INTENT(IN)       :: time
699
700  REAL :: field_glo(klon_glo,field_size)
701  REAL :: field_glo_tmp(klon_glo,field_size)
702  INTEGER :: ind_cell_glo_glo(klon_glo) ! cell indexes on global grid
703
704  INTEGER                        :: ierr
705  INTEGER                        :: nvarid
706  INTEGER                        :: idim
707  INTEGER :: i
708
709    ! gather indexes on global grid
710    CALL gather(ind_cell_glo,ind_cell_glo_glo)
711    ! gather field on master
712    CALL gather(field,field_glo_tmp)
713
714    IF (is_master) THEN
715      ! reorder columns
716      DO i=1,klon_glo
717        field_glo(ind_cell_glo_glo(i),:)=field_glo_tmp(i,:)
718      ENDDO
719
720      IF (field_size==1) THEN
721        ! input is a 1D "surface field" array
722        ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
723        if (ierr/=NF90_NOERR.or..not.present(time)) then ! for a time-independent field
724          ierr=NF90_REDEF(nid_restart)
725#ifdef NC_DOUBLE
726          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
727                            (/idim2/),nvarid)
728#else
729          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
730                            (/idim2/),nvarid)
731#endif
732          if (ierr.ne.NF90_NOERR) then
733            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
734            write(*,*)trim(nf90_strerror(ierr))
735          endif
736          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
737          ierr=NF90_ENDDEF(nid_restart)
738          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
739        else
740          ! check if the variable has already been defined:
741          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
742          if (ierr/=NF90_NOERR) then ! variable not found, define it
743            ierr=NF90_REDEF(nid_restart)
744#ifdef NC_DOUBLE
745            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
746                              (/idim2,idim7/),nvarid)
747#else
748            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
749                              (/idim2,idim7/),nvarid)
750#endif
751            if (ierr.ne.NF90_NOERR) then
752              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
753              write(*,*)trim(nf90_strerror(ierr))
754            endif
755            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
756            ierr=NF90_ENDDEF(nid_restart)
757          endif
758          ! Write the variable
759          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
760                            start=(/1,timeindex/))
761        endif ! of if (.not.present(timeindex))
762
763      ELSE IF (field_size==klev) THEN
764        ! input is a 2D "atmospheric field" array
765        ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
766        if (ierr/=NF90_NOERR.or..not.present(time)) then ! for a time-independent field
767          ierr=NF90_REDEF(nid_restart)
768#ifdef NC_DOUBLE
769          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
770                            (/idim2,idim6/),nvarid)
771#else
772          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
773                            (/idim2,idim6/),nvarid)
774#endif
775          if (ierr.ne.NF90_NOERR) then
776            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
777            write(*,*)trim(nf90_strerror(ierr))
778          endif
779          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
780          ierr=NF90_ENDDEF(nid_restart)
781          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
782        else
783          ! check if the variable has already been defined:
784          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
785          if (ierr/=NF90_NOERR) then ! variable not found, define it
786            ierr=NF90_REDEF(nid_restart)
787#ifdef NC_DOUBLE
788            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
789                              (/idim2,idim6,idim7/),nvarid)
790#else
791            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
792                              (/idim2,idim6,idim7/),nvarid)
793#endif
794            if (ierr.ne.NF90_NOERR) then
795              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
796              write(*,*)trim(nf90_strerror(ierr))
797            endif
798            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
799            ierr=NF90_ENDDEF(nid_restart)
800          endif
801          ! Write the variable
802          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
803                            start=(/1,1,timeindex/))
804        endif ! of if (.not.present(time))
805
806      ELSE IF (field_size==klevp1) THEN
807        ! input is a 2D "interlayer atmospheric field" array
808        ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
809        if (ierr/=NF90_NOERR.or..not.present(time)) then ! for a time-independent field
810          ierr=NF90_REDEF(nid_restart)
811#ifdef NC_DOUBLE
812          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
813                            (/idim2,idim4/),nvarid)
814#else
815          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
816                            (/idim2,idim4/),nvarid)
817#endif
818          if (ierr.ne.NF90_NOERR) then
819            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
820            write(*,*)trim(nf90_strerror(ierr))
821          endif
822          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
823          ierr = NF90_ENDDEF(nid_restart)
824          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
825        else
826          ! check if the variable has already been defined:
827          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
828          if (ierr/=NF90_NOERR) then ! variable not found, define it
829            ierr=NF90_REDEF(nid_restart)
830#ifdef NC_DOUBLE
831            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
832                              (/idim2,idim4,idim7/),nvarid)
833#else
834            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
835                              (/idim2,idim4,idim7/),nvarid)
836#endif
837            if (ierr.ne.NF90_NOERR) then
838              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
839              write(*,*)trim(nf90_strerror(ierr))
840            endif
841            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
842            ierr=NF90_ENDDEF(nid_restart)
843          endif
844          ! Write the variable
845          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
846                            start=(/1,1,timeindex/))
847        endif ! of if (.not.present(timeindex))
848
849      ELSE IF (field_size==nsoilmx) THEN
850        ! input is a 2D "subsurface field" array
851        ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
852        if (ierr/=NF90_NOERR.or..not.present(time)) then ! for a time-independent field
853          ierr = NF90_REDEF(nid_restart)
854#ifdef NC_DOUBLE
855          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
856                            (/idim2,idim3/),nvarid)
857#else
858          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
859                            (/idim2,idim3/),nvarid)
860#endif
861          if (ierr.ne.NF90_NOERR) then
862            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
863            write(*,*)trim(nf90_strerror(ierr))
864          endif
865          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
866          ierr = NF90_ENDDEF(nid_restart)
867          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
868        else
869          ! check if the variable has already been defined:
870          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
871          if (ierr/=NF90_NOERR) then ! variable not found, define it
872            ierr=NF90_REDEF(nid_restart)
873#ifdef NC_DOUBLE
874            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
875                              (/idim2,idim3,idim7/),nvarid)
876#else
877            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
878                              (/idim2,idim3,idim7/),nvarid)
879#endif
880           if (ierr.ne.NF90_NOERR) then
881              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
882              write(*,*)trim(nf90_strerror(ierr))
883            endif
884            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
885            ierr=NF90_ENDDEF(nid_restart)
886          endif
887          ! Write the variable
888          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
889                            start=(/1,1,timeindex/))
890
891        endif ! of if (.not.present(time))
892
893      ELSE
894        PRINT *, "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name)
895        write(*,*) "  field_size =",field_size
896        CALL abort_physic("put_field_rgen","wrong field dimension",1)
897      ENDIF
898
899      ! Check the writing of field to file went OK
900      if (ierr.ne.NF90_NOERR) then
901        write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name)
902        write(*,*)trim(nf90_strerror(ierr))
903        CALL abort_physic("put_field_rgen","Failed writing field",1)
904      endif
905
906    ENDIF ! of IF (is_master)
907
908  END SUBROUTINE put_field_rgen
909
910  SUBROUTINE put_var_r0(nid_restart,var_name,title,var)
911  ! Put a scalar in file
912   IMPLICIT NONE
913     INTEGER,INTENT(IN)          :: nid_restart
914     CHARACTER(LEN=*),INTENT(IN) :: var_name
915     CHARACTER(LEN=*),INTENT(IN) :: title
916     REAL,INTENT(IN)             :: var
917     REAL                        :: varin(1)
918
919     varin(1)=var
920
921     CALL put_var_rgen(nid_restart,var_name,title,varin,size(varin))
922
923  END SUBROUTINE put_var_r0
924
925
926  SUBROUTINE put_var_r1(nid_restart,var_name,title,var)
927  ! Put a vector in file
928   IMPLICIT NONE
929     INTEGER,INTENT(IN)          :: nid_restart
930     CHARACTER(LEN=*),INTENT(IN) :: var_name
931     CHARACTER(LEN=*),INTENT(IN) :: title
932     REAL,INTENT(IN)             :: var(:)
933
934     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
935
936  END SUBROUTINE put_var_r1
937
938  SUBROUTINE put_var_r2(nid_restart,var_name,title,var)
939  ! Put a 2D field in file
940   IMPLICIT NONE
941     INTEGER,INTENT(IN)          :: nid_restart
942     CHARACTER(LEN=*),INTENT(IN) :: var_name
943     CHARACTER(LEN=*),INTENT(IN) :: title
944     REAL,INTENT(IN)             :: var(:,:)
945
946     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
947
948  END SUBROUTINE put_var_r2
949
950  SUBROUTINE put_var_r3(nid_restart,var_name,title,var)
951  ! Put a 3D field in file
952   IMPLICIT NONE
953     INTEGER,INTENT(IN)          :: nid_restart
954     CHARACTER(LEN=*),INTENT(IN) :: var_name
955     CHARACTER(LEN=*),INTENT(IN) :: title
956     REAL,INTENT(IN)             :: var(:,:,:)
957
958     CALL put_var_rgen(nid_restart,var_name,title,var,size(var))
959
960  END SUBROUTINE put_var_r3
961
962  SUBROUTINE put_var_rgen(nid_restart,var_name,title,var,var_size)
963  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
964                    NF90_FLOAT, NF90_DOUBLE, &
965                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
966                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
967  USE comsoil_h, only: nsoilmx
968  USE mod_phys_lmdz_para, only: is_master
969!  USE slab_ice_h, only: noceanmx
970  IMPLICIT NONE
971    INTEGER,INTENT(IN)          :: nid_restart
972     CHARACTER(LEN=*),INTENT(IN) :: var_name
973     CHARACTER(LEN=*),INTENT(IN) :: title
974     INTEGER,INTENT(IN)          :: var_size
975     REAL,INTENT(IN)             :: var(var_size)
976
977     INTEGER :: ierr
978     INTEGER :: nvarid
979     INTEGER :: idim1d
980     logical,save :: firsttime=.true.
981!$OMP THREADPRIVATE(firsttime)
982
983    IF (is_master) THEN
984
985      IF (var_name=="Time") THEN
986        ! Very specific case of "Time" variable
987        if (firsttime) then
988          ! Create the "Time variable"
989          ierr=NF90_REDEF(nid_restart)
990#ifdef NC_DOUBLE
991          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,&
992                            (/idim7/),nvarid)
993#else
994          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,&
995                            (/idim7/),nvarid)
996#endif
997          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
998          ierr=NF90_ENDDEF(nid_restart)
999
1000          firsttime=.false.
1001        endif
1002        ! Append "time" value
1003        ! get current length of "Time" dimension
1004        ierr=nf90_inq_dimid(nid_restart,var_name,idim1d)
1005        ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex)
1006        timeindex=timeindex+1
1007        ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid)
1008        ierr=NF90_PUT_VAR(nid_restart,nvarid,var,&
1009                           start=(/timeindex/))
1010        IF (ierr/=NF90_NOERR) THEN
1011          write(*,*)'put_var_rgen: problem writing Time'
1012          write(*,*)trim(nf90_strerror(ierr))
1013          CALL abort_physic("put_var_rgen","Failed to write Time",1)
1014        ENDIF
1015        return ! nothing left to do
1016      ELSEIF (var_size==length) THEN
1017        ! We know it is a "controle" kind of 1D array
1018        idim1d=idim1
1019      ELSEIF (var_size==nsoilmx) THEN
1020        ! We know it is an  "mlayer" kind of 1D array
1021        idim1d=idim3
1022      ELSE
1023        PRINT *, "put_var_rgen error : wrong dimension"
1024        write(*,*) "  var_size =",var_size
1025        CALL abort_physic("put_var_rgen","Wrong variable dimension",1)
1026
1027      ENDIF ! of IF (var_size==length) THEN
1028
1029      ! Swich to NetCDF define mode
1030      ierr=NF90_REDEF (nid_restart)
1031      ! Define the variable
1032#ifdef NC_DOUBLE
1033      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid)
1034#else
1035      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid)
1036#endif
1037      ! Add a "title" attribute
1038      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1039      ! Swich out of define mode
1040      ierr=NF90_ENDDEF(nid_restart)
1041      ! Write variable to file
1042      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
1043      IF (ierr/=NF90_NOERR) THEN
1044        write(*,*)'put_var_rgen: problem writing '//trim(var_name)
1045        write(*,*)trim(nf90_strerror(ierr))
1046        CALL abort_physic("put_var_rgen","Failed writing variable",1)
1047      ENDIF
1048    ENDIF ! of IF (is_master)
1049
1050  END SUBROUTINE put_var_rgen
1051
1052  SUBROUTINE put_var_c1(nid_restart,var_name,title,var)
1053  ! Put a vector of characters in file
1054
1055  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
1056                    NF90_CHAR, &
1057                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
1058                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
1059  USE mod_phys_lmdz_para, only: is_master
1060
1061   IMPLICIT NONE
1062     INTEGER,INTENT(IN)          :: nid_restart
1063     CHARACTER(LEN=*),INTENT(IN) :: var_name
1064     CHARACTER(LEN=*),INTENT(IN) :: title
1065     CHARACTER(LEN=*),INTENT(IN) :: var(:)
1066
1067     INTEGER :: ierr
1068     INTEGER :: nvarid
1069     INTEGER :: idim1d_1, idim1d_2
1070     INTEGER :: var_size
1071
1072    IF (is_master) THEN
1073
1074      var_size = size(var)
1075      IF (var_size==ldscrpt) THEN
1076        ! We know it is a "controle descriptor" kind of 1D array
1077        idim1d_1=idim11
1078        idim1d_2=idim10
1079      ELSE
1080        PRINT *, "put_var_cgen error : wrong dimension"
1081        write(*,*) "  var_size =",var_size
1082        CALL abort_physic("put_var_cgen","Wrong variable dimension",1)
1083
1084      ENDIF ! of IF (var_size==length) THEN
1085
1086      ! Swich to NetCDF define mode
1087      ierr=NF90_REDEF (nid_restart)
1088      ! Define the variable
1089      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_CHAR,(/idim1d_1,idim1d_2/),nvarid)
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_cgen: problem writing '//trim(var_name)
1098        write(*,*)trim(nf90_strerror(ierr))
1099        CALL abort_physic("put_var_cgen","Failed writing variable",1)
1100      ENDIF
1101    ENDIF ! of IF (is_master)
1102
1103  END SUBROUTINE put_var_c1
1104
1105!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1106
1107  SUBROUTINE create_restart1D(filename,nid_restart)
1108    USE netcdf, only: NF90_CREATE, NF90_CLOBBER, NF90_64BIT_OFFSET, &
1109                      NF90_NOERR, nf90_strerror, &
1110                      NF90_PUT_ATT, NF90_GLOBAL, NF90_DEF_DIM, &
1111                      NF90_UNLIMITED, NF90_ENDDEF
1112    USE mod_phys_lmdz_para, only: is_master
1113    USE mod_grid_phy_lmdz, only: klon_glo
1114    USE dimphy, only: klev, klevp1
1115    USE tracer_h, only: nqtot
1116    USE comsoil_h, only: nsoilmx
1117  !  USE slab_ice_h, only: noceanmx
1118
1119    IMPLICIT NONE
1120      CHARACTER(LEN=*),INTENT(IN) :: filename
1121      INTEGER,INTENT(INOUT)       :: nid_restart
1122      INTEGER                     :: ierr
1123
1124      IF (is_master) THEN
1125
1126          ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), &
1127                            nid_restart)
1128          IF (ierr/=NF90_NOERR) THEN
1129            write(*,*)'create_restart1D: problem creating file '//trim(filename)
1130            write(*,*)trim(nf90_strerror(ierr))
1131            CALL abort_physic("create_restart1D","Failed creating file",1)
1132          ENDIF
1133
1134        ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",&
1135                          "Physics start file")
1136        IF (ierr/=NF90_NOERR) THEN
1137          write(*,*)'create_restart1D: problem writing title '
1138          write(*,*)trim(nf90_strerror(ierr))
1139        ENDIF
1140
1141        ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2)
1142        IF (ierr/=NF90_NOERR) THEN
1143          write(*,*)'create_restart1D: problem defining physical_points dimension '
1144          write(*,*)trim(nf90_strerror(ierr))
1145          CALL abort_physic("create_restart1D","Failed defining physical_points",1)
1146        ENDIF
1147
1148        ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
1149        IF (ierr/=NF90_NOERR) THEN
1150          write(*,*)'create_restart1D: problem defining nlayer dimension '
1151          write(*,*)trim(nf90_strerror(ierr))
1152          CALL abort_physic("create_restart1D","Failed defining nlayer",1)
1153        ENDIF
1154
1155        ierr=NF90_ENDDEF(nid_restart)
1156        IF (ierr/=NF90_NOERR) THEN
1157          write(*,*)'create_restart1D: problem ending definition mode '
1158          write(*,*)trim(nf90_strerror(ierr))
1159          CALL abort_physic("create_restart1D","Failed ending definition mode",1)
1160        ENDIF
1161      ENDIF
1162
1163    END SUBROUTINE create_restart1D
1164
1165END MODULE iostart
Note: See TracBrowser for help on using the repository browser.