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

Last change on this file since 3985 was 3983, checked in by jbclement, 6 days ago

PEM:

  • Removing completely the ice metamorphism computed by a threshold at the end of the PCM (which was commented).
  • Addition of a module "metamorphism" to compute the PCM frost at the PEM beginning and give it back to the PCM at the PEM end. The frost is considered as the ice given by the PCM "startfi.nc" which is above the yearly minimum. Thereby, metamorphism is performed through this operation.
  • Ice reservoirs representation in the PEM is modernized.

JBC

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