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

Last change on this file since 2850 was 2842, checked in by romain.vande, 3 years ago

Mars PEM:
PEM is now adapted to run with XIOS diurnal averages (when they will work properly)
RV

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