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

Last change on this file since 2284 was 2167, checked in by musat, 10 years ago

Correction sorties fichier 6heures en "inst" au lieu de "ave".
Changement noms par defaut des fichiers de sorties:
histmth, histday, histhf6h, histhf3h, histhf3hm
FH/IM

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