source: LMDZ6/trunk/libf/phylmd/phys_output_mod.F90 @ 3651

Last change on this file since 3651 was 3630, checked in by Laurent Fairhead, 4 years ago

Parameter new_aod is not needed anymore as it is assumed to be true
all the time. This means that we cannot replay AR4 simulations with new
LMDZ sources (we probably couldn't anyway)
LF, OB

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