source: trunk/LMDZ.TITAN/libf/phytitan/iostart.F90 @ 1871

Last change on this file since 1871 was 1871, checked in by jvatant, 7 years ago

Making chemistry handling more flexible - Step 1
+ Enable upper_chemistry_dimension in the startfi files
+ Added a comchem_h for all the chemistry related stuff
+ start2archive adapted to these modifs : next step newstart !
--JVO

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