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

Last change on this file since 3523 was 3515, checked in by emillour, 11 days ago

Generic PCM:
Add a controle_descriptor to restartfi.nc file, describing contents of the
controle array stored therein. Remove obsolete unused variables "lmixmin"
and "emin_turb" in the process.
EM

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