source: trunk/LMDZ.MARS/libf/phymars/iostart.F90 @ 1242

Last change on this file since 1242 was 1130, checked in by emillour, 11 years ago

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

File size: 31.3 KB
Line 
1MODULE iostart
2
3    IMPLICIT NONE
4PRIVATE
5    INTEGER,SAVE :: nid_start ! NetCDF file identifier for startfi.nc file
6    INTEGER,SAVE :: nid_restart ! NetCDF file identifier for restartfi.nc file
7   
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      ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",nqtot,idim5)
531      IF (ierr/=NF90_NOERR) THEN
532        write(*,*)'open_restartphy: problem defining number_of_advected_fields dimension '
533        write(*,*)trim(nf90_strerror(ierr))
534        CALL ABORT
535      ENDIF
536
537      ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
538      IF (ierr/=NF90_NOERR) THEN
539        write(*,*)'open_restartphy: problem defining nlayer dimension '
540        write(*,*)trim(nf90_strerror(ierr))
541        CALL ABORT
542      ENDIF
543     
544      ierr=NF90_DEF_DIM(nid_restart,"Time",NF90_UNLIMITED,idim7)
545      IF (ierr/=NF90_NOERR) THEN
546        write(*,*)'open_restartphy: problem defining Time dimension '
547        write(*,*)trim(nf90_strerror(ierr))
548        CALL ABORT
549      ENDIF
550
551      ierr=NF90_ENDDEF(nid_restart)
552      IF (ierr/=NF90_NOERR) THEN
553        write(*,*)'open_restartphy: problem ending definition mode '
554        write(*,*)trim(nf90_strerror(ierr))
555        CALL ABORT
556      ENDIF
557    ENDIF
558
559  END SUBROUTINE open_restartphy
560
561  SUBROUTINE close_restartphy
562  USE netcdf, only: NF90_CLOSE
563  USE mod_phys_lmdz_para, only: is_master
564  IMPLICIT NONE
565    INTEGER          :: ierr
566
567    IF (is_master) THEN
568      ierr = NF90_CLOSE (nid_restart)
569    ENDIF
570 
571  END SUBROUTINE close_restartphy
572
573  SUBROUTINE put_field_r1(field_name,title,field,time)
574  ! For a surface field
575  IMPLICIT NONE
576  CHARACTER(LEN=*),INTENT(IN)    :: field_name
577  CHARACTER(LEN=*),INTENT(IN)    :: title
578  REAL,INTENT(IN)                :: field(:)
579  REAL,OPTIONAL,INTENT(IN)       :: time
580 
581  IF (present(time)) THEN
582    ! if timeindex is present, it is a time-dependent variable
583    CALL put_field_rgen(field_name,title,field,1,time)
584  ELSE
585    CALL put_field_rgen(field_name,title,field,1)
586  ENDIF
587 
588  END SUBROUTINE put_field_r1
589
590  SUBROUTINE put_field_r2(field_name,title,field,time)
591  ! For a "3D" horizontal-vertical field
592  IMPLICIT NONE
593  CHARACTER(LEN=*),INTENT(IN)    :: field_name
594  CHARACTER(LEN=*),INTENT(IN)    :: title
595  REAL,INTENT(IN)                :: field(:,:)
596  REAL,OPTIONAL,INTENT(IN)       :: time
597 
598  IF (present(time)) THEN
599    ! if timeindex is present, it is a time-dependent variable
600    CALL put_field_rgen(field_name,title,field,size(field,2),time)
601  ELSE
602    CALL put_field_rgen(field_name,title,field,size(field,2))
603  ENDIF
604 
605  END SUBROUTINE put_field_r2
606
607  SUBROUTINE put_field_r3(field_name,title,field,time)
608  ! For a "4D" field surf/alt/??
609  IMPLICIT NONE
610  CHARACTER(LEN=*),INTENT(IN)    :: field_name
611  CHARACTER(LEN=*),INTENT(IN)    :: title
612  REAL,INTENT(IN)                :: field(:,:,:)
613  REAL,OPTIONAL,INTENT(IN)       :: time
614 
615  IF (present(time)) THEN
616    ! if timeindex is present, it is a time-dependent variable
617    CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3),&
618                        time)
619  ELSE 
620    CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3))
621  ENDIF
622 
623  END SUBROUTINE put_field_r3
624 
625  SUBROUTINE put_field_rgen(field_name,title,field,field_size,time)
626  USE netcdf
627  USE dimphy
628  USE comsoil_h, only: nsoilmx
629  USE mod_grid_phy_lmdz
630  USE mod_phys_lmdz_para
631  IMPLICIT NONE
632  CHARACTER(LEN=*),INTENT(IN)    :: field_name
633  CHARACTER(LEN=*),INTENT(IN)    :: title
634  INTEGER,INTENT(IN)             :: field_size
635  REAL,INTENT(IN)                :: field(klon,field_size)
636  REAL,OPTIONAL,INTENT(IN)       :: time
637 
638  REAL                           :: field_glo(klon_glo,field_size)
639  INTEGER                        :: ierr
640  INTEGER                        :: nvarid
641  INTEGER                        :: idim
642   
643    CALL gather(field,field_glo)
644   
645    IF (is_master) THEN
646
647      IF (field_size==1) THEN
648        ! input is a 1D "surface field" array
649        if (.not.present(time)) then ! for a time-independent field
650          ierr=NF90_REDEF(nid_restart)
651#ifdef NC_DOUBLE
652          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
653                            (/idim2/),nvarid)
654#else
655          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
656                            (/idim2/),nvarid)
657#endif
658          if (ierr.ne.NF90_NOERR) then
659            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
660            write(*,*)trim(nf90_strerror(ierr))
661          endif
662          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
663          ierr=NF90_ENDDEF(nid_restart)
664          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
665        else
666          ! check if the variable has already been defined:
667          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
668          if (ierr/=NF90_NOERR) then ! variable not found, define it
669            ierr=NF90_REDEF(nid_restart)
670#ifdef NC_DOUBLE
671            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
672                              (/idim2,idim7/),nvarid)
673#else
674            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
675                              (/idim2,idim7/),nvarid)
676#endif
677            if (ierr.ne.NF90_NOERR) then
678              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
679              write(*,*)trim(nf90_strerror(ierr))
680            endif
681            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
682            ierr=NF90_ENDDEF(nid_restart)
683          endif
684          ! Write the variable
685          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
686                            start=(/1,timeindex/))
687        endif ! of if (.not.present(timeindex))
688
689      ELSE IF (field_size==klev) THEN
690        ! input is a 2D "atmospheric field" array
691        if (.not.present(time)) then ! for a time-independent field
692          ierr=NF90_REDEF(nid_restart)
693#ifdef NC_DOUBLE
694          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
695                            (/idim2,idim6/),nvarid)
696#else
697          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
698                            (/idim2,idim6/),nvarid)
699#endif
700          if (ierr.ne.NF90_NOERR) then
701            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
702            write(*,*)trim(nf90_strerror(ierr))
703          endif
704          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
705          ierr=NF90_ENDDEF(nid_restart)
706          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
707        else
708          ! check if the variable has already been defined:
709          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
710          if (ierr/=NF90_NOERR) then ! variable not found, define it
711            ierr=NF90_REDEF(nid_restart)
712#ifdef NC_DOUBLE
713            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
714                              (/idim2,idim6,idim7/),nvarid)
715#else
716            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
717                              (/idim2,idim6,idim7/),nvarid)
718#endif
719            if (ierr.ne.NF90_NOERR) then
720              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
721              write(*,*)trim(nf90_strerror(ierr))
722            endif
723            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
724            ierr=NF90_ENDDEF(nid_restart)
725          endif
726          ! Write the variable
727          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
728                            start=(/1,1,timeindex/))
729        endif ! of if (.not.present(time))
730
731      ELSE IF (field_size==klevp1) THEN
732        ! input is a 2D "interlayer atmospheric field" array
733        if (.not.present(time)) then ! for a time-independent field
734          ierr=NF90_REDEF(nid_restart)
735#ifdef NC_DOUBLE
736          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
737                            (/idim2,idim4/),nvarid)
738#else
739          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
740                            (/idim2,idim4/),nvarid)
741#endif
742          if (ierr.ne.NF90_NOERR) then
743            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
744            write(*,*)trim(nf90_strerror(ierr))
745          endif
746          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
747          ierr = NF90_ENDDEF(nid_restart)
748          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
749        else
750          ! check if the variable has already been defined:
751          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
752          if (ierr/=NF90_NOERR) then ! variable not found, define it
753            ierr=NF90_REDEF(nid_restart)
754#ifdef NC_DOUBLE
755            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
756                              (/idim2,idim4,idim7/),nvarid)
757#else
758            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
759                              (/idim2,idim4,idim7/),nvarid)
760#endif
761            if (ierr.ne.NF90_NOERR) then
762              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
763              write(*,*)trim(nf90_strerror(ierr))
764            endif
765            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
766            ierr=NF90_ENDDEF(nid_restart)
767          endif
768          ! Write the variable
769          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
770                            start=(/1,1,timeindex/))
771        endif ! of if (.not.present(timeindex))
772
773      ELSE IF (field_size==nsoilmx) THEN
774        ! input is a 2D "subsurface field" array
775        if (.not.present(time)) then ! for a time-independent field
776          ierr = NF90_REDEF(nid_restart)
777#ifdef NC_DOUBLE
778          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
779                            (/idim2,idim3/),nvarid)
780#else
781          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
782                            (/idim2,idim3/),nvarid)
783#endif
784          if (ierr.ne.NF90_NOERR) then
785            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
786            write(*,*)trim(nf90_strerror(ierr))
787          endif
788          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
789          ierr = NF90_ENDDEF(nid_restart)
790          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
791        else
792          ! check if the variable has already been defined:
793          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
794          if (ierr/=NF90_NOERR) then ! variable not found, define it
795            ierr=NF90_REDEF(nid_restart)
796#ifdef NC_DOUBLE
797            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
798                              (/idim2,idim3,idim7/),nvarid)
799#else
800            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
801                              (/idim2,idim3,idim7/),nvarid)
802#endif
803           if (ierr.ne.NF90_NOERR) then
804              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
805              write(*,*)trim(nf90_strerror(ierr))
806            endif
807            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
808            ierr=NF90_ENDDEF(nid_restart)
809          endif
810          ! Write the variable
811          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
812                            start=(/1,1,timeindex/))
813
814        endif ! of if (.not.present(time))
815
816      ELSE
817        PRINT *, "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name)
818        write(*,*) "  field_size =",field_size
819        CALL ABORT
820      ENDIF
821
822      ! Check the writting of field to file went OK
823      if (ierr.ne.NF90_NOERR) then
824        write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name)
825        write(*,*)trim(nf90_strerror(ierr))
826        call abort
827      endif
828
829    ENDIF ! of IF (is_master)
830   
831  END SUBROUTINE put_field_rgen 
832 
833  SUBROUTINE put_var_r0(var_name,title,var)
834  ! Put a scalar in file
835   IMPLICIT NONE
836     CHARACTER(LEN=*),INTENT(IN) :: var_name
837     CHARACTER(LEN=*),INTENT(IN) :: title
838     REAL,INTENT(IN)             :: var
839     REAL                        :: varin(1)
840     
841     varin(1)=var
842     
843     CALL put_var_rgen(var_name,title,varin,size(varin))
844
845  END SUBROUTINE put_var_r0
846
847
848  SUBROUTINE put_var_r1(var_name,title,var)
849  ! Put a vector in file
850   IMPLICIT NONE
851     CHARACTER(LEN=*),INTENT(IN) :: var_name
852     CHARACTER(LEN=*),INTENT(IN) :: title
853     REAL,INTENT(IN)             :: var(:)
854     
855     CALL put_var_rgen(var_name,title,var,size(var))
856
857  END SUBROUTINE put_var_r1
858 
859  SUBROUTINE put_var_r2(var_name,title,var)
860  ! Put a 2D field in file
861   IMPLICIT NONE
862     CHARACTER(LEN=*),INTENT(IN) :: var_name
863     CHARACTER(LEN=*),INTENT(IN) :: title
864     REAL,INTENT(IN)             :: var(:,:)
865     
866     CALL put_var_rgen(var_name,title,var,size(var))
867
868  END SUBROUTINE put_var_r2     
869 
870  SUBROUTINE put_var_r3(var_name,title,var)
871  ! Put a 3D field in file
872   IMPLICIT NONE
873     CHARACTER(LEN=*),INTENT(IN) :: var_name
874     CHARACTER(LEN=*),INTENT(IN) :: title
875     REAL,INTENT(IN)             :: var(:,:,:)
876     
877     CALL put_var_rgen(var_name,title,var,size(var))
878
879  END SUBROUTINE put_var_r3
880
881  SUBROUTINE put_var_rgen(var_name,title,var,var_size)
882  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
883                    NF90_FLOAT, NF90_DOUBLE, NF90_DEF_VAR, &
884                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
885                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
886  USE comsoil_h, only: nsoilmx
887  USE mod_phys_lmdz_para, only: is_master
888  IMPLICIT NONE
889     CHARACTER(LEN=*),INTENT(IN) :: var_name
890     CHARACTER(LEN=*),INTENT(IN) :: title
891     INTEGER,INTENT(IN)          :: var_size
892     REAL,INTENT(IN)             :: var(var_size)
893     
894     INTEGER :: ierr
895     INTEGER :: nvarid
896     INTEGER :: idim1d
897     logical,save :: firsttime=.true.
898         
899    IF (is_master) THEN
900
901      IF (var_name=="Time") THEN
902        ! Very specific case of "Time" variable
903        if (firsttime) then
904          ! Create the "Time variable"
905          ierr=NF90_REDEF(nid_restart)
906#ifdef NC_DOUBLE
907          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,&
908                            (/idim7/),nvarid)
909#else
910          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,&
911                            (/idim7/),nvarid)
912#endif
913          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
914          ierr=NF90_ENDDEF(nid_restart)
915         
916          firsttime=.false.
917        endif
918        ! Append "time" value
919        ! get current length of "Time" dimension
920        ierr=nf90_inq_dimid(nid_restart,var_name,idim1d)
921        ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex)
922        timeindex=timeindex+1
923        ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid)
924        ierr=NF90_PUT_VAR(nid_restart,nvarid,var,&
925                           start=(/timeindex/))
926        IF (ierr/=NF90_NOERR) THEN
927          write(*,*)'put_var_rgen: problem writing Time'
928          write(*,*)trim(nf90_strerror(ierr))
929          CALL ABORT
930        ENDIF
931        return ! nothing left to do
932      ELSEIF (var_size==length) THEN
933        ! We know it is a "controle" kind of 1D array
934        idim1d=idim1
935      ELSEIF (var_size==nsoilmx) THEN
936        ! We know it is an  "mlayer" kind of 1D array
937        idim1d=idim3
938      ELSE
939        PRINT *, "put_var_rgen error : wrong dimension"
940        write(*,*) "  var_size =",var_size
941        CALL abort
942
943      ENDIF ! of IF (var_size==length) THEN
944
945      ! Swich to NetCDF define mode
946      ierr=NF90_REDEF (nid_restart)
947      ! Define the variable
948#ifdef NC_DOUBLE
949      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid)
950#else
951      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid)
952#endif
953      ! Add a "title" attribute
954      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
955      ! Swich out of define mode
956      ierr=NF90_ENDDEF(nid_restart)
957      ! Write variable to file
958      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
959      IF (ierr/=NF90_NOERR) THEN
960        write(*,*)'put_var_rgen: problem writing '//trim(var_name)
961        write(*,*)trim(nf90_strerror(ierr))
962        CALL ABORT
963      ENDIF
964    ENDIF ! of IF (is_master)
965   
966  END SUBROUTINE put_var_rgen     
967
968END MODULE iostart
Note: See TracBrowser for help on using the repository browser.