source: trunk/LMDZ.COMMON/libf/evolution/iostart_PEM.F90 @ 3495

Last change on this file since 3495 was 3457, checked in by jbclement, 6 weeks ago

PCM:
Modification of dimension detection for the variables written in "diagpem.nc": in particular for 'nb_str_max' which can evolve and match the value of other dimensions.
JBC

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