source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/iostart.F90 @ 1057

Last change on this file since 1057 was 313, checked in by ymipsl, 10 years ago
  • implement splitting of XIOS file for lmdz physics
  • Termination is done properly in parallel by calling MPI_ABORT instead of abort or stop

YM

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