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

Last change on this file since 3509 was 3509, checked in by jbclement, 13 days ago

Dynamic + Mars PCM:
Addition of the description for the 'controle' array in the "start.nc" and "startfi.nc" files. It is given by the variable 'controle_descriptor' whose the element 'controle_descriptor(i)' explains 'controle(i)'.
JBC

File size: 40.0 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    INTEGER,PARAMETER :: length = 100 ! size of tab_cntrl array
23    INTEGER,PARAMETER :: ldscrpt = 35 ! size of dscrpt_tab_cntrl array
24    INTEGER,PARAMETER :: ndscrpt = 50 ! size of characters in dscrpt_tab_cntrl array
25
26    INTERFACE get_field
27      MODULE PROCEDURE Get_field_r1,Get_field_r2,Get_field_r3
28    END INTERFACE get_field
29
30    INTERFACE get_var
31      MODULE PROCEDURE get_var_r0,Get_var_r1,Get_var_r2,Get_var_r3
32    END INTERFACE get_var
33
34    INTERFACE put_field
35      MODULE PROCEDURE put_field_r1,put_field_r2,put_field_r3
36    END INTERFACE put_field
37
38    INTERFACE put_var
39      MODULE PROCEDURE put_var_r0,put_var_r1,put_var_r2,put_var_r3, put_var_c1
40    END INTERFACE put_var
41
42    PUBLIC nid_start, length, ldscrpt, ndscrpt
43    PUBLIC get_field,get_var,put_field,put_var
44    PUBLIC inquire_dimension, inquire_dimension_length
45    PUBLIC inquire_field, inquire_field_ndims
46    PUBLIC open_startphy,close_startphy,open_restartphy,close_restartphy
47
48CONTAINS
49
50  SUBROUTINE open_startphy(filename)
51  USE netcdf, only: NF90_OPEN, NF90_NOERR, NF90_NOWRITE, nf90_strerror
52  USE mod_phys_lmdz_para, only: is_master, bcast
53  IMPLICIT NONE
54    CHARACTER(LEN=*) :: filename
55    INTEGER          :: ierr
56
57    IF (is_master) THEN
58      ierr = NF90_OPEN (filename, NF90_NOWRITE, nid_start)
59      IF (ierr.NE.NF90_NOERR) THEN
60        write(*,*)'open_startphy: problem opening file '//trim(filename)
61        write(*,*)trim(nf90_strerror(ierr))
62        CALL abort_physic("open_startphy","Cannot open file",1)
63      ENDIF
64    ENDIF
65   
66    CALL bcast(nid_start) ! tell all procs about nid_start
67 
68  END SUBROUTINE open_startphy
69
70  SUBROUTINE close_startphy
71  USE netcdf, only: NF90_CLOSE
72  USE mod_phys_lmdz_para, only: is_master
73  IMPLICIT NONE
74    INTEGER          :: ierr
75
76    IF (is_master) THEN
77        ierr = NF90_CLOSE (nid_start)
78    ENDIF
79
80  END SUBROUTINE close_startphy
81
82
83  FUNCTION inquire_field(Field_name)
84  ! check if a given field is present in the input file
85  USE netcdf, only: NF90_INQ_VARID, NF90_NOERR
86  USE mod_phys_lmdz_para, only: is_master, bcast
87  IMPLICIT NONE
88    CHARACTER(LEN=*),INTENT(IN) :: Field_name
89    LOGICAL :: inquire_field
90    INTEGER :: varid
91    INTEGER :: ierr
92   
93    IF (is_master) THEN
94      ierr=NF90_INQ_VARID(nid_start,Field_name,varid)
95      IF (ierr==NF90_NOERR) THEN
96        Inquire_field=.TRUE.
97      ELSE
98        Inquire_field=.FALSE.
99      ENDIF
100    ENDIF
101
102    CALL bcast(inquire_field)
103
104  END FUNCTION inquire_field
105
106
107  FUNCTION inquire_field_ndims(Field_name)
108  ! give the number of dimensions of "Field_name" stored in the input file
109  USE netcdf, only: nf90_inq_varid, nf90_inquire_variable, &
110                    NF90_NOERR, nf90_strerror
111  USE mod_phys_lmdz_para, only: is_master, bcast
112  IMPLICIT NONE
113    CHARACTER(LEN=*),INTENT(IN) :: Field_name
114    INTEGER :: inquire_field_ndims
115    INTEGER :: varid
116    INTEGER :: ierr
117   
118    IF (is_master) THEN
119      ierr=nf90_inq_varid(nid_start,Field_name,varid)
120      ierr=nf90_inquire_variable(nid_start,varid,&
121                                  ndims=inquire_field_ndims)
122      IF (ierr.NE.NF90_NOERR) THEN
123        write(*,*)'inquire_field_ndims: problem obtaining ndims of '&
124                  //trim(field_name)
125        write(*,*)trim(nf90_strerror(ierr))
126        CALL abort_physic("inquire_field_ndims","Failed to get ndims",1)
127      ENDIF
128    ENDIF
129
130    CALL bcast(inquire_field_ndims)
131
132  END FUNCTION inquire_field_ndims
133
134
135  FUNCTION inquire_dimension(Field_name)
136  ! check if a given dimension is present in the input file
137  USE netcdf, only: nf90_inq_dimid, NF90_NOERR
138  USE mod_phys_lmdz_para, only: is_master, bcast
139  IMPLICIT NONE
140    CHARACTER(LEN=*),INTENT(IN) :: Field_name
141    LOGICAL :: inquire_dimension
142    INTEGER :: varid
143    INTEGER :: ierr
144   
145    IF (is_master) THEN
146      ierr=NF90_INQ_DIMID(nid_start,Field_name,varid)
147      IF (ierr==NF90_NOERR) THEN
148        Inquire_dimension=.TRUE.
149      ELSE
150        Inquire_dimension=.FALSE.
151      ENDIF
152    ENDIF
153
154    CALL bcast(inquire_dimension)
155
156  END FUNCTION inquire_dimension
157
158  FUNCTION inquire_dimension_length(Field_name)
159  ! give the length of the "Field_name" dimension stored in the input file
160  USE netcdf, only: nf90_inquire_dimension, nf90_inq_dimid, &
161                    NF90_NOERR, nf90_strerror
162  USE mod_phys_lmdz_para, only: is_master, bcast
163  IMPLICIT NONE
164    CHARACTER(LEN=*),INTENT(IN) :: Field_name
165    INTEGER :: inquire_dimension_length
166    INTEGER :: varid
167    INTEGER :: ierr
168   
169    IF (is_master) THEN
170      ierr=nf90_inq_dimid(nid_start,Field_name,varid)
171      ierr=nf90_inquire_dimension(nid_start,varid,&
172                                  len=inquire_dimension_length)
173      IF (ierr.NE.NF90_NOERR) THEN
174        write(*,*)'inquire_field_length: problem obtaining length of '&
175                  //trim(field_name)
176        write(*,*)trim(nf90_strerror(ierr))
177        CALL abort_physic("inquire_field_ndims","Failed to get length",1)
178      ENDIF
179    ENDIF
180
181    CALL bcast(inquire_dimension_length)
182
183  END FUNCTION inquire_dimension_length
184
185  SUBROUTINE Get_Field_r1(field_name,field,found,timeindex)
186  ! For a surface field
187  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
188  IMPLICIT NONE
189    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
190    REAL,INTENT(INOUT)               :: Field(:)
191    LOGICAL,INTENT(OUT),OPTIONAL   :: found
192    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
193
194    integer :: edges(4), corners(4)
195
196    edges(1)=klon_glo
197    edges(2:4)=1
198    corners(1:4)=1
199    if (PRESENT(timeindex)) then
200      corners(2)=timeindex
201    endif
202
203    IF (PRESENT(found)) THEN
204      CALL Get_field_rgen(field_name,field,1,corners,edges,found)
205    ELSE
206      CALL Get_field_rgen(field_name,field,1,corners,edges)
207    ENDIF
208     
209  END SUBROUTINE Get_Field_r1
210 
211  SUBROUTINE Get_Field_r2(field_name,field,found,timeindex)
212  ! For a "3D" horizontal-vertical field
213  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
214  IMPLICIT NONE
215    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
216    REAL,INTENT(INOUT)               :: Field(:,:)
217    LOGICAL,INTENT(OUT),OPTIONAL   :: found
218    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
219
220    integer :: edges(4), corners(4)
221
222    edges(1)=klon_glo
223    edges(2)=size(field,2)
224    edges(3:4)=1
225    corners(1:4)=1
226    if (PRESENT(timeindex)) then
227      corners(3)=timeindex
228    endif
229   
230    IF (PRESENT(found)) THEN
231      CALL Get_field_rgen(field_name,field,size(field,2),&
232                          corners,edges,found)
233    ELSE
234      CALL Get_field_rgen(field_name,field,size(field,2),&
235                          corners,edges)
236    ENDIF
237
238     
239  END SUBROUTINE Get_Field_r2
240 
241  SUBROUTINE Get_Field_r3(field_name,field,found,timeindex)
242  ! for a "4D" field surf/alt/??
243  use mod_grid_phy_lmdz, only: klon_glo ! number of atmospheric columns (full grid)
244  IMPLICIT NONE
245    CHARACTER(LEN=*),INTENT(IN)    :: Field_name
246    REAL,INTENT(INOUT)               :: Field(:,:,:)
247    LOGICAL,INTENT(OUT),OPTIONAL   :: found
248    INTEGER,INTENT(IN),OPTIONAL    :: timeindex ! time index of sought data
249
250    integer :: edges(4), corners(4)
251
252    edges(1)=klon_glo
253    edges(2)=size(field,2)
254    edges(3)=size(field,3)
255    edges(4)=1
256    corners(1:4)=1
257    if (PRESENT(timeindex)) then
258      corners(4)=timeindex
259    endif
260   
261    IF (PRESENT(found)) THEN
262      CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),&
263                          corners,edges,found)
264    ELSE
265      CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),&
266                          corners,edges)
267    ENDIF
268     
269  END SUBROUTINE Get_Field_r3
270 
271  SUBROUTINE Get_field_rgen(field_name,field,field_size, &
272                            corners,edges,found)
273  USE netcdf
274  USE dimphy
275  USE mod_grid_phy_lmdz
276  USE mod_phys_lmdz_para
277  IMPLICIT NONE
278    CHARACTER(LEN=*) :: Field_name
279    INTEGER          :: field_size
280    REAL             :: field(klon,field_size)
281    INTEGER,INTENT(IN) :: corners(4)
282    INTEGER,INTENT(IN) :: edges(4)
283    LOGICAL,OPTIONAL :: found
284   
285    REAL    :: field_glo(klon_glo,field_size)
286    LOGICAL :: tmp_found
287    INTEGER :: varid
288    INTEGER :: ierr
289   
290    IF (is_master) THEN
291 
292      ierr=NF90_INQ_VARID(nid_start,Field_name,varid)
293     
294      IF (ierr==NF90_NOERR) THEN
295        CALL body(field_glo)
296        tmp_found=.TRUE.
297      ELSE
298        tmp_found=.FALSE.
299      ENDIF
300   
301    ENDIF
302   
303    CALL bcast(tmp_found)
304
305    IF (tmp_found) THEN
306      CALL scatter(field_glo,field)
307    ENDIF
308   
309    IF (PRESENT(found)) THEN
310      found=tmp_found
311    ELSE
312      IF (.NOT. tmp_found) THEN
313        PRINT*, 'get_field_rgen: Field <'//field_name//'> not found'
314        CALL abort_physic("get_field_rgen","Failed to get field",1)
315      ENDIF
316    ENDIF
317 
318   
319    CONTAINS
320     
321     SUBROUTINE body(field_glo)
322       REAL :: field_glo(klon_glo*field_size)
323         ierr=NF90_GET_VAR(nid_start,varid,field_glo,corners,edges)
324         IF (ierr/=NF90_NOERR) THEN
325           ! La variable exist dans le fichier mais la lecture a echouee.
326           PRINT*, 'get_field_rgen: Failed reading <'//field_name//'>'
327
328!           IF (field_name=='CLWCON' .OR. field_name=='RNEBCON' .OR. field_name=='RATQS') THEN
329!              ! Essaye de lire le variable sur surface uniqument, comme fait avant
330!              field_glo(:)=0.
331!              ierr=NF90_GET_VAR(nid_start,varid,field_glo(1:klon_glo))
332!              IF (ierr/=NF90_NOERR) THEN
333!                 PRINT*, 'phyetat0: Lecture echouee aussi en 2D pour <'//field_name//'>'
334!                 CALL abort
335!              ELSE
336!                 PRINT*, 'phyetat0: La variable <'//field_name//'> lu sur surface seulement'!, selon ancien format, le reste mis a zero'
337!              END IF
338!           ELSE
339              CALL abort_physic("get_field_rgen","Failed to read field",1)
340!           ENDIF
341         ENDIF
342
343     END SUBROUTINE body
344
345  END SUBROUTINE Get_field_rgen
346
347
348  SUBROUTINE get_var_r0(var_name,var,found)
349  ! Get a scalar from input file
350  IMPLICIT NONE 
351    CHARACTER(LEN=*),INTENT(IN)  :: var_name
352    REAL,INTENT(INOUT)             :: var
353    LOGICAL,OPTIONAL,INTENT(OUT) :: found
354
355    REAL                         :: varout(1)
356   
357    IF (PRESENT(found)) THEN
358      CALL Get_var_rgen(var_name,varout,size(varout),found)
359    ELSE
360      CALL Get_var_rgen(var_name,varout,size(varout))
361    ENDIF
362    var=varout(1)
363 
364  END SUBROUTINE get_var_r0
365
366  SUBROUTINE get_var_r1(var_name,var,found)
367  ! Get a vector from input file
368  IMPLICIT NONE 
369    CHARACTER(LEN=*),INTENT(IN)  :: var_name
370    REAL,INTENT(INOUT)             :: var(:)
371    LOGICAL,OPTIONAL,INTENT(OUT) :: found
372   
373    IF (PRESENT(found)) THEN
374      CALL Get_var_rgen(var_name,var,size(var),found)
375    ELSE
376      CALL Get_var_rgen(var_name,var,size(var))
377    ENDIF
378 
379  END SUBROUTINE get_var_r1
380
381  SUBROUTINE get_var_r2(var_name,var,found)
382  ! Get a 2D field from input file
383  IMPLICIT NONE 
384    CHARACTER(LEN=*),INTENT(IN)  :: var_name
385    REAL,INTENT(OUT)             :: var(:,:)
386    LOGICAL,OPTIONAL,INTENT(OUT) :: found
387   
388    IF (PRESENT(found)) THEN
389      CALL Get_var_rgen(var_name,var,size(var),found)
390    ELSE
391      CALL Get_var_rgen(var_name,var,size(var))
392    ENDIF
393 
394  END SUBROUTINE get_var_r2
395
396  SUBROUTINE get_var_r3(var_name,var,found)
397  ! Get a 3D field frominput file
398  IMPLICIT NONE 
399    CHARACTER(LEN=*),INTENT(IN)  :: var_name
400    REAL,INTENT(INOUT)             :: var(:,:,:)
401    LOGICAL,OPTIONAL,INTENT(OUT) :: found
402   
403    IF (PRESENT(found)) THEN
404      CALL Get_var_rgen(var_name,var,size(var),found)
405    ELSE
406      CALL Get_var_rgen(var_name,var,size(var))
407    ENDIF
408 
409  END SUBROUTINE get_var_r3
410
411  SUBROUTINE Get_var_rgen(var_name,var,var_size,found)
412  USE netcdf
413  USE dimphy
414  USE mod_grid_phy_lmdz
415  USE mod_phys_lmdz_para
416  IMPLICIT NONE
417    CHARACTER(LEN=*) :: var_name
418    INTEGER          :: var_size
419    REAL             :: var(var_size)
420    LOGICAL,OPTIONAL :: found
421   
422    LOGICAL :: tmp_found
423    INTEGER :: varid
424    INTEGER :: ierr
425   
426    IF (is_mpi_root .AND. is_omp_root) 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
672  USE dimphy
673  USE comsoil_h, only: nsoilmx
674  USE comslope_mod, ONLY: nslope
675  USE mod_grid_phy_lmdz
676  USE mod_phys_lmdz_para
677  IMPLICIT NONE
678  CHARACTER(LEN=*),INTENT(IN)    :: field_name
679  CHARACTER(LEN=*),INTENT(IN)    :: title
680  INTEGER,INTENT(IN)             :: field_size
681  REAL,INTENT(IN)                :: field(klon,field_size)
682  REAL,OPTIONAL,INTENT(IN)       :: time
683 
684  REAL                           :: field_glo(klon_glo,field_size)
685  REAL                           :: field_glo_reshape(klon_glo,nsoilmx,nslope,timeindex)
686  INTEGER                        :: ierr
687  INTEGER                        :: nvarid
688  INTEGER                        :: idim
689   
690    CALL gather(field,field_glo)
691   
692    IF (is_master) THEN
693
694      IF (field_size==1) THEN
695        ! input is a 1D "surface field" array
696        if (.not.present(time)) then ! for a time-independent field
697          ierr=NF90_REDEF(nid_restart)
698#ifdef NC_DOUBLE
699          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
700                            (/idim2/),nvarid)
701#else
702          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
703                            (/idim2/),nvarid)
704#endif
705          if (ierr.ne.NF90_NOERR) then
706            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
707            write(*,*)trim(nf90_strerror(ierr))
708          endif
709          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
710          ierr=NF90_ENDDEF(nid_restart)
711          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
712        else
713          ! check if the variable has already been defined:
714          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
715          if (ierr/=NF90_NOERR) then ! variable not found, define it
716            ierr=NF90_REDEF(nid_restart)
717#ifdef NC_DOUBLE
718            if(field_name.eq. "tauscaling" .or. field_name.eq. "totcloudfrac" .or. field_name.eq. "wstar" ) then
719            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
720                              (/idim2,idim7/),nvarid)
721            else
722            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
723                              (/idim2,idim8,idim7/),nvarid)
724            endif
725#else
726            if(field_name.eq. "tauscaling" .or. field_name.eq. "totcloudfrac" .or. field_name.eq. "wstar" ) then
727            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
728                              (/idim2,idim7/),nvarid)
729            else
730            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
731                              (/idim2,idim8,idim7/),nvarid)
732            endif
733#endif
734            if (ierr.ne.NF90_NOERR) then
735              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
736              write(*,*)trim(nf90_strerror(ierr))
737            endif
738            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
739            ierr=NF90_ENDDEF(nid_restart)
740          endif
741          ! Write the variable
742          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
743                            start=(/1,timeindex/))
744        endif ! of if (.not.present(timeindex))
745
746      ELSE IF (field_size==klev) THEN
747        ! input is a 2D "atmospheric field" array
748        if (.not.present(time)) then ! for a time-independent field
749          ierr=NF90_REDEF(nid_restart)
750#ifdef NC_DOUBLE
751          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
752                            (/idim2,idim6/),nvarid)
753#else
754          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
755                            (/idim2,idim6/),nvarid)
756#endif
757          if (ierr.ne.NF90_NOERR) then
758            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
759            write(*,*)trim(nf90_strerror(ierr))
760          endif
761          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
762          ierr=NF90_ENDDEF(nid_restart)
763          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
764        else
765          ! check if the variable has already been defined:
766          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
767          if (ierr/=NF90_NOERR) then ! variable not found, define it
768            ierr=NF90_REDEF(nid_restart)
769#ifdef NC_DOUBLE
770            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
771                              (/idim2,idim6,idim7/),nvarid)
772#else
773            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
774                              (/idim2,idim6,idim7/),nvarid)
775#endif
776            if (ierr.ne.NF90_NOERR) then
777              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
778              write(*,*)trim(nf90_strerror(ierr))
779            endif
780            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
781            ierr=NF90_ENDDEF(nid_restart)
782          endif
783          ! Write the variable
784          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
785                            start=(/1,1,timeindex/))
786        endif ! of if (.not.present(time))
787
788      ELSE IF (field_size==klevp1) THEN
789        ! input is a 2D "interlayer atmospheric field" array
790        if (.not.present(time)) then ! for a time-independent field
791          ierr=NF90_REDEF(nid_restart)
792#ifdef NC_DOUBLE
793          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
794                            (/idim2,idim4/),nvarid)
795#else
796          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
797                            (/idim2,idim4/),nvarid)
798#endif
799          if (ierr.ne.NF90_NOERR) then
800            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
801            write(*,*)trim(nf90_strerror(ierr))
802          endif
803          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
804          ierr = NF90_ENDDEF(nid_restart)
805          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
806        else
807          ! check if the variable has already been defined:
808          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
809          if (ierr/=NF90_NOERR) then ! variable not found, define it
810            ierr=NF90_REDEF(nid_restart)
811#ifdef NC_DOUBLE
812            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
813                              (/idim2,idim4,idim7/),nvarid)
814#else
815            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
816                              (/idim2,idim4,idim7/),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          endif
825          ! Write the variable
826          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
827                            start=(/1,1,timeindex/))
828        endif ! of if (.not.present(timeindex))
829
830      ELSE IF (field_size==nsoilmx .or. field_size==nsoilmx*nslope) THEN
831        ! input is a 2D "subsurface field" array
832        if (.not.present(time)) then ! for a time-independent field
833          ierr = NF90_REDEF(nid_restart)
834#ifdef NC_DOUBLE
835          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
836                            (/idim2,idim3/),nvarid)
837#else
838          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
839                            (/idim2,idim3/),nvarid)
840#endif
841          if (ierr.ne.NF90_NOERR) then
842            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
843            write(*,*)trim(nf90_strerror(ierr))
844          endif
845          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
846          ierr = NF90_ENDDEF(nid_restart)
847          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
848        else
849          ! check if the variable has already been defined:
850          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
851          if (ierr/=NF90_NOERR) then ! variable not found, define it
852            ierr=NF90_REDEF(nid_restart)
853#ifdef NC_DOUBLE
854            if(field_name.eq. "tsoil".or. field_name.eq. "inertiesoil") then
855            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
856                              (/idim2,idim3,idim8,idim7/),nvarid)
857            else
858            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
859                              (/idim2,idim3,idim7/),nvarid)
860            endif
861#else
862            if(field_name.eq. "tsoil".or. field_name.eq. "inertiesoil") then
863            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
864                              (/idim2,idim3,idim8,idim7/),nvarid)
865            else
866            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
867                              (/idim2,idim3,idim7/),nvarid)
868            endif
869#endif
870           if (ierr.ne.NF90_NOERR) then
871              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
872              write(*,*)trim(nf90_strerror(ierr))
873            endif
874            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
875            ierr=NF90_ENDDEF(nid_restart)
876          endif
877          ! Write the variable
878          if(field_name.eq. "tsoil" .or. field_name.eq. "inertiesoil") then
879
880          field_glo_reshape=RESHAPE(field_glo, (/klon_glo, nsoilmx, nslope, timeindex/))
881
882          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo_reshape,&
883                            start=(/1,1,1,timeindex/))
884          else
885          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
886                            start=(/1,1,timeindex/))
887          endif
888        endif ! of if (.not.present(time))
889
890      ELSE IF (field_size==nslope) THEN
891        ! input is a 2D "subsurface field" array
892        if (.not.present(time)) then ! for a time-independent field
893          ierr = NF90_REDEF(nid_restart)
894#ifdef NC_DOUBLE
895          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
896                            (/idim2,idim8/),nvarid)
897#else
898          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
899                            (/idim2,idim8/),nvarid)
900#endif
901          if (ierr.ne.NF90_NOERR) then
902            write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
903            write(*,*)trim(nf90_strerror(ierr))
904          endif
905          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
906          ierr = NF90_ENDDEF(nid_restart)
907          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
908        else
909          ! check if the variable has already been defined:
910          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
911          if (ierr/=NF90_NOERR) then ! variable not found, define it
912            ierr=NF90_REDEF(nid_restart)
913#ifdef NC_DOUBLE
914            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
915                              (/idim2,idim8,idim7/),nvarid)
916#else
917            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
918                              (/idim2,idim8,idim7/),nvarid)
919#endif
920           if (ierr.ne.NF90_NOERR) then
921              write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
922              write(*,*)trim(nf90_strerror(ierr))
923            endif
924            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
925            ierr=NF90_ENDDEF(nid_restart)
926          endif
927          ! Write the variable
928          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
929                            start=(/1,1,timeindex/))
930
931        endif ! of if (.not.present(time))
932
933      ELSE
934        PRINT *, "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name)
935        write(*,*) "  field_size =",field_size
936        CALL abort_physic("put_field_rgen","wrong field dimension",1)
937      ENDIF
938
939      ! Check the writting of field to file went OK
940      if (ierr.ne.NF90_NOERR) then
941        write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name)
942        write(*,*)trim(nf90_strerror(ierr))
943        CALL abort_physic("put_field_rgen","Failed writing field",1)
944      endif
945
946    ENDIF ! of IF (is_master)
947   
948  END SUBROUTINE put_field_rgen 
949 
950  SUBROUTINE put_var_r0(var_name,title,var)
951  ! Put a scalar in file
952   IMPLICIT NONE
953     CHARACTER(LEN=*),INTENT(IN) :: var_name
954     CHARACTER(LEN=*),INTENT(IN) :: title
955     REAL,INTENT(IN)             :: var
956     REAL                        :: varin(1)
957     
958     varin(1)=var
959     
960     CALL put_var_rgen(var_name,title,varin,size(varin))
961
962  END SUBROUTINE put_var_r0
963
964  SUBROUTINE put_var_r1(var_name,title,var)
965  ! Put a vector in file
966   IMPLICIT NONE
967     CHARACTER(LEN=*),INTENT(IN) :: var_name
968     CHARACTER(LEN=*),INTENT(IN) :: title
969     REAL,INTENT(IN)             :: var(:)
970     
971     CALL put_var_rgen(var_name,title,var,size(var))
972
973  END SUBROUTINE put_var_r1
974 
975  SUBROUTINE put_var_r2(var_name,title,var)
976  ! Put a 2D field in file
977   IMPLICIT NONE
978     CHARACTER(LEN=*),INTENT(IN) :: var_name
979     CHARACTER(LEN=*),INTENT(IN) :: title
980     REAL,INTENT(IN)             :: var(:,:)
981     
982     CALL put_var_rgen(var_name,title,var,size(var))
983
984  END SUBROUTINE put_var_r2
985 
986  SUBROUTINE put_var_r3(var_name,title,var)
987  ! Put a 3D field in file
988   IMPLICIT NONE
989     CHARACTER(LEN=*),INTENT(IN) :: var_name
990     CHARACTER(LEN=*),INTENT(IN) :: title
991     REAL,INTENT(IN)             :: var(:,:,:)
992     
993     CALL put_var_rgen(var_name,title,var,size(var))
994
995  END SUBROUTINE put_var_r3
996
997  SUBROUTINE put_var_rgen(var_name,title,var,var_size)
998  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
999                    NF90_FLOAT, NF90_DOUBLE, &
1000                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
1001                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
1002  USE comsoil_h, only: nsoilmx
1003  USE comslope_mod, only: nslope
1004  USE mod_phys_lmdz_para, only: is_master
1005  IMPLICIT NONE
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.
1015         
1016    IF (is_master) THEN
1017
1018      IF (var_name=="Time") THEN
1019        ! Very specific case of "Time" variable
1020        if (firsttime .or. first_time1d) then
1021          ! Create the "Time variable"
1022          ierr=NF90_REDEF(nid_restart)
1023#ifdef NC_DOUBLE
1024          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,&
1025                            (/idim7/),nvarid)
1026#else
1027          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,&
1028                            (/idim7/),nvarid)
1029#endif
1030          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1031          ierr=NF90_ENDDEF(nid_restart)
1032          firsttime=.false.
1033          if(first_time1d) first_time1d=.false.
1034        endif
1035        ! Append "time" value
1036        ! get current length of "Time" dimension
1037        ierr=nf90_inq_dimid(nid_restart,var_name,idim1d)
1038        ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex)
1039        timeindex=timeindex+1
1040        ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid)
1041        ierr=NF90_PUT_VAR(nid_restart,nvarid,var,&
1042                           start=(/timeindex/))
1043        IF (ierr/=NF90_NOERR) THEN
1044          write(*,*)'put_var_rgen: problem writing Time'
1045          write(*,*)trim(nf90_strerror(ierr))
1046          CALL abort_physic("put_var_rgen","Failed to write Time",1)
1047        ENDIF
1048        return ! nothing left to do
1049      ELSEIF (var_size==length) THEN
1050        ! We know it is a "controle" kind of 1D array
1051        idim1d=idim1
1052      ELSEIF (var_size==nsoilmx) THEN
1053        ! We know it is an "mlayer" kind of 1D array
1054        idim1d=idim3
1055      ELSEIF (var_size==nslope+1) THEN
1056        ! We know it is an "inter slope" kind of 1D array
1057        idim1d=idim9
1058      ELSE
1059        PRINT *, "put_var_rgen error : wrong dimension"
1060        write(*,*) "  var_size =",var_size
1061        CALL abort_physic("put_var_rgen","Wrong variable dimension",1)
1062
1063      ENDIF ! of IF (var_size==length) THEN
1064
1065      ! Swich to NetCDF define mode
1066      ierr=NF90_REDEF (nid_restart)
1067      ! Define the variable
1068#ifdef NC_DOUBLE
1069      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid)
1070#else
1071      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid)
1072#endif
1073     
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))
1083        CALL abort_physic("put_var_rgen","Failed writing variable",1)
1084      ENDIF
1085    ENDIF ! of IF (is_master)
1086   
1087  END SUBROUTINE put_var_rgen
1088
1089  SUBROUTINE put_var_c1(var_name,title,var)
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 comsoil_h, only: nsoilmx
1097  USE comslope_mod, only: nslope
1098  USE mod_phys_lmdz_para, only: is_master
1099
1100   IMPLICIT NONE
1101     CHARACTER(LEN=*),INTENT(IN) :: var_name
1102     CHARACTER(LEN=*),INTENT(IN) :: title
1103     CHARACTER(LEN=*),INTENT(IN) :: var(:)
1104
1105     INTEGER :: ierr
1106     INTEGER :: nvarid
1107     INTEGER :: idim1d_1, idim1d_2
1108     INTEGER :: var_size
1109     logical,save :: firsttime=.true.
1110
1111    IF (is_master) THEN
1112
1113      var_size = size(var)
1114      IF (var_size==ldscrpt) THEN
1115        ! We know it is a "controle descriptor" kind of 1D array
1116        idim1d_1=idim11
1117        idim1d_2=idim10
1118      ELSE
1119        PRINT *, "put_var_cgen error : wrong dimension"
1120        write(*,*) "  var_size =",var_size
1121        CALL abort_physic("put_var_cgen","Wrong variable dimension",1)
1122
1123      ENDIF ! of IF (var_size==length) THEN
1124
1125      ! Swich to NetCDF define mode
1126      ierr=NF90_REDEF (nid_restart)
1127      ! Define the variable
1128      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_CHAR,(/idim1d_1,idim1d_2/),nvarid)
1129      ! Add a "title" attribute
1130      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1131      ! Swich out of define mode
1132      ierr=NF90_ENDDEF(nid_restart)
1133      ! Write variable to file
1134      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
1135      IF (ierr/=NF90_NOERR) THEN
1136        write(*,*)'put_var_cgen: problem writing '//trim(var_name)
1137        write(*,*)trim(nf90_strerror(ierr))
1138        CALL abort_physic("put_var_cgen","Failed writing variable",1)
1139      ENDIF
1140    ENDIF ! of IF (is_master)
1141
1142  END SUBROUTINE put_var_c1
1143
1144END MODULE iostart
Note: See TracBrowser for help on using the repository browser.