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

Last change on this file since 1255 was 1217, checked in by milmd, 11 years ago

Bug fix for when there is zero tracer.
MI

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