source: LMDZ5/trunk/libf/phylmd/phys_output_mod.F90 @ 3003

Last change on this file since 3003 was 3003, checked in by Laurent Fairhead, 7 years ago

Modifications to the code and xml files to output Ap and B, the coefficients
of the hybrid coordinates as requested by the CMIP6 DataRequest?
LF (with guidance from A. Caubel and S. Senesi)

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 25.0 KB
Line 
1! $Id: phys_output_mod.F90 3003 2017-10-03 14:45:26Z fairhead $
2!
3
4MODULE phys_output_mod
5  USE indice_sol_mod
6  USE phys_output_var_mod
7  USE phys_output_write_mod, ONLY : phys_output_write
8  REAL, DIMENSION(nfiles),SAVE :: ecrit_files
9
10! Abderrahmane 12 2007
11!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
12!!! Ecreture des Sorties du modele dans les fichiers Netcdf :
13! histmth.nc : moyennes mensuelles
14! histday.nc : moyennes journalieres
15! histhf.nc  : moyennes toutes les 3 heures
16! histins.nc : valeurs instantanees
17!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
18
19CONTAINS
20
21!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
22!!!!!!!!! Ouverture des fichier et definition des variable de sortie !!!!!!!!
23  !! histbeg, histvert et histdef
24!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
25
26  SUBROUTINE phys_output_open(rlon,rlat,pim,tabij,ipt,jpt,plon,plat, &
27       jjmp1,nlevSTD,clevSTD,rlevSTD, dtime, ok_veget, &
28       type_ocean, iflag_pbl,iflag_pbl_split,ok_mensuel,ok_journe, &
29       ok_hf,ok_instan,ok_LES,ok_ade,ok_aie, read_climoz, &
30       phys_out_filestations, &
31       new_aod, aerosol_couple, flag_aerosol_strat, &
32       pdtphys, paprs, pphis, pplay, lmax_th, ptconv, ptconvth, ivap, &
33       d_u, d_t, qx, d_qx, zmasse, ok_sync)   
34
35    USE iophy
36    USE dimphy
37    USE infotrac_phy, ONLY: nqtot, nqo, niadv, tname, ttext
38    USE ioipsl
39    USE phys_cal_mod, only : hour, calend
40    USE mod_phys_lmdz_para
41    !Martin
42    USE surface_data, ONLY : ok_snow
43    USE phys_output_ctrlout_mod
44    USE mod_grid_phy_lmdz, only: klon_glo,nbp_lon,nbp_lat
45    USE print_control_mod, ONLY: prt_level,lunout
46    USE vertical_layers_mod, ONLY: ap,bp,preff,presnivs, aps, bps, pseudoalt
47    USE time_phylmdz_mod, ONLY: day_ini, itau_phy, start_time, annee_ref, day_ref
48#ifdef CPP_XIOS
49    ! ug Pour les sorties XIOS
50    USE wxios
51#endif
52
53    IMPLICIT NONE
54    include "clesphys.h"
55    include "thermcell.h"
56    include "YOMCST.h"
57
58    ! ug Nouveaux arguments n\'ecessaires au histwrite_mod:
59    INTEGER, INTENT(IN)                         :: ivap
60    INTEGER, DIMENSION(klon), INTENT(IN)        :: lmax_th
61    LOGICAL, INTENT(IN)                         :: ok_sync
62    LOGICAL, DIMENSION(klon, klev), INTENT(IN)  :: ptconv, ptconvth
63    REAL, INTENT(IN)                            :: pdtphys
64    REAL, DIMENSION(klon), INTENT(IN)           :: pphis
65    REAL, DIMENSION(klon, klev), INTENT(IN)     :: pplay, d_u, d_t
66    REAL, DIMENSION(klon, klev+1), INTENT(IN)   :: paprs
67    REAL, DIMENSION(klon,klev,nqtot), INTENT(IN):: qx, d_qx
68    REAL, DIMENSION(klon, klev), INTENT(IN)      :: zmasse
69
70
71    REAL,DIMENSION(klon),INTENT(IN) :: rlon
72    REAL,DIMENSION(klon),INTENT(IN) :: rlat
73    INTEGER, INTENT(IN)             :: pim
74    INTEGER, DIMENSION(pim)            :: tabij
75    INTEGER,DIMENSION(pim), INTENT(IN) :: ipt, jpt
76    REAL,DIMENSION(pim), INTENT(IN) :: plat, plon
77    REAL,DIMENSION(pim,2) :: plat_bounds, plon_bounds
78
79    INTEGER                               :: jjmp1
80    INTEGER                               :: nlevSTD, radpas
81    LOGICAL                               :: ok_mensuel, ok_journe, ok_hf, ok_instan
82    LOGICAL                               :: ok_LES,ok_ade,ok_aie
83    INTEGER                               :: flag_aerosol_strat
84    LOGICAL                               :: new_aod, aerosol_couple
85    INTEGER, INTENT(IN)::  read_climoz ! read ozone climatology
86    !     Allowed values are 0, 1 and 2
87    !     0: do not read an ozone climatology
88    !     1: read a single ozone climatology that will be used day and night
89    !     2: read two ozone climatologies, the average day and night
90    !     climatology and the daylight climatology
91
92    REAL                                  :: dtime
93    INTEGER                               :: idayref
94    REAL                                  :: zjulian_start, zjulian
95    CHARACTER(LEN=4), DIMENSION(nlevSTD)  :: clevSTD
96    REAL, DIMENSION(nlevSTD)              :: rlevSTD
97    INTEGER                               :: nsrf, k, iq, iiq, iff, i, j, ilev
98    INTEGER                               :: naero
99    LOGICAL                               :: ok_veget
100    INTEGER                               :: iflag_pbl
101    INTEGER                               :: iflag_pbl_split
102    CHARACTER(LEN=4)                      :: bb2
103    CHARACTER(LEN=2)                      :: bb3
104    CHARACTER(LEN=6)                      :: type_ocean
105    INTEGER, DIMENSION(nbp_lon*jjmp1)         ::  ndex2d
106    INTEGER, DIMENSION(nbp_lon*jjmp1*klev)    :: ndex3d
107    INTEGER                               :: imin_ins, imax_ins
108    INTEGER                               :: jmin_ins, jmax_ins
109    INTEGER, DIMENSION(nfiles)            :: phys_out_levmin, phys_out_levmax
110    INTEGER, DIMENSION(nfiles)            :: phys_out_filelevels
111    CHARACTER(LEN=20), DIMENSION(nfiles)  :: chtimestep = (/ 'Default', 'Default', 'Default', 'Default', 'Default', &
112                                                             'Default', 'Default', 'Default', 'Default', 'Default' /)
113    LOGICAL, DIMENSION(nfiles)            :: phys_out_filekeys
114    LOGICAL, DIMENSION(nfiles)            :: phys_out_filestations
115
116!!!!!!!!!! stockage dans une region limitee pour chaque fichier !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
117    !                 entre [phys_out_lonmin,phys_out_lonmax] et [phys_out_latmin,phys_out_latmax]
118
119    LOGICAL, DIMENSION(nfiles), SAVE  :: phys_out_regfkey       = (/ .FALSE., .FALSE., .FALSE., .FALSE., .FALSE., &
120                                                                     .FALSE., .FALSE., .FALSE., .FALSE., .FALSE. /)
121    REAL, DIMENSION(nfiles), SAVE     :: phys_out_lonmin        = (/ -180., -180., -180., -180., -180., &
122                                                                     -180., -180., -180., -180., -180. /)
123    REAL, DIMENSION(nfiles), SAVE     :: phys_out_lonmax        = (/  180.,  180.,  180.,  180.,  180., &
124                                                                      180.,  180.,  180.,  180.,  180. /)
125    REAL, DIMENSION(nfiles), SAVE     :: phys_out_latmin        = (/  -90.,  -90.,  -90.,  -90.,  -90., &
126                                                                      -90.,  -90.,  -90.,  -90.,  -90. /)
127    REAL, DIMENSION(nfiles), SAVE     :: phys_out_latmax        = (/   90.,   90.,   90.,   90.,   90., &
128                                                                       90.,   90.,   90.,   90.,   90. /)
129    REAL, DIMENSION(klev,2) :: Ahyb_bounds, Bhyb_bounds
130    REAL, DIMENSION(klev)   :: lev_index
131               
132#ifdef CPP_XIOS
133    ! ug Variables utilis\'ees pour r\'ecup\'erer le calendrier pour xios
134    INTEGER :: x_an, x_mois, x_jour
135    REAL :: x_heure
136    INTEGER :: ini_an, ini_mois, ini_jour
137    REAL :: ini_heure
138#endif
139
140    WRITE(lunout,*) 'Debut phys_output_mod.F90'
141    ! Initialisations (Valeurs par defaut
142
143    DO ilev=1,klev
144      Ahyb_bounds(ilev,1) = ap(ilev)
145      Ahyb_bounds(ilev,2) = ap(ilev+1)
146      Bhyb_bounds(ilev,1) = bp(ilev)
147      Bhyb_bounds(ilev,2) = bp(ilev+1)
148      lev_index(ilev) = REAL(ilev)
149    END DO
150
151    IF (.NOT. ALLOCATED(o_trac)) ALLOCATE(o_trac(nqtot))
152    IF (.NOT. ALLOCATED(o_trac_cum)) ALLOCATE(o_trac_cum(nqtot))
153    ALLOCATE(o_dtr_the(nqtot),o_dtr_con(nqtot),o_dtr_lessi_impa(nqtot))
154    ALLOCATE(o_dtr_lessi_nucl(nqtot),o_dtr_insc(nqtot),o_dtr_bcscav(nqtot))
155    ALLOCATE(o_dtr_evapls(nqtot),o_dtr_ls(nqtot),o_dtr_trsp(nqtot))
156    ALLOCATE(o_dtr_sscav(nqtot),o_dtr_sat(nqtot),o_dtr_uscav(nqtot))
157    ALLOCATE(o_dtr_dry(nqtot),o_dtr_vdf(nqtot))
158
159    levmax = (/ klev, klev, klev, klev, klev, klev, nlevSTD, nlevSTD, nlevSTD, klev /)
160
161    phys_out_filenames(1) = 'histmth'
162    phys_out_filenames(2) = 'histday'
163    phys_out_filenames(3) = 'histhf6h'
164    phys_out_filenames(4) = 'histhf3h'
165    phys_out_filenames(5) = 'histhf3hm'
166    phys_out_filenames(6) = 'histstn'
167    phys_out_filenames(7) = 'histmthNMC'
168    phys_out_filenames(8) = 'histdayNMC'
169    phys_out_filenames(9) = 'histhfNMC'
170    phys_out_filenames(10)= 'histstrataer'
171
172    type_ecri(1) = 'ave(X)'
173    type_ecri(2) = 'ave(X)'
174    type_ecri(3) = 'inst(X)'
175    type_ecri(4) = 'inst(X)'
176    type_ecri(5) = 'ave(X)'
177    type_ecri(6) = 'inst(X)'
178    type_ecri(7) = 'inst(X)'
179    type_ecri(8) = 'inst(X)'
180    type_ecri(9) = 'inst(X)'
181    type_ecri(10)= 'ave(X)'
182
183    clef_files(1) = ok_mensuel
184    clef_files(2) = ok_journe
185    clef_files(3) = ok_hf
186    clef_files(4) = ok_instan
187    clef_files(5) = ok_LES
188    clef_files(6) = ok_instan
189    clef_files(7) = ok_histNMC(1)
190    clef_files(8) = ok_histNMC(2)
191    clef_files(9) = ok_histNMC(3)
192#ifdef CPP_StratAer
193    clef_files(10)= .TRUE.
194#else
195    clef_files(10)= .FALSE.
196#endif
197
198    !sortir des fichiers "stations" si clef_stations(:)=.TRUE.
199    clef_stations(1) = .FALSE.
200    clef_stations(2) = .FALSE.
201    clef_stations(3) = .FALSE.
202    clef_stations(4) = .FALSE.
203    clef_stations(5) = .FALSE.
204    clef_stations(6) = .FALSE.
205    clef_stations(7) = .FALSE.
206    clef_stations(8) = .FALSE.
207    clef_stations(9) = .FALSE.
208    clef_stations(10)= .FALSE.
209
210    lev_files(1) = lev_histmth
211    lev_files(2) = lev_histday
212    lev_files(3) = lev_histhf
213    lev_files(4) = lev_histins
214    lev_files(5) = lev_histLES
215    lev_files(6) = lev_histins
216    lev_files(7) = levout_histNMC(1)
217    lev_files(8) = levout_histNMC(2)
218    lev_files(9) = levout_histNMC(3)
219    lev_files(10)= 5
220
221    ecrit_files(1) = ecrit_mth
222    ecrit_files(2) = ecrit_day
223    ecrit_files(3) = ecrit_hf
224    ecrit_files(4) = ecrit_ins
225    ecrit_files(5) = ecrit_LES
226    ecrit_files(6) = ecrit_ins
227    ecrit_files(7) = freq_outNMC(1)
228    ecrit_files(8) = freq_outNMC(2)
229    ecrit_files(9) = freq_outNMC(3)
230    ecrit_files(10)= ecrit_mth
231
232    !! Lectures des parametres de sorties dans physiq.def
233
234    CALL getin('phys_out_regfkey',phys_out_regfkey)
235    CALL getin('phys_out_lonmin',phys_out_lonmin)
236    CALL getin('phys_out_lonmax',phys_out_lonmax)
237    CALL getin('phys_out_latmin',phys_out_latmin)
238    CALL getin('phys_out_latmax',phys_out_latmax)
239    phys_out_levmin(:)=levmin(:)
240    CALL getin('phys_out_levmin',levmin)
241    phys_out_levmax(:)=levmax(:)
242    CALL getin('phys_out_levmax',levmax)
243    CALL getin('phys_out_filenames',phys_out_filenames)
244    phys_out_filekeys(:)=clef_files(:)
245    CALL getin('phys_out_filekeys',clef_files)
246    phys_out_filestations(:)=clef_stations(:)
247    CALL getin('phys_out_filestations',clef_stations)
248    phys_out_filelevels(:)=lev_files(:)
249    CALL getin('phys_out_filelevels',lev_files)
250    CALL getin('phys_out_filetimesteps',chtimestep)
251    phys_out_filetypes(:)=type_ecri(:)
252    CALL getin('phys_out_filetypes',type_ecri)
253
254    type_ecri_files(:)=type_ecri(:)
255
256!    if (ok_all_xml) phys_out_filelevels = 999
257
258    WRITE(lunout,*)'phys_out_lonmin=',phys_out_lonmin
259    WRITE(lunout,*)'phys_out_lonmax=',phys_out_lonmax
260    WRITE(lunout,*)'phys_out_latmin=',phys_out_latmin
261    WRITE(lunout,*)'phys_out_latmax=',phys_out_latmax
262    WRITE(lunout,*)'phys_out_filenames=',phys_out_filenames
263    WRITE(lunout,*)'phys_out_filetypes=',type_ecri
264    WRITE(lunout,*)'phys_out_filekeys=',clef_files
265    WRITE(lunout,*)'phys_out_filestations=',clef_stations
266    WRITE(lunout,*)'phys_out_filelevels=',lev_files
267    WRITE(lunout,*)'phys_out_regfkey=',phys_out_regfkey
268
269! A noter pour
270! l heure initiale - dans les fichiers histoire hist* - on met comme 
271! heure de debut soit la vraie heure (pour le 1D) soit 0h (pour le 3D)
272! afin d avoir une seule sortie mensuelle par mois lorsque l on tourne
273! par annee (IM).
274!
275     idayref = day_ref
276     IF (klon_glo==1) THEN
277       ! current_time (used to compute hour) is updated at the begining of
278       ! the physics; to set the correct outputs "initial time" we thus
279       ! have to use (hour-dtphys).
280         CALL ymds2ju(annee_ref, 1, idayref, hour-pdtphys, zjulian)
281         print *,'phys_output_mod: annee,iday,hour,zjulian=',annee_ref,idayref, hour, zjulian
282     ELSE
283         CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
284         CALL ymds2ju(annee_ref, 1, day_ini, start_time*rday, zjulian_start)
285     ENDIF
286
287#ifdef CPP_XIOS
288    ! ug R\'eglage du calendrier xios
289    !Temps julian => an, mois, jour, heure
290    CALL ju2ymds(zjulian, x_an, x_mois, x_jour, x_heure)
291    CALL ju2ymds(zjulian_start, ini_an, ini_mois, ini_jour, ini_heure)
292    CALL wxios_set_cal(dtime, calend, x_an, x_mois, x_jour, x_heure, ini_an, &
293                       ini_mois, ini_jour, ini_heure )
294#endif
295
296!!!!!!!!!!!!!!!!!!!!!!! Boucle sur les fichiers !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
297    ! Appel de histbeg et histvert pour creer le fichier et les niveaux verticaux !!
298    ! Appel des histbeg pour definir les variables (nom, moy ou inst, freq de sortie ..
299!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
300
301    zdtime_moy = dtime         ! Frequence ou l on moyenne
302
303
304  ecrit_files(7) = ecrit_files(1)
305  ecrit_files(8) = ecrit_files(2)
306  ecrit_files(9) = ecrit_files(3)
307
308  DO iff=1,nfiles
309
310       ! Calculate ecrit_files for all files
311      IF ( chtimestep(iff).eq.'Default' ) THEN
312          ! Par defaut ecrit_files = (ecrit_mensuel ecrit_jour ecrit_hf
313          ! ...)*86400.
314          ecrit_files(iff)=ecrit_files(iff)*86400.
315      ELSE IF (chtimestep(iff).eq.'-1') THEN
316          PRINT*,'ecrit_files(',iff,') < 0 so IOIPSL work on different'
317          PRINT*,'months length'
318          ecrit_files(iff)=-1.
319      ELSE
320       CALL convers_timesteps(chtimestep(iff),dtime,ecrit_files(iff))
321      ENDIF
322
323       WRITE(lunout,*)'ecrit_files(',iff,')= ',ecrit_files(iff)
324       zoutm(iff) = ecrit_files(iff) ! Frequence ou l on ecrit en seconde
325
326
327#ifdef CPP_XIOS
328!!! Ouverture de chaque fichier XIOS !!!!!!!!!!!
329    IF (.not. ok_all_xml) THEN
330      IF (prt_level >= 10) THEN
331        print*,'phys_output_open: call wxios_add_file with phys_out_filenames(iff)=',trim(phys_out_filenames(iff))   
332      ENDIF
333      CALL wxios_add_file(phys_out_filenames(iff),chtimestep(iff),lev_files(iff)) 
334    ENDIF
335
336!!! Declaration des axes verticaux de chaque fichier:
337    IF (prt_level >= 10) THEN
338      print*,'phys_output_open: Declare vertical axes for each file'
339    ENDIF
340
341   IF (iff.LE.6.OR.iff.EQ.10) THEN
342    CALL wxios_add_vaxis("presnivs", &
343            levmax(iff) - levmin(iff) + 1, presnivs(levmin(iff):levmax(iff)))
344    CALL wxios_add_vaxis("Ahyb", &
345            levmax(iff) - levmin(iff) + 1, aps(levmin(iff):levmax(iff)), positif='down', &
346            bnds=Ahyb_bounds(levmin(iff):levmax(iff),:))
347    CALL wxios_add_vaxis("Bhyb", &
348            levmax(iff) - levmin(iff) + 1, bps(levmin(iff):levmax(iff)), positif='down', &
349            bnds=Bhyb_bounds(levmin(iff):levmax(iff),:))
350    CALL wxios_add_vaxis("klev", levmax(iff) - levmin(iff) + 1, &
351                          lev_index(levmin(iff):levmax(iff)))
352    CALL wxios_add_vaxis("bnds", 2, (/1.,2./))
353
354     CALL wxios_add_vaxis("Alt", &
355            levmax(iff) - levmin(iff) + 1, pseudoalt)
356   ELSE
357    ! NMC files
358    CALL wxios_add_vaxis("plev", &
359            levmax(iff) - levmin(iff) + 1, rlevSTD(levmin(iff):levmax(iff)))
360   ENDIF
361#endif
362
363        IF (clef_files(iff)) THEN
364!!!!!!!!!!!!!!!!! Traitement dans le cas ou l'on veut stocker sur un domaine limite !!
365!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
366          IF (phys_out_regfkey(iff)) THEN
367             imin_ins=1
368             imax_ins=nbp_lon
369             jmin_ins=1
370             jmax_ins=jjmp1
371
372             ! correction abderr       
373             DO i=1,nbp_lon
374                WRITE(lunout,*)'io_lon(i)=',io_lon(i)
375                IF (io_lon(i).le.phys_out_lonmin(iff)) imin_ins=i
376                IF (io_lon(i).le.phys_out_lonmax(iff)) imax_ins=i+1
377             ENDDO
378
379             DO j=1,jjmp1
380                WRITE(lunout,*)'io_lat(j)=',io_lat(j)
381                IF (io_lat(j).ge.phys_out_latmin(iff)) jmax_ins=j+1
382                IF (io_lat(j).ge.phys_out_latmax(iff)) jmin_ins=j
383             ENDDO
384
385             WRITE(lunout,*)'On stoke le fichier histoire numero ',iff,' sur ', &
386                  imin_ins,imax_ins,jmin_ins,jmax_ins
387             WRITE(lunout,*)'longitudes : ', &
388                  io_lon(imin_ins),io_lon(imax_ins), &
389                  'latitudes : ', &
390                  io_lat(jmax_ins),io_lat(jmin_ins)
391
392             CALL histbeg(phys_out_filenames(iff),nbp_lon,io_lon,jjmp1,io_lat, &
393                  imin_ins,imax_ins-imin_ins+1, &
394                  jmin_ins,jmax_ins-jmin_ins+1, &
395                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
396!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
397             !IM fichiers stations
398          ELSE IF (clef_stations(iff)) THEN
399
400             IF (prt_level >= 10) THEN
401             WRITE(lunout,*)'phys_output_open: iff=',iff,'  phys_out_filenames(iff)=',phys_out_filenames(iff)
402             ENDIF
403             
404             CALL histbeg_phy_all(rlon,rlat,pim,tabij,ipt,jpt,plon,plat,plon_bounds,plat_bounds, &
405                  phys_out_filenames(iff), &
406                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
407          ELSE
408
409             IF (prt_level >= 10) THEN
410             WRITE(lunout,*)'phys_output_open: iff=',iff,'  phys_out_filenames(iff)=',phys_out_filenames(iff)
411             ENDIF
412
413             CALL histbeg_phy_all(phys_out_filenames(iff),itau_phy,zjulian,&
414                 dtime,nhorim(iff),nid_files(iff))
415          ENDIF
416
417#ifndef CPP_IOIPSL_NO_OUTPUT
418          IF (iff.LE.6.OR.iff.EQ.10) THEN
419             CALL histvert(nid_files(iff), "presnivs", "Vertical levels", "Pa", & 
420               levmax(iff) - levmin(iff) + 1, &
421               presnivs(levmin(iff):levmax(iff)), nvertm(iff),"down")
422!!!! Composantes de la coordonnee sigma-hybride
423          CALL histvert(nid_files(iff), "Ahyb","Ahyb comp of Hyb Cord ", "Pa", &
424               levmax(iff) - levmin(iff) + 1,aps,nvertap(iff))
425
426          CALL histvert(nid_files(iff), "Bhyb","Bhyb comp of Hyb Cord", " ", &
427               levmax(iff) - levmin(iff) + 1,bps,nvertbp(iff))
428
429          CALL histvert(nid_files(iff), "Alt","Height approx for scale heigh of 8km at levels", "Km", &                       
430               levmax(iff) - levmin(iff) + 1,pseudoalt,nvertAlt(iff))
431
432          ELSE
433          ! NMC files
434             CALL histvert(nid_files(iff), "plev", "pressure", "Pa", &
435               levmax(iff) - levmin(iff) + 1, &
436              rlevSTD(levmin(iff):levmax(iff)), nvertm(iff), "down")
437          ENDIF
438#endif
439
440     ENDIF ! clef_files
441
442       IF (nqtot>=nqo+1) THEN
443!
444            DO iq=nqo+1,nqtot
445            iiq=niadv(iq)
446            o_trac(iq-nqo) = ctrl_out((/ 1, 5, 5, 5, 10, 10, 11, 11, 11, 11 /), &
447                           tname(iiq),'Tracer '//ttext(iiq), "-",  &
448                           (/ '', '', '', '', '', '', '', '', '', '' /))
449            o_dtr_vdf(iq-nqo) = ctrl_out((/ 4, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
450                              'd'//trim(tname(iq))//'_vdf',  &
451                              'Tendance tracer '//ttext(iiq), "-" , &
452                              (/ '', '', '', '', '', '', '', '', '', '' /))
453
454            o_dtr_the(iq-nqo) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
455                              'd'//trim(tname(iq))//'_the', &
456                              'Tendance tracer '//ttext(iiq), "-", &
457                              (/ '', '', '', '', '', '', '', '', '', '' /))
458
459            o_dtr_con(iq-nqo) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
460                              'd'//trim(tname(iq))//'_con', &
461                              'Tendance tracer '//ttext(iiq), "-", &
462                              (/ '', '', '', '', '', '', '', '', '', '' /))
463
464            o_dtr_lessi_impa(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
465                                     'd'//trim(tname(iq))//'_lessi_impa', &
466                                     'Tendance tracer '//ttext(iiq), "-", &
467                                     (/ '', '', '', '', '', '', '', '', '', '' /))
468
469            o_dtr_lessi_nucl(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
470                                     'd'//trim(tname(iq))//'_lessi_nucl', &
471                                     'Tendance tracer '//ttext(iiq), "-", &
472                                     (/ '', '', '', '', '', '', '', '', '', '' /))
473
474            o_dtr_insc(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
475                               'd'//trim(tname(iq))//'_insc', &
476                               'Tendance tracer '//ttext(iiq), "-", &
477                               (/ '', '', '', '', '', '', '', '', '', '' /))
478
479            o_dtr_bcscav(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
480                                 'd'//trim(tname(iq))//'_bcscav', &
481                                 'Tendance tracer '//ttext(iiq), "-", &
482                                 (/ '', '', '', '', '', '', '', '', '', '' /))
483
484            o_dtr_evapls(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
485                                 'd'//trim(tname(iq))//'_evapls', &
486                                 'Tendance tracer '//ttext(iiq), "-", &
487                                 (/ '', '', '', '', '', '', '', '', '', '' /))
488
489            o_dtr_ls(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
490                             'd'//trim(tname(iq))//'_ls', &
491                             'Tendance tracer '//ttext(iiq), "-", &
492                             (/ '', '', '', '', '', '', '', '', '', '' /))
493
494            o_dtr_trsp(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
495                               'd'//trim(tname(iq))//'_trsp', &
496                               'Tendance tracer '//ttext(iiq), "-", &
497                               (/ '', '', '', '', '', '', '', '', '', '' /))
498
499            o_dtr_sscav(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
500                                'd'//trim(tname(iq))//'_sscav', &
501                                'Tendance tracer '//ttext(iiq), "-", &
502                                (/ '', '', '', '', '', '', '', '', '', '' /))
503
504            o_dtr_sat(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
505                               'd'//trim(tname(iq))//'_sat', &
506                               'Tendance tracer '//ttext(iiq), "-", &
507                               (/ '', '', '', '', '', '', '', '', '', '' /))
508
509            o_dtr_uscav(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
510                                'd'//trim(tname(iq))//'_uscav', &
511                                'Tendance tracer '//ttext(iiq), "-", &
512                                 (/ '', '', '', '', '', '', '', '', '', '' /))
513
514            o_dtr_dry(iq-nqo) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
515                              'cum'//'d'//trim(tname(iq))//'_dry', &
516                              'tracer tendency dry deposition'//ttext(iiq), "-", &
517                              (/ '', '', '', '', '', '', '', '', '', '' /))
518
519            o_trac_cum(iq-nqo) = ctrl_out((/ 1, 4, 10, 10, 10, 10, 11, 11, 11, 11 /), &
520                               'cum'//tname(iiq),&
521                               'Cumulated tracer '//ttext(iiq), "-", &
522                               (/ '', '', '', '', '', '', '', '', '', '' /))
523            ENDDO
524      ENDIF
525
526   ENDDO !  iff
527
528    ! Updated write frequencies due to phys_out_filetimesteps.
529    ! Write frequencies are now in seconds. 
530    ecrit_mth = ecrit_files(1)
531    ecrit_day = ecrit_files(2)
532    ecrit_hf  = ecrit_files(3)
533    ecrit_ins = ecrit_files(4)
534    ecrit_LES = ecrit_files(5)
535    ecrit_ins = ecrit_files(6)
536
537    IF (prt_level >= 10) THEN
538      WRITE(lunout,*)'swaero_diag=',swaero_diag
539      WRITE(lunout,*)'dryaod_diag=',dryaod_diag
540      WRITE(lunout,*)'ok_4xCO2atm=',ok_4xCO2atm
541      WRITE(lunout,*)'phys_output_open: ends here'
542    ENDIF
543
544  END SUBROUTINE phys_output_open
545
546  SUBROUTINE convers_timesteps(str,dtime,timestep)
547
548    use ioipsl
549    USE phys_cal_mod
550    USE time_phylmdz_mod, ONLY: day_ref, annee_ref
551    USE print_control_mod, ONLY: lunout
552
553    IMPLICIT NONE
554
555    CHARACTER(LEN=20)   :: str
556    CHARACTER(LEN=10)   :: type
557    INTEGER             :: ipos,il
558    real                :: ttt,xxx,timestep,dayseconde,dtime
559    parameter (dayseconde=86400.)
560
561    ipos=scan(str,'0123456789.',.TRUE.)
562    ! 
563    il=len_trim(str)
564    WRITE(lunout,*) "ipos = ", ipos
565    WRITE(lunout,*) "il = ", il
566    IF (ipos == 0) CALL abort_physic("convers_timesteps", "bad str", 1)
567    read(str(1:ipos),*) ttt
568    WRITE(lunout,*)ttt
569    type=str(ipos+1:il)
570
571    IF ( il == ipos ) THEN
572       type='day'
573    ENDIF
574
575    IF ( type == 'day'.or.type == 'days'.or.type == 'jours'.or.type == 'jour' ) timestep = ttt * dayseconde
576    IF ( type == 'mounths'.or.type == 'mth'.or.type == 'mois' ) THEN
577       WRITE(lunout,*)'annee_ref,day_ref mon_len',annee_ref,day_ref,mth_len
578       timestep = ttt * dayseconde * mth_len
579    ENDIF
580    IF ( type == 'hours'.or.type == 'hr'.or.type == 'heurs') timestep = ttt * dayseconde / 24.
581    IF ( type == 'mn'.or.type == 'minutes'  ) timestep = ttt * 60.
582    IF ( type == 's'.or.type == 'sec'.or.type == 'secondes'   ) timestep = ttt
583    IF ( type == 'TS' ) timestep = ttt * dtime
584
585    WRITE(lunout,*)'type =      ',type
586    WRITE(lunout,*)'nb j/h/m =  ',ttt
587    WRITE(lunout,*)'timestep(s)=',timestep
588
589  END SUBROUTINE convers_timesteps
590
591END MODULE phys_output_mod
Note: See TracBrowser for help on using the repository browser.