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

Last change on this file since 2114 was 2114, checked in by Laurent Fairhead, 10 years ago

Suite (et fin?) des modifications pour permettre un controle 'tout xml'
des fichiers de sorties XIOS.
Un nouveau paramètre logique est introduit, ok_all_xml, false par défaut, et lu
dans le run.def qui permet de faire du 'tout xml'


Follow-up modifications to ensure total xlm control over the output files
from XIOS.
A new logical parameter, ok_all_wml, is introduced. False by default, it is
read from the run.def file and, if true, will give over control to the
XIOS xml files

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