source: trunk/LMDZ.MARS/libf/phymars/iostart.F90 @ 3702

Last change on this file since 3702 was 3702, checked in by emillour, 3 months ago

Mars PCM:
Add reindexing of columns when reading/writing (re)startfi files. This is not
necessary with the lon-lat (LMDZ.COMMON) dynamical core, but required when
using DYNAMICO (where correspondance between dynamics and physics column
indexes changes with number of computing cores).
EM

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