source: LMDZ6/trunk/libf/phylmdiso/phys_output_mod.F90 @ 4038

Last change on this file since 4038 was 3940, checked in by crisi, 3 years ago

replace files by symbloic liks from phylmdiso towards phylmd.
Many files at once

File size: 32.0 KB
Line 
1! $Id: phys_output_mod.F90 3792 2021-01-04 17:01:25Z evignon $
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       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, type_trac, &
38        nqtottr,itr_indice ! C Risi
39    USE ioipsl
40    USE phys_cal_mod, only : hour, calend
41    USE mod_phys_lmdz_para
42    !Martin
43    USE surface_data, ONLY : landice_opt
44    USE phys_output_ctrlout_mod
45    USE mod_grid_phy_lmdz, only: klon_glo,nbp_lon,nbp_lat
46    USE print_control_mod, ONLY: prt_level,lunout
47    USE vertical_layers_mod, ONLY: ap,bp,preff,presnivs, aps, bps, pseudoalt
48    USE time_phylmdz_mod, ONLY: day_ini, itau_phy, start_time, annee_ref, day_ref
49#ifdef CPP_XIOS
50    ! ug Pour les sorties XIOS
51    USE wxios
52#endif
53#ifdef ISO
54    USE infotrac_phy,ONLY: niso, ntraciso
55    USE isotopes_mod, ONLY: striso,iso_HTO
56#ifdef ISOTRAC
57    use isotrac_mod, only: index_zone,index_iso,strtrac
58#endif
59#endif
60
61    IMPLICIT NONE
62    include "clesphys.h"
63    include "thermcell.h"
64    include "YOMCST.h"
65
66    ! ug Nouveaux arguments n\'ecessaires au histwrite_mod:
67    INTEGER, INTENT(IN)                         :: ivap
68    INTEGER, DIMENSION(klon), INTENT(IN)        :: lmax_th
69    LOGICAL, INTENT(IN)                         :: ok_sync
70    LOGICAL, DIMENSION(klon, klev), INTENT(IN)  :: ptconv, ptconvth
71    REAL, INTENT(IN)                            :: pdtphys
72    REAL, DIMENSION(klon), INTENT(IN)           :: pphis
73    REAL, DIMENSION(klon, klev), INTENT(IN)     :: pplay, d_u, d_t
74    REAL, DIMENSION(klon, klev+1), INTENT(IN)   :: paprs
75    REAL, DIMENSION(klon,klev,nqtot), INTENT(IN):: qx, d_qx
76    REAL, DIMENSION(klon, klev), INTENT(IN)      :: zmasse
77
78
79    REAL,DIMENSION(klon),INTENT(IN) :: rlon
80    REAL,DIMENSION(klon),INTENT(IN) :: rlat
81    INTEGER, INTENT(IN)             :: pim
82    INTEGER, DIMENSION(pim)            :: tabij
83    INTEGER,DIMENSION(pim), INTENT(IN) :: ipt, jpt
84    REAL,DIMENSION(pim), INTENT(IN) :: plat, plon
85    REAL,DIMENSION(pim,2) :: plat_bounds, plon_bounds
86
87    INTEGER                               :: jjmp1
88    INTEGER                               :: nlevSTD, radpas
89    LOGICAL                               :: ok_mensuel, ok_journe, ok_hf, ok_instan
90    LOGICAL                               :: ok_LES,ok_ade,ok_aie
91    INTEGER                               :: flag_aerosol_strat
92    LOGICAL                               :: aerosol_couple
93    INTEGER, INTENT(IN)::  read_climoz ! read ozone climatology
94    !     Allowed values are 0, 1 and 2
95    !     0: do not read an ozone climatology
96    !     1: read a single ozone climatology that will be used day and night
97    !     2: read two ozone climatologies, the average day and night
98    !     climatology and the daylight climatology
99
100    REAL                                  :: dtime
101    INTEGER                               :: idayref
102    REAL                                  :: zjulian_start, zjulian
103    CHARACTER(LEN=4), DIMENSION(nlevSTD)  :: clevSTD
104    REAL, DIMENSION(nlevSTD)              :: rlevSTD
105    INTEGER                               :: nsrf, k, iq, iiq, iff, i, j, ilev
106    INTEGER                               :: itr ! C Risi
107    INTEGER                               :: naero
108    LOGICAL                               :: ok_veget
109    INTEGER                               :: iflag_pbl
110    INTEGER                               :: iflag_pbl_split
111    CHARACTER(LEN=4)                      :: bb2
112    CHARACTER(LEN=2)                      :: bb3
113    CHARACTER(LEN=6)                      :: type_ocean
114    INTEGER, DIMENSION(nbp_lon*jjmp1)         ::  ndex2d
115    INTEGER, DIMENSION(nbp_lon*jjmp1*klev)    :: ndex3d
116    INTEGER                               :: imin_ins, imax_ins
117    INTEGER                               :: jmin_ins, jmax_ins
118    INTEGER, DIMENSION(nfiles)            :: phys_out_levmin, phys_out_levmax
119    INTEGER, DIMENSION(nfiles)            :: phys_out_filelevels
120    CHARACTER(LEN=20), DIMENSION(nfiles)  :: chtimestep = (/ 'Default', 'Default', 'Default', 'Default', 'Default', &
121                                                             'Default', 'Default', 'Default', 'Default', 'Default' /)
122    LOGICAL, DIMENSION(nfiles)            :: phys_out_filekeys
123    LOGICAL, DIMENSION(nfiles)            :: phys_out_filestations
124
125#ifdef ISO
126      INTEGER  :: ixt,iiso,izone
127      CHARACTER*50 :: striso_sortie
128      integer :: lnblnk
129#endif
130
131!!!!!!!!!! stockage dans une region limitee pour chaque fichier !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
132    !                 entre [phys_out_lonmin,phys_out_lonmax] et [phys_out_latmin,phys_out_latmax]
133
134    LOGICAL, DIMENSION(nfiles), SAVE  :: phys_out_regfkey       = (/ .FALSE., .FALSE., .FALSE., .FALSE., .FALSE., &
135                                                                     .FALSE., .FALSE., .FALSE., .FALSE., .FALSE. /)
136    REAL, DIMENSION(nfiles), SAVE     :: phys_out_lonmin        = (/ -180., -180., -180., -180., -180., &
137                                                                     -180., -180., -180., -180., -180. /)
138    REAL, DIMENSION(nfiles), SAVE     :: phys_out_lonmax        = (/  180.,  180.,  180.,  180.,  180., &
139                                                                      180.,  180.,  180.,  180.,  180. /)
140    REAL, DIMENSION(nfiles), SAVE     :: phys_out_latmin        = (/  -90.,  -90.,  -90.,  -90.,  -90., &
141                                                                      -90.,  -90.,  -90.,  -90.,  -90. /)
142    REAL, DIMENSION(nfiles), SAVE     :: phys_out_latmax        = (/   90.,   90.,   90.,   90.,   90., &
143                                                                       90.,   90.,   90.,   90.,   90. /)
144    REAL, DIMENSION(klev,2) :: Ahyb_bounds, Bhyb_bounds
145    REAL, DIMENSION(klev+1)   :: lev_index
146               
147#ifdef CPP_XIOS
148    ! ug Variables utilis\'ees pour r\'ecup\'erer le calendrier pour xios
149    INTEGER :: x_an, x_mois, x_jour
150    REAL :: x_heure
151    INTEGER :: ini_an, ini_mois, ini_jour
152    REAL :: ini_heure
153#endif
154    INTEGER                         :: ISW
155    REAL, DIMENSION(NSW)            :: wl1_sun, wl2_sun !wavelength bounds (in um) for SW
156    REAL, DIMENSION(NSW)            :: wn1_sun, wn2_sun !wavenumber bounds (in m-1) for SW
157    REAL, DIMENSION(NSW)            :: spectband  !mean wavenumb. of each sp.band
158    REAL, DIMENSION(NSW,2)          :: spbnds_sun !bounds of spectband
159
160    WRITE(lunout,*) 'Debut phys_output_mod.F90'
161    ! Initialisations (Valeurs par defaut
162
163    DO ilev=1,klev
164      Ahyb_bounds(ilev,1) = ap(ilev)
165      Ahyb_bounds(ilev,2) = ap(ilev+1)
166      Bhyb_bounds(ilev,1) = bp(ilev)
167      Bhyb_bounds(ilev,2) = bp(ilev+1)
168      lev_index(ilev) = REAL(ilev)
169    END DO
170      lev_index(klev+1) = REAL(klev+1)
171
172    IF (.NOT. ALLOCATED(o_trac)) ALLOCATE(o_trac(nqtot))
173    IF (.NOT. ALLOCATED(o_trac_cum)) ALLOCATE(o_trac_cum(nqtot))
174    ALLOCATE(o_dtr_the(nqtot),o_dtr_con(nqtot),o_dtr_lessi_impa(nqtot))
175    ALLOCATE(o_dtr_lessi_nucl(nqtot),o_dtr_insc(nqtot),o_dtr_bcscav(nqtot))
176    ALLOCATE(o_dtr_evapls(nqtot),o_dtr_ls(nqtot),o_dtr_trsp(nqtot))
177    ALLOCATE(o_dtr_sscav(nqtot),o_dtr_sat(nqtot),o_dtr_uscav(nqtot))
178    ALLOCATE(o_dtr_dry(nqtot),o_dtr_vdf(nqtot))
179
180#ifdef ISO
181    ALLOCATE(o_xtprecip(ntraciso))
182    ALLOCATE(o_xtplul(ntraciso))
183    ALLOCATE(o_xtpluc(ntraciso))
184    ALLOCATE(o_xtevap(ntraciso))
185    ALLOCATE(o_xtovap(ntraciso))
186    ALLOCATE(o_xtoliq(ntraciso))
187    ALLOCATE(o_xtcond(ntraciso))
188    ALLOCATE(o_dxtdyn(ntraciso))
189    ALLOCATE(o_dxtldyn(ntraciso))
190    ALLOCATE(o_dxtcon(ntraciso))
191    ALLOCATE(o_dxtlsc(ntraciso))
192    ALLOCATE(o_dxteva(ntraciso))
193    ALLOCATE(o_dxtajs(ntraciso))
194    ALLOCATE(o_dxtvdf(ntraciso))
195    ALLOCATE(o_dxtthe(ntraciso))
196    ALLOCATE(o_dxtch4(ntraciso))
197    if (iso_HTO.gt.0) then
198      ALLOCATE(o_dxtprod_nucl(ntraciso))
199      ALLOCATE(o_dxtcosmo(ntraciso))
200      ALLOCATE(o_dxtdecroiss(ntraciso))
201    endif
202#endif
203
204    levmax = (/ klev, klev, klev, klev, klev, klev, nlevSTD, nlevSTD, nlevSTD, klev /)
205
206    phys_out_filenames(1) = 'histmth'
207    phys_out_filenames(2) = 'histday'
208    phys_out_filenames(3) = 'histhf6h'
209    phys_out_filenames(4) = 'histhf3h'
210    phys_out_filenames(5) = 'histhf3hm'
211    phys_out_filenames(6) = 'histstn'
212    phys_out_filenames(7) = 'histmthNMC'
213    phys_out_filenames(8) = 'histdayNMC'
214    phys_out_filenames(9) = 'histhfNMC'
215    phys_out_filenames(10)= 'histstrataer'
216
217    type_ecri(1) = 'ave(X)'
218    type_ecri(2) = 'ave(X)'
219    type_ecri(3) = 'inst(X)'
220    type_ecri(4) = 'inst(X)'
221    type_ecri(5) = 'ave(X)'
222    type_ecri(6) = 'inst(X)'
223    type_ecri(7) = 'inst(X)'
224    type_ecri(8) = 'inst(X)'
225    type_ecri(9) = 'inst(X)'
226    type_ecri(10)= 'ave(X)'
227
228    clef_files(1) = ok_mensuel
229    clef_files(2) = ok_journe
230    clef_files(3) = ok_hf
231    clef_files(4) = ok_instan
232    clef_files(5) = ok_LES
233    clef_files(6) = ok_instan
234    clef_files(7) = ok_histNMC(1)
235    clef_files(8) = ok_histNMC(2)
236    clef_files(9) = ok_histNMC(3)
237#ifdef CPP_StratAer
238    clef_files(10)= .TRUE.
239#else
240    clef_files(10)= .FALSE.
241#endif
242
243    !sortir des fichiers "stations" si clef_stations(:)=.TRUE.
244    clef_stations(1) = .FALSE.
245    clef_stations(2) = .FALSE.
246    clef_stations(3) = .FALSE.
247    clef_stations(4) = .FALSE.
248    clef_stations(5) = .FALSE.
249    clef_stations(6) = .FALSE.
250    clef_stations(7) = .FALSE.
251    clef_stations(8) = .FALSE.
252    clef_stations(9) = .FALSE.
253    clef_stations(10)= .FALSE.
254
255    lev_files(1) = lev_histmth
256    lev_files(2) = lev_histday
257    lev_files(3) = lev_histhf
258    lev_files(4) = lev_histins
259    lev_files(5) = lev_histLES
260    lev_files(6) = lev_histins
261    lev_files(7) = levout_histNMC(1)
262    lev_files(8) = levout_histNMC(2)
263    lev_files(9) = levout_histNMC(3)
264    lev_files(10)= 5
265
266    ecrit_files(1) = ecrit_mth
267    ecrit_files(2) = ecrit_day
268    ecrit_files(3) = ecrit_hf
269    ecrit_files(4) = ecrit_ins
270    ecrit_files(5) = ecrit_LES
271    ecrit_files(6) = ecrit_ins
272    ecrit_files(7) = freq_outNMC(1)
273    ecrit_files(8) = freq_outNMC(2)
274    ecrit_files(9) = freq_outNMC(3)
275    ecrit_files(10)= ecrit_mth
276
277    !! Lectures des parametres de sorties dans physiq.def
278
279    CALL getin('phys_out_regfkey',phys_out_regfkey)
280    CALL getin('phys_out_lonmin',phys_out_lonmin)
281    CALL getin('phys_out_lonmax',phys_out_lonmax)
282    CALL getin('phys_out_latmin',phys_out_latmin)
283    CALL getin('phys_out_latmax',phys_out_latmax)
284    phys_out_levmin(:)=levmin(:)
285    CALL getin('phys_out_levmin',levmin)
286    phys_out_levmax(:)=levmax(:)
287    CALL getin('phys_out_levmax',levmax)
288    CALL getin('phys_out_filenames',phys_out_filenames)
289    phys_out_filekeys(:)=clef_files(:)
290    CALL getin('phys_out_filekeys',clef_files)
291    phys_out_filestations(:)=clef_stations(:)
292    CALL getin('phys_out_filestations',clef_stations)
293    phys_out_filelevels(:)=lev_files(:)
294    CALL getin('phys_out_filelevels',lev_files)
295    CALL getin('phys_out_filetimesteps',chtimestep)
296    phys_out_filetypes(:)=type_ecri(:)
297    CALL getin('phys_out_filetypes',type_ecri)
298
299    type_ecri_files(:)=type_ecri(:)
300
301!    if (ok_all_xml) phys_out_filelevels = 999
302
303    WRITE(lunout,*)'phys_out_lonmin=',phys_out_lonmin
304    WRITE(lunout,*)'phys_out_lonmax=',phys_out_lonmax
305    WRITE(lunout,*)'phys_out_latmin=',phys_out_latmin
306    WRITE(lunout,*)'phys_out_latmax=',phys_out_latmax
307    WRITE(lunout,*)'phys_out_filenames=',phys_out_filenames
308    WRITE(lunout,*)'phys_out_filetypes=',type_ecri
309    WRITE(lunout,*)'phys_out_filekeys=',clef_files
310    WRITE(lunout,*)'phys_out_filestations=',clef_stations
311    WRITE(lunout,*)'phys_out_filelevels=',lev_files
312    WRITE(lunout,*)'phys_out_regfkey=',phys_out_regfkey
313
314! A noter pour
315! l heure initiale - dans les fichiers histoire hist* - on met comme 
316! heure de debut soit la vraie heure (pour le 1D) soit 0h (pour le 3D)
317! afin d avoir une seule sortie mensuelle par mois lorsque l on tourne
318! par annee (IM).
319!
320     idayref = day_ref
321     IF (klon_glo==1) THEN
322       ! current_time (used to compute hour) is updated at the begining of
323       ! the physics; to set the correct outputs "initial time" we thus
324       ! have to use (hour-dtphys).
325         CALL ymds2ju(annee_ref, 1, idayref, hour-pdtphys, zjulian)
326         print *,'phys_output_mod: annee,iday,hour,zjulian=',annee_ref,idayref, hour, zjulian
327     ELSE
328         CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
329         CALL ymds2ju(annee_ref, 1, day_ini, start_time*rday, zjulian_start)
330     ENDIF
331
332#ifdef CPP_XIOS
333    ! ug R\'eglage du calendrier xios
334    !Temps julian => an, mois, jour, heure
335    CALL ju2ymds(zjulian, x_an, x_mois, x_jour, x_heure)
336    CALL ju2ymds(zjulian_start, ini_an, ini_mois, ini_jour, ini_heure)
337    CALL wxios_set_cal(dtime, calend, x_an, x_mois, x_jour, x_heure, ini_an, &
338                       ini_mois, ini_jour, ini_heure )
339#endif
340
341!!!!!!!!!!!!!!!!!!!!!!! Boucle sur les fichiers !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
342    ! Appel de histbeg et histvert pour creer le fichier et les niveaux verticaux !!
343    ! Appel des histbeg pour definir les variables (nom, moy ou inst, freq de sortie ..
344!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345
346    zdtime_moy = dtime         ! Frequence ou l on moyenne
347
348
349  ecrit_files(7) = ecrit_files(1)
350  ecrit_files(8) = ecrit_files(2)
351  ecrit_files(9) = ecrit_files(3)
352
353  DO iff=1,nfiles
354
355       ! Calculate ecrit_files for all files
356      IF ( chtimestep(iff).eq.'Default' ) THEN
357          ! Par defaut ecrit_files = (ecrit_mensuel ecrit_jour ecrit_hf
358          ! ...)*86400.
359          ecrit_files(iff)=ecrit_files(iff)*86400.
360      ELSE IF (chtimestep(iff).eq.'-1') THEN
361          PRINT*,'ecrit_files(',iff,') < 0 so IOIPSL work on different'
362          PRINT*,'months length'
363          ecrit_files(iff)=-1.
364      ELSE
365       CALL convers_timesteps(chtimestep(iff),dtime,ecrit_files(iff))
366      ENDIF
367
368       WRITE(lunout,*)'ecrit_files(',iff,')= ',ecrit_files(iff)
369       zoutm(iff) = ecrit_files(iff) ! Frequence ou l on ecrit en seconde
370
371
372#ifdef CPP_XIOS
373!!! Ouverture de chaque fichier XIOS !!!!!!!!!!!
374    IF (.not. ok_all_xml) THEN
375      IF (prt_level >= 10) THEN
376        print*,'phys_output_open: call wxios_add_file with phys_out_filenames(iff)=',trim(phys_out_filenames(iff))   
377      ENDIF
378      CALL wxios_add_file(phys_out_filenames(iff),chtimestep(iff),lev_files(iff)) 
379    ENDIF
380
381!!! Declaration des axes verticaux de chaque fichier:
382    IF (prt_level >= 10) THEN
383      print*,'phys_output_open: Declare vertical axes for each file'
384    ENDIF
385
386   IF (iff.LE.6.OR.iff.EQ.10) THEN
387    CALL wxios_add_vaxis("presnivs", &
388            levmax(iff) - levmin(iff) + 1, presnivs(levmin(iff):levmax(iff)))
389    CALL wxios_add_vaxis("Ahyb", &
390            levmax(iff) - levmin(iff) + 1, aps(levmin(iff):levmax(iff)), positif='down', &
391            bnds=Ahyb_bounds(levmin(iff):levmax(iff),:))
392    CALL wxios_add_vaxis("Bhyb", &
393            levmax(iff) - levmin(iff) + 1, bps(levmin(iff):levmax(iff)), positif='down', &
394            bnds=Bhyb_bounds(levmin(iff):levmax(iff),:))
395    CALL wxios_add_vaxis("klev", levmax(iff) - levmin(iff) + 1, &
396                          lev_index(levmin(iff):levmax(iff)))
397    CALL wxios_add_vaxis("klevp1", klev+1, &
398                          lev_index(1:klev+1))
399    CALL wxios_add_vaxis("bnds", 2, (/1.,2./))
400
401     CALL wxios_add_vaxis("Alt", &
402            levmax(iff) - levmin(iff) + 1, pseudoalt)
403
404    IF (NSW.EQ.6) THEN
405!
406!wl1_sun: minimum bound of wavelength (in um)
407!
408      wl1_sun(1)=0.180
409      wl1_sun(2)=0.250
410      wl1_sun(3)=0.440
411      wl1_sun(4)=0.690
412      wl1_sun(5)=1.190
413      wl1_sun(6)=2.380
414!
415!wl2_sun: maximum bound of wavelength (in um)
416!
417      wl2_sun(1)=0.250
418      wl2_sun(2)=0.440
419      wl2_sun(3)=0.690
420      wl2_sun(4)=1.190
421      wl2_sun(5)=2.380
422      wl2_sun(6)=4.000
423!
424    ELSE IF(NSW.EQ.2) THEN
425!
426!wl1_sun: minimum bound of wavelength (in um)
427!
428      wl1_sun(1)=0.250
429      wl1_sun(2)=0.690
430!
431!wl2_sun: maximum bound of wavelength (in um)
432!
433      wl2_sun(1)=0.690
434      wl2_sun(2)=4.000
435    ENDIF
436
437    DO ISW=1, NSW
438     wn1_sun(ISW)=1.e+6/wl1_sun(ISW)
439     wn2_sun(ISW)=1.e+6/wl2_sun(ISW)
440     spbnds_sun(ISW,1)=wn2_sun(ISW)
441     spbnds_sun(ISW,2)=wn1_sun(ISW)
442     spectband(ISW)=(wn1_sun(ISW)+wn2_sun(ISW))/2
443    ENDDO
444!
445!!! ajout axe vertical spectband : solar band number
446    CALL wxios_add_vaxis("spectband", NSW, spectband, positif='down')
447   ELSE
448    ! NMC files
449    CALL wxios_add_vaxis("plev", &
450            levmax(iff) - levmin(iff) + 1, rlevSTD(levmin(iff):levmax(iff)))
451   ENDIF
452#endif
453
454        IF (clef_files(iff)) THEN
455!!!!!!!!!!!!!!!!! Traitement dans le cas ou l'on veut stocker sur un domaine limite !!
456!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
457          IF (phys_out_regfkey(iff)) THEN
458             imin_ins=1
459             imax_ins=nbp_lon
460             jmin_ins=1
461             jmax_ins=jjmp1
462
463             ! correction abderr       
464             DO i=1,nbp_lon
465                WRITE(lunout,*)'io_lon(i)=',io_lon(i)
466                IF (io_lon(i).le.phys_out_lonmin(iff)) imin_ins=i
467                IF (io_lon(i).le.phys_out_lonmax(iff)) imax_ins=i+1
468             ENDDO
469
470             DO j=1,jjmp1
471                WRITE(lunout,*)'io_lat(j)=',io_lat(j)
472                IF (io_lat(j).ge.phys_out_latmin(iff)) jmax_ins=j+1
473                IF (io_lat(j).ge.phys_out_latmax(iff)) jmin_ins=j
474             ENDDO
475
476             WRITE(lunout,*)'On stoke le fichier histoire numero ',iff,' sur ', &
477                  imin_ins,imax_ins,jmin_ins,jmax_ins
478             WRITE(lunout,*)'longitudes : ', &
479                  io_lon(imin_ins),io_lon(imax_ins), &
480                  'latitudes : ', &
481                  io_lat(jmax_ins),io_lat(jmin_ins)
482
483             CALL histbeg(phys_out_filenames(iff),nbp_lon,io_lon,jjmp1,io_lat, &
484                  imin_ins,imax_ins-imin_ins+1, &
485                  jmin_ins,jmax_ins-jmin_ins+1, &
486                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
487!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
488             !IM fichiers stations
489          ELSE IF (clef_stations(iff)) THEN
490
491             IF (prt_level >= 10) THEN
492             WRITE(lunout,*)'phys_output_open: iff=',iff,'  phys_out_filenames(iff)=',phys_out_filenames(iff)
493             ENDIF
494             
495             CALL histbeg_phy_all(rlon,rlat,pim,tabij,ipt,jpt,plon,plat,plon_bounds,plat_bounds, &
496                  phys_out_filenames(iff), &
497                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
498          ELSE
499
500             IF (prt_level >= 10) THEN
501             WRITE(lunout,*)'phys_output_open: iff=',iff,'  phys_out_filenames(iff)=',phys_out_filenames(iff)
502             ENDIF
503
504             CALL histbeg_phy_all(phys_out_filenames(iff),itau_phy,zjulian,&
505                 dtime,nhorim(iff),nid_files(iff))
506          ENDIF
507
508#ifndef CPP_IOIPSL_NO_OUTPUT
509          IF (iff.LE.6.OR.iff.EQ.10) THEN
510             CALL histvert(nid_files(iff), "presnivs", "Vertical levels", "Pa", & 
511               levmax(iff) - levmin(iff) + 1, &
512               presnivs(levmin(iff):levmax(iff)), nvertm(iff),"down")
513!!!! Composantes de la coordonnee sigma-hybride
514          CALL histvert(nid_files(iff), "Ahyb","Ahyb comp of Hyb Cord ", "Pa", &
515               levmax(iff) - levmin(iff) + 1,aps,nvertap(iff))
516
517          CALL histvert(nid_files(iff), "Bhyb","Bhyb comp of Hyb Cord", " ", &
518               levmax(iff) - levmin(iff) + 1,bps,nvertbp(iff))
519
520          CALL histvert(nid_files(iff), "Alt","Height approx for scale heigh of 8km at levels", "Km", &                       
521               levmax(iff) - levmin(iff) + 1,pseudoalt,nvertAlt(iff))
522
523          ELSE
524          ! NMC files
525             CALL histvert(nid_files(iff), "plev", "pressure", "Pa", &
526               levmax(iff) - levmin(iff) + 1, &
527              rlevSTD(levmin(iff):levmax(iff)), nvertm(iff), "down")
528          ENDIF
529#endif
530
531     ENDIF ! clef_files
532
533        write(lunout,*) 'phys_output_mid 496: nqtottr=',nqtottr
534        write(lunout,*) 'itr_indice=',itr_indice
535!       IF (nqtot>=nqo+1) THEN
536        IF (nqtottr>=1) THEN
537!
538            !DO iq=nqo+1,nqtot
539            ! C Risi: on modifie la boucle
540            do itr=1,nqtottr ! C Risi
541            iq=itr_indice(itr)  ! C Risi
542            write(*,*) 'phys_output_mid 503: itr=',itr
543 
544            iiq=niadv(iq)
545            o_trac(itr) = ctrl_out((/ 1, 5, 5, 5, 10, 10, 11, 11, 11, 11 /), &
546                           tname(iiq),'Tracer '//ttext(iiq), "-",  &
547                           (/ '', '', '', '', '', '', '', '', '', '' /))
548            o_dtr_vdf(itr) = ctrl_out((/ 4, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
549                              'd'//trim(tname(iq))//'_vdf',  &
550                              'Tendance tracer '//ttext(iiq), "-" , &
551                              (/ '', '', '', '', '', '', '', '', '', '' /))
552
553            o_dtr_the(itr) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
554                              'd'//trim(tname(iq))//'_the', &
555                              'Tendance tracer '//ttext(iiq), "-", &
556                              (/ '', '', '', '', '', '', '', '', '', '' /))
557
558            o_dtr_con(itr) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
559                              'd'//trim(tname(iq))//'_con', &
560                              'Tendance tracer '//ttext(iiq), "-", &
561                              (/ '', '', '', '', '', '', '', '', '', '' /))
562
563            o_dtr_lessi_impa(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
564                                     'd'//trim(tname(iq))//'_lessi_impa', &
565                                     'Tendance tracer '//ttext(iiq), "-", &
566                                     (/ '', '', '', '', '', '', '', '', '', '' /))
567
568            o_dtr_lessi_nucl(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
569                                     'd'//trim(tname(iq))//'_lessi_nucl', &
570                                     'Tendance tracer '//ttext(iiq), "-", &
571                                     (/ '', '', '', '', '', '', '', '', '', '' /))
572
573            o_dtr_insc(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
574                               'd'//trim(tname(iq))//'_insc', &
575                               'Tendance tracer '//ttext(iiq), "-", &
576                               (/ '', '', '', '', '', '', '', '', '', '' /))
577
578            o_dtr_bcscav(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
579                                 'd'//trim(tname(iq))//'_bcscav', &
580                                 'Tendance tracer '//ttext(iiq), "-", &
581                                 (/ '', '', '', '', '', '', '', '', '', '' /))
582
583            o_dtr_evapls(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
584                                 'd'//trim(tname(iq))//'_evapls', &
585                                 'Tendance tracer '//ttext(iiq), "-", &
586                                 (/ '', '', '', '', '', '', '', '', '', '' /))
587
588            o_dtr_ls(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
589                             'd'//trim(tname(iq))//'_ls', &
590                             'Tendance tracer '//ttext(iiq), "-", &
591                             (/ '', '', '', '', '', '', '', '', '', '' /))
592
593            o_dtr_trsp(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
594                               'd'//trim(tname(iq))//'_trsp', &
595                               'Tendance tracer '//ttext(iiq), "-", &
596                               (/ '', '', '', '', '', '', '', '', '', '' /))
597
598            o_dtr_sscav(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
599                                'd'//trim(tname(iq))//'_sscav', &
600                                'Tendance tracer '//ttext(iiq), "-", &
601                                (/ '', '', '', '', '', '', '', '', '', '' /))
602
603            o_dtr_sat(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
604                               'd'//trim(tname(iq))//'_sat', &
605                               'Tendance tracer '//ttext(iiq), "-", &
606                               (/ '', '', '', '', '', '', '', '', '', '' /))
607
608            o_dtr_uscav(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
609                                'd'//trim(tname(iq))//'_uscav', &
610                                'Tendance tracer '//ttext(iiq), "-", &
611                                 (/ '', '', '', '', '', '', '', '', '', '' /))
612
613            o_dtr_dry(itr) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11, 11 /), &
614                              'cum'//'d'//trim(tname(iq))//'_dry', &
615                              'tracer tendency dry deposition'//ttext(iiq), "-", &
616                              (/ '', '', '', '', '', '', '', '', '', '' /))
617
618            o_trac_cum(itr) = ctrl_out((/ 1, 4, 10, 10, 10, 10, 11, 11, 11, 11 /), &
619                               'cum'//tname(iiq),&
620                               'Cumulated tracer '//ttext(iiq), "-", &
621                               (/ '', '', '', '', '', '', '', '', '', '' /))
622            ENDDO
623      ENDIF
624
625   ENDDO !  iff
626
627        write(*,*) 'phys_output_mid 589'
628#ifdef ISO
629  do ixt=1,ntraciso
630     if (ixt.le.niso) then
631        striso_sortie=striso(ixt)
632     else
633#ifdef ISOTRAC
634        iiso=index_iso(ixt)
635        izone=index_zone(ixt)       
636        striso_sortie=striso(iiso)//strtrac(izone)
637#else
638        write(*,*) 'phys_output_mod 546: ixt,ntraciso=', ixt,ntraciso
639        stop
640#endif
641     endif
642
643   o_xtprecip(ixt)=ctrl_out((/ 1, 1, 1, 10, 5, 10, 11, 11, 11, 11/), &
644    'precip'//striso_sortie(1:lnblnk(striso_sortie)),  &
645    'Precip Totale liq+sol', 'kg/(s*m2)', (/ ('', i=1, 10) /))   
646   o_xtplul(ixt) = ctrl_out((/ 1, 1, 1, 10, 10, 10, 11, 11, 11, 11/), &
647    'plul'//striso_sortie(1:lnblnk(striso_sortie)),  &
648    'Large-scale Precip.', 'kg/(s*m2)', (/ ('', i=1, 10) /))
649   o_xtpluc(ixt) = ctrl_out((/ 1, 1, 1, 10, 5, 10, 11, 11, 11, 11/), &
650    'pluc'//striso_sortie(1:lnblnk(striso_sortie)),  &
651    'Convective Precip.', 'kg/(s*m2)', (/ ('', i=1, 10) /))
652   o_xtevap(ixt) = ctrl_out((/ 1, 1, 10, 10, 10, 10, 11, 11, 11, 11/), &
653    'evap'//striso_sortie(1:lnblnk(striso_sortie)),  &
654    'Evaporat', 'kg/(s*m2)', (/ ('', i=1, 10) /))
655   o_xtovap(ixt) = ctrl_out((/ 2, 3, 4, 10, 10, 10, 11, 11, 11, 11/), &
656    'ovap'//striso_sortie(1:lnblnk(striso_sortie)),  &
657    'Specific humidity', 'kg/kg', (/ ('', i=1, 10) /))
658   o_xtoliq(ixt) = ctrl_out((/ 2, 3, 4, 10, 10, 10, 11, 11, 11, 11/), &
659    'oliq'//striso_sortie(1:lnblnk(striso_sortie)),  &
660    'Liquid water', 'kg/kg', (/ ('', i=1, 10) /))
661   o_xtcond(ixt) = ctrl_out((/ 2, 3, 4, 10, 10, 10, 11, 11, 11, 11/), &
662    'ocond'//striso_sortie(1:lnblnk(striso_sortie)),  &
663    'Condensed water', 'kg/kg', (/ ('', i=1, 10) /))     
664   o_dxtdyn(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
665    'dqdyn'//striso_sortie(1:lnblnk(striso_sortie)),  &
666    'Dynamics dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
667   o_dxtldyn(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
668    'dqldyn'//striso_sortie(1:lnblnk(striso_sortie)),  &
669    'Dynamics dQL', '(kg/kg)/s', (/ ('', i=1, 10) /))
670   o_dxtcon(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
671    'dqcon'//striso_sortie(1:lnblnk(striso_sortie)),  &
672    'Convection dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
673   o_dxteva(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
674    'dqeva'//striso_sortie(1:lnblnk(striso_sortie)),  &
675    'Reevaporation dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
676   o_dxtlsc(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
677    'dqlsc'//striso_sortie(1:lnblnk(striso_sortie)),  &
678    'Condensation dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
679   o_dxtajs(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
680    'dqajs'//striso_sortie(1:lnblnk(striso_sortie)),  &
681    'Dry adjust. dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
682   o_dxtvdf(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
683    'dqvdf'//striso_sortie(1:lnblnk(striso_sortie)),  &
684    'Boundary-layer dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
685   o_dxtthe(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
686    'dqthe'//striso_sortie(1:lnblnk(striso_sortie)),  &
687    'Thermal dQ', '(kg/kg)/s', (/ ('', i=1, 10) /))
688
689   IF (ok_qch4) then
690     o_dxtch4(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
691        'dqch4'//striso_sortie(1:lnblnk(striso_sortie)),  &
692    'H2O due to CH4 oxidation & photolysis', '(kg/kg)/s', (/ ('', i=1, 10) /))
693   endif ! IF (ok_qch4) then
694
695   if (ixt.eq.iso_HTO) then
696      o_dxtprod_nucl(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
697        'dqprodnucl'//striso_sortie(1:lnblnk(striso_sortie)),  &
698        'dHTO/dt due to nuclear production', '(kg/kg)/s', (/ ('', i=1, 10) /))
699      o_dxtcosmo(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
700        'dqcosmo'//striso_sortie(1:lnblnk(striso_sortie)),  &
701        'dHTO/dt due to cosmogenic production', '(kg/kg)/s', (/ ('', i=1, 10) /))
702      o_dxtdecroiss(ixt) = ctrl_out((/ 4, 10, 10, 10, 10, 10, 11, 11, 11, 11/), &
703        'dqdecroiss'//striso_sortie(1:lnblnk(striso_sortie)),  &
704        'dHTO/dt due to radiative destruction', '(kg/kg)/s', (/ ('', i=1, 10) /))
705   endif !if (ixt.eq.iso_HTO) then
706  enddo !do ixt=1,niso
707#endif
708        write(*,*) 'phys_output_mid 596'
709
710
711    ! Updated write frequencies due to phys_out_filetimesteps.
712    ! Write frequencies are now in seconds. 
713    ecrit_mth = ecrit_files(1)
714    ecrit_day = ecrit_files(2)
715    ecrit_hf  = ecrit_files(3)
716    ecrit_ins = ecrit_files(4)
717    ecrit_LES = ecrit_files(5)
718    ecrit_ins = ecrit_files(6)
719
720    IF (prt_level >= 10) THEN
721      WRITE(lunout,*)'swaerofree_diag=',swaerofree_diag
722      WRITE(lunout,*)'swaero_diag=',swaero_diag
723      WRITE(lunout,*)'dryaod_diag=',dryaod_diag
724      WRITE(lunout,*)'ok_4xCO2atm=',ok_4xCO2atm
725      WRITE(lunout,*)'phys_output_open: ends here'
726    ENDIF
727
728  END SUBROUTINE phys_output_open
729
730  SUBROUTINE convers_timesteps(str,dtime,timestep)
731
732    use ioipsl
733    USE phys_cal_mod
734    USE time_phylmdz_mod, ONLY: day_ref, annee_ref
735    USE print_control_mod, ONLY: lunout
736
737    IMPLICIT NONE
738
739    CHARACTER(LEN=20)   :: str
740    CHARACTER(LEN=10)   :: type
741    INTEGER             :: ipos,il
742    real                :: ttt,xxx,timestep,dayseconde,dtime
743    parameter (dayseconde=86400.)
744
745    ipos=scan(str,'0123456789.',.TRUE.)
746    ! 
747    il=len_trim(str)
748    WRITE(lunout,*) "ipos = ", ipos
749    WRITE(lunout,*) "il = ", il
750    IF (ipos == 0) CALL abort_physic("convers_timesteps", "bad str", 1)
751    read(str(1:ipos),*) ttt
752    WRITE(lunout,*)ttt
753    type=str(ipos+1:il)
754
755    IF ( il == ipos ) THEN
756       type='day'
757    ENDIF
758
759    IF ( type == 'day'.or.type == 'days'.or.type == 'jours'.or.type == 'jour' ) timestep = ttt * dayseconde
760    IF ( type == 'mounths'.or.type == 'mth'.or.type == 'mois' ) THEN
761       WRITE(lunout,*)'annee_ref,day_ref mon_len',annee_ref,day_ref,mth_len
762       timestep = ttt * dayseconde * mth_len
763    ENDIF
764    IF ( type == 'hours'.or.type == 'hr'.or.type == 'heurs') timestep = ttt * dayseconde / 24.
765    IF ( type == 'mn'.or.type == 'minutes'  ) timestep = ttt * 60.
766    IF ( type == 's'.or.type == 'sec'.or.type == 'secondes'   ) timestep = ttt
767    IF ( type == 'TS' ) timestep = ttt * dtime
768
769    WRITE(lunout,*)'type =      ',type
770    WRITE(lunout,*)'nb j/h/m =  ',ttt
771    WRITE(lunout,*)'timestep(s)=',timestep
772
773  END SUBROUTINE convers_timesteps
774
775END MODULE phys_output_mod
Note: See TracBrowser for help on using the repository browser.