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

Last change on this file since 1303 was 1297, checked in by bclmd, 10 years ago

LMDZ.GENERIC: slab ocean added

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