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

Last change on this file since 1910 was 1907, checked in by lguez, 11 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
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

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

  • 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: 22.6 KB
RevLine 
[1279]1! $Id: phys_output_mod.F90 1907 2013-11-26 13:10:46Z fairhead $
2!
[1813]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
[907]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
[1562]19CONTAINS
[907]20
21!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
22!!!!!!!!! Ouverture des fichier et definition des variable de sortie !!!!!!!!
[1562]23  !! histbeg, histvert et histdef
[907]24!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1562]25
[1539]26  SUBROUTINE phys_output_open(rlon,rlat,pim,tabij,ipt,jpt,plon,plat, &
[1828]27       jjmp1,nlevSTD,clevSTD,rlevSTD,nbteta, &
[1539]28       ctetaSTD, dtime, ok_veget, &
[1279]29       type_ocean, iflag_pbl,ok_mensuel,ok_journe, &
30       ok_hf,ok_instan,ok_LES,ok_ade,ok_aie, read_climoz, &
[1539]31       phys_out_filestations, &
[1807]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)   
[907]35
[1562]36    USE iophy
37    USE dimphy
38    USE infotrac
39    USE ioipsl
[1742]40    USE phys_cal_mod, only : hour
[1562]41    USE mod_phys_lmdz_para
42    USE aero_mod, only : naero_spc,name_aero
[1865]43        !Martin
44    USE surface_data, ONLY : ok_snow
[1807]45    USE phys_output_ctrlout_mod
[907]46
[1825]47#ifdef CPP_XIOS
48    ! ug Pour les sorties XIOS
49    USE wxios
50#endif
51
[1562]52    IMPLICIT NONE
53    include "dimensions.h"
54    include "temps.h"
55    include "clesphys.h"
56    include "thermcell.h"
57    include "comvert.h"
[1575]58    include "iniprint.h"
[907]59
[1807]60    ! ug Nouveaux arguments nécessaires au histwrite_mod:
61    INTEGER, INTENT(IN)                         :: ivap
62    INTEGER, DIMENSION(klon), INTENT(IN)        :: lmax_th
63    LOGICAL, INTENT(IN)                         :: ok_sync
64    LOGICAL, DIMENSION(klon, klev), INTENT(IN)  :: ptconv, ptconvth
65    REAL, INTENT(IN)                            :: pdtphys
66    REAL, DIMENSION(klon), INTENT(IN)           :: pphis
67    REAL, DIMENSION(klon, klev), INTENT(IN)     :: pplay, d_t
68    REAL, DIMENSION(klon, klev+1), INTENT(IN)   :: paprs
69    REAL, DIMENSION(klon,klev,nqtot), INTENT(IN):: qx, d_qx
70    REAL, DIMENSION(klon, llm), INTENT(IN)      :: zmasse
71
72
73    REAL,DIMENSION(klon),INTENT(IN) :: rlon
74    REAL,DIMENSION(klon),INTENT(IN) :: rlat
75    INTEGER, INTENT(IN)             :: pim
[1539]76    INTEGER, DIMENSION(pim)            :: tabij
[1807]77    INTEGER,DIMENSION(pim), INTENT(IN) :: ipt, jpt
78    REAL,DIMENSION(pim), INTENT(IN) :: plat, plon
79    REAL,DIMENSION(pim,2) :: plat_bounds, plon_bounds
[1539]80
[1791]81    INTEGER                               :: jjmp1
82    INTEGER                               :: nbteta, nlevSTD, radpas
83    LOGICAL                               :: ok_mensuel, ok_journe, ok_hf, ok_instan
84    LOGICAL                               :: ok_LES,ok_ade,ok_aie,flag_aerosol_strat
85    LOGICAL                               :: new_aod, aerosol_couple
[1807]86    INTEGER, INTENT(IN)::  read_climoz ! read ozone climatology
[1562]87    !     Allowed values are 0, 1 and 2
88    !     0: do not read an ozone climatology
89    !     1: read a single ozone climatology that will be used day and night
90    !     2: read two ozone climatologies, the average day and night
91    !     climatology and the daylight climatology
[1279]92
[1791]93    REAL                                  :: dtime
94    INTEGER                               :: idayref
95    REAL                                  :: zjulian
96    REAL, DIMENSION(klev)                 :: Ahyb, Bhyb, Alt
97    CHARACTER(LEN=4), DIMENSION(nlevSTD)  :: clevSTD
[1828]98    REAL, DIMENSION(nlevSTD)              :: rlevSTD
[1791]99    INTEGER                               :: nsrf, k, iq, iiq, iff, i, j, ilev
100    INTEGER                               :: naero
101    LOGICAL                               :: ok_veget
102    INTEGER                               :: iflag_pbl
103    CHARACTER(LEN=4)                      :: bb2
104    CHARACTER(LEN=2)                      :: bb3
105    CHARACTER(LEN=6)                      :: type_ocean
106    CHARACTER(LEN=3)                      :: ctetaSTD(nbteta)
107    REAL, DIMENSION(nfiles)               :: ecrit_files
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
[1884]114    CHARACTER(LEN=20), DIMENSION(nfiles)  :: chtimestep = (/ 'Default', 'Default', 'Default', 'Default', 'Default', &
115                                                             'Default', 'Default', 'Default', 'Default' /)
[1791]116    LOGICAL, DIMENSION(nfiles)            :: phys_out_filekeys
117    LOGICAL, DIMENSION(nfiles)            :: phys_out_filestations
[907]118
[1065]119!!!!!!!!!! stockage dans une region limitee pour chaque fichier !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1562]120    !                 entre [phys_out_lonmin,phys_out_lonmax] et [phys_out_latmin,phys_out_latmax]
[1054]121
[1828]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. /)                       
[1825]132#ifdef CPP_XIOS
133    ! ug Variables utilisées pour récupérer le calendrier pour xios
134    INTEGER :: x_an, x_mois, x_jour
135    REAL :: x_heure
136#endif
137
[1791]138    WRITE(lunout,*) 'Debut phys_output_mod.F90'
[1562]139    ! Initialisations (Valeurs par defaut
[1403]140
[1791]141    IF (.NOT. ALLOCATED(o_trac)) ALLOCATE(o_trac(nqtot))
142    IF (.NOT. ALLOCATED(o_trac_cum)) ALLOCATE(o_trac_cum(nqtot))
[1813]143    ALLOCATE(o_dtr_the(nqtot),o_dtr_con(nqtot),o_dtr_lessi_impa(nqtot))
144    ALLOCATE(o_dtr_lessi_nucl(nqtot),o_dtr_insc(nqtot),o_dtr_bcscav(nqtot))
145    ALLOCATE(o_dtr_evapls(nqtot),o_dtr_ls(nqtot),o_dtr_trsp(nqtot))
146    ALLOCATE(o_dtr_sscav(nqtot),o_dtr_sat(nqtot),o_dtr_uscav(nqtot))
147    ALLOCATE(o_dtr_dry(nqtot),o_dtr_vdf(nqtot))
[1403]148
[1813]149
[1828]150    levmax = (/ klev, klev, klev, klev, klev, klev, nlevSTD, nlevSTD, nlevSTD /)
[1065]151
[1562]152    phys_out_filenames(1) = 'histmth'
153    phys_out_filenames(2) = 'histday'
154    phys_out_filenames(3) = 'histhf'
155    phys_out_filenames(4) = 'histins'
156    phys_out_filenames(5) = 'histLES'
157    phys_out_filenames(6) = 'histstn'
[1828]158    phys_out_filenames(7) = 'histmthNMC'
159    phys_out_filenames(8) = 'histdayNMC'
160    phys_out_filenames(9) = 'histhfNMC.nc'
[907]161
[1562]162    type_ecri(1) = 'ave(X)'
163    type_ecri(2) = 'ave(X)'
164    type_ecri(3) = 'ave(X)'
165    type_ecri(4) = 'inst(X)'
166    type_ecri(5) = 'ave(X)'
167    type_ecri(6) = 'inst(X)'
[1828]168    type_ecri(7) = 'inst(X)'
169    type_ecri(8) = 'inst(X)'
170    type_ecri(9) = 'inst(X)'
[907]171
[1562]172    clef_files(1) = ok_mensuel
173    clef_files(2) = ok_journe
174    clef_files(3) = ok_hf
175    clef_files(4) = ok_instan
176    clef_files(5) = ok_LES
177    clef_files(6) = ok_instan
[1828]178    clef_files(7) = ok_histNMC(1)
179    clef_files(8) = ok_histNMC(2)
180    clef_files(9) = ok_histNMC(3)
[907]181
[1562]182    !sortir des fichiers "stations" si clef_stations(:)=.TRUE.
183    clef_stations(1) = .FALSE.
184    clef_stations(2) = .FALSE.
185    clef_stations(3) = .FALSE.
186    clef_stations(4) = .FALSE.
187    clef_stations(5) = .FALSE.
188    clef_stations(6) = .FALSE.
[1828]189    clef_stations(7) = .FALSE.
190    clef_stations(8) = .FALSE.
191    clef_stations(9) = .FALSE.
[1539]192
[1562]193    lev_files(1) = lev_histmth
194    lev_files(2) = lev_histday
195    lev_files(3) = lev_histhf
196    lev_files(4) = lev_histins
197    lev_files(5) = lev_histLES
198    lev_files(6) = lev_histins
[1828]199    lev_files(7) = levout_histNMC(1)
200    lev_files(8) = levout_histNMC(2)
201    lev_files(9) = levout_histNMC(3)
[907]202
[1562]203    ecrit_files(1) = ecrit_mth
204    ecrit_files(2) = ecrit_day
205    ecrit_files(3) = ecrit_hf
206    ecrit_files(4) = ecrit_ins
207    ecrit_files(5) = ecrit_LES
208    ecrit_files(6) = ecrit_ins
[1843]209    ecrit_files(7) = freq_outNMC(1)/86400.
210    ecrit_files(8) = freq_outNMC(2)/86400.
211    ecrit_files(9) = freq_outNMC(3)/86400.
[1279]212
[1562]213    !! Lectures des parametres de sorties dans physiq.def
[1279]214
[1791]215    CALL getin('phys_out_regfkey',phys_out_regfkey)
216    CALL getin('phys_out_lonmin',phys_out_lonmin)
217    CALL getin('phys_out_lonmax',phys_out_lonmax)
218    CALL getin('phys_out_latmin',phys_out_latmin)
219    CALL getin('phys_out_latmax',phys_out_latmax)
[1562]220    phys_out_levmin(:)=levmin(:)
[1791]221    CALL getin('phys_out_levmin',levmin)
[1562]222    phys_out_levmax(:)=levmax(:)
[1791]223    CALL getin('phys_out_levmax',levmax)
224    CALL getin('phys_out_filenames',phys_out_filenames)
[1562]225    phys_out_filekeys(:)=clef_files(:)
[1791]226    CALL getin('phys_out_filekeys',clef_files)
[1562]227    phys_out_filestations(:)=clef_stations(:)
[1791]228    CALL getin('phys_out_filestations',clef_stations)
[1562]229    phys_out_filelevels(:)=lev_files(:)
[1791]230    CALL getin('phys_out_filelevels',lev_files)
231    CALL getin('phys_out_filetimesteps',chtimestep)
[1562]232    phys_out_filetypes(:)=type_ecri(:)
[1791]233    CALL getin('phys_out_filetypes',type_ecri)
[1279]234
[1562]235    type_ecri_files(:)=type_ecri(:)
[1279]236
[1791]237    WRITE(lunout,*)'phys_out_lonmin=',phys_out_lonmin
238    WRITE(lunout,*)'phys_out_lonmax=',phys_out_lonmax
239    WRITE(lunout,*)'phys_out_latmin=',phys_out_latmin
240    WRITE(lunout,*)'phys_out_latmax=',phys_out_latmax
241    WRITE(lunout,*)'phys_out_filenames=',phys_out_filenames
242    WRITE(lunout,*)'phys_out_filetypes=',type_ecri
243    WRITE(lunout,*)'phys_out_filekeys=',clef_files
244    WRITE(lunout,*)'phys_out_filestations=',clef_stations
245    WRITE(lunout,*)'phys_out_filelevels=',lev_files
[1562]246
[1825]247#ifdef CPP_XIOS
248    ! ug Réglage du calendrier xios
249    !Temps julian => an, mois, jour, heure
250    CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian)
251    CALL ju2ymds(zjulian, x_an, x_mois, x_jour, x_heure)
252    CALL wxios_set_cal(dtime, calend, x_an, x_mois, x_jour, x_heure)
253#endif
254
[907]255!!!!!!!!!!!!!!!!!!!!!!! Boucle sur les fichiers !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1562]256    ! Appel de histbeg et histvert pour creer le fichier et les niveaux verticaux !!
257    ! Appel des histbeg pour definir les variables (nom, moy ou inst, freq de sortie ..
[907]258!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
259
[1807]260    zdtime_moy = dtime         ! Frequence ou l on moyenne
[1279]261
[1562]262    ! Calcul des Ahyb, Bhyb et Alt
[1791]263    DO k=1,klev
[1562]264       Ahyb(k)=(ap(k)+ap(k+1))/2.
265       Bhyb(k)=(bp(k)+bp(k+1))/2.
266       Alt(k)=log(preff/presnivs(k))*8.
[1791]267    ENDDO
[1562]268    !          if(prt_level.ge.1) then
[1791]269    WRITE(lunout,*)'Ap Hybrid = ',Ahyb(1:klev)
270    WRITE(lunout,*)'Bp Hybrid = ',Bhyb(1:klev)
271    WRITE(lunout,*)'Alt approx des couches pour une haut d echelle de 8km = ',Alt(1:klev)
[1562]272    !          endif
273    DO iff=1,nfiles
[929]274
[1641]275       ! Calculate ecrit_files for all files
[1884]276      IF ( chtimestep(iff).eq.'Default' ) then
277          ! Par defaut ecrit_files = (ecrit_mensuel ecrit_jour ecrit_hf
278          ! ...)*86400.
[1641]279          ecrit_files(iff)=ecrit_files(iff)*86400.
[1884]280      else
281       CALL convers_timesteps(chtimestep(iff),dtime,ecrit_files(iff))
282      ENDIF
283
[1791]284       WRITE(lunout,*)'ecrit_files(',iff,')= ',ecrit_files(iff)
[1641]285       zoutm(iff) = ecrit_files(iff) ! Frequence ou l on ecrit en seconde
[1279]286
[1852]287#ifdef CPP_XIOS
288    ! ug déclaration des axes verticaux de chaque fichier:
289    CALL wxios_add_vaxis("presnivs", phys_out_filenames(iff), &
290            levmax(iff) - levmin(iff) + 1, presnivs(levmin(iff):levmax(iff)))
291    CALL wxios_add_vaxis("Ahyb", phys_out_filenames(iff), &
292            levmax(iff) - levmin(iff) + 1, Ahyb)
293    CALL wxios_add_vaxis("Bhyb", phys_out_filenames(iff), &
294            levmax(iff) - levmin(iff) + 1, Bhyb)
295    CALL wxios_add_vaxis("Ahyb", phys_out_filenames(iff), &
296            levmax(iff) - levmin(iff) + 1, Alt)
297
298    IF (clef_files(iff)) THEN
299         CALL wxios_add_file(phys_out_filenames(iff),chtimestep(iff),lev_files(iff))
300    END IF
301#endif
302
[1641]303       IF (clef_files(iff)) THEN
[1279]304
[1562]305          idayref = day_ref
[1742]306!          CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)       
[1638]307! correction pour l heure initiale                               !jyg
308!                                                                !jyg
[1742]309          CALL ymds2ju(annee_ref, 1, idayref, hour, zjulian)         !jyg
310! correction pour l heure initiale                               !jyg
311!                                                                !jyg
312!!!      CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)       !jyg
313! correction pour l heure initiale                               !jyg
314!                                                                !jyg
[1638]315!      CALL ymds2ju(annee_ref, 1, idayref, hour, zjulian)         !jyg
[907]316
317!!!!!!!!!!!!!!!!! Traitement dans le cas ou l'on veut stocker sur un domaine limite !!
318!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1791]319          IF (phys_out_regfkey(iff)) then
[1562]320             imin_ins=1
321             imax_ins=iim
322             jmin_ins=1
323             jmax_ins=jjmp1
[907]324
[1562]325             ! correction abderr       
326             do i=1,iim
[1791]327                WRITE(lunout,*)'io_lon(i)=',io_lon(i)
328                IF (io_lon(i).le.phys_out_lonmin(iff)) imin_ins=i
329                IF (io_lon(i).le.phys_out_lonmax(iff)) imax_ins=i+1
[1562]330             enddo
[907]331
[1562]332             do j=1,jjmp1
[1791]333                WRITE(lunout,*)'io_lat(j)=',io_lat(j)
334                IF (io_lat(j).ge.phys_out_latmin(iff)) jmax_ins=j+1
335                IF (io_lat(j).ge.phys_out_latmax(iff)) jmin_ins=j
[1562]336             enddo
[907]337
[1791]338             WRITE(lunout,*)'On stoke le fichier histoire numero ',iff,' sur ', &
[1562]339                  imin_ins,imax_ins,jmin_ins,jmax_ins
[1791]340             WRITE(lunout,*)'longitudes : ', &
[1562]341                  io_lon(imin_ins),io_lon(imax_ins), &
342                  'latitudes : ', &
343                  io_lat(jmax_ins),io_lat(jmin_ins)
[907]344
[1562]345             CALL histbeg(phys_out_filenames(iff),iim,io_lon,jjmp1,io_lat, &
346                  imin_ins,imax_ins-imin_ins+1, &
347                  jmin_ins,jmax_ins-jmin_ins+1, &
348                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
[907]349!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[1562]350             !IM fichiers stations
[1791]351          else IF (clef_stations(iff)) THEN
[1539]352
[1791]353             WRITE(lunout,*)'phys_output_mod phys_out_filenames=',phys_out_filenames(iff)
[1539]354
[1791]355             CALL histbeg_phy_all(rlon,rlat,pim,tabij,ipt,jpt,plon,plat,plon_bounds,plat_bounds, &
[1562]356                  phys_out_filenames(iff), &
357                  itau_phy,zjulian,dtime,nhorim(iff),nid_files(iff))
358          else
[1825]359             CALL histbeg_phy_all(phys_out_filenames(iff),itau_phy,zjulian,&
[1852]360                 dtime,nhorim(iff),nid_files(iff))
[1562]361          endif
[907]362
[1852]363#ifndef CPP_NO_IOIPSL
[1828]364          if (iff.le.6) then
365             CALL histvert(nid_files(iff), "presnivs", "Vertical levels", "Pa", & 
[1562]366               levmax(iff) - levmin(iff) + 1, &
367               presnivs(levmin(iff):levmax(iff)), nvertm(iff),"down")
[1807]368!!!! Composantes de la coordonnee sigma-hybride
[1562]369          CALL histvert(nid_files(iff), "Ahyb","Ahyb comp of Hyb Cord ", "Pa", &
370               levmax(iff) - levmin(iff) + 1,Ahyb,nvertap(iff))
[907]371
[1562]372          CALL histvert(nid_files(iff), "Bhyb","Bhyb comp of Hyb Cord", " ", &
373               levmax(iff) - levmin(iff) + 1,Bhyb,nvertbp(iff))
[1279]374
[1828]375          CALL histvert(nid_files(iff), "Alt","Height approx for scale heigh of 8km at levels", "Km", &                       
[1562]376               levmax(iff) - levmin(iff) + 1,Alt,nvertAlt(iff))
[1279]377
[1828]378          else
379             CALL histvert(nid_files(iff), "plev", "pressure", "Pa", &
380               levmax(iff) - levmin(iff) + 1, &
381              rlevSTD(levmin(iff):levmax(iff)), nvertm(iff), "down")
382          endif
383#endif
[907]384
[1791]385      IF (nqtot>=3) THEN
386            DO iq=3,nqtot 
387            iiq=niadv(iq)
[1828]388            o_trac(iq-2) = ctrl_out((/ 4, 5, 1, 1, 1, 10, 11, 11, 11 /), &
389                           tname(iiq),'Tracer '//ttext(iiq), "-",  &
390                           (/ '', '', '', '', '', '', '', '', '' /))
[1813]391
[1828]392            o_dtr_vdf(iq-2) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11 /), &
393                              'd'//trim(tname(iq))//'_vdf',  &
394                              'Tendance tracer '//ttext(iiq), "-" , &
395                              (/ '', '', '', '', '', '', '', '', '' /))
[1813]396
[1828]397            o_dtr_the(iq-2) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11 /), &
398                              'd'//trim(tname(iq))//'_the', &
399                              'Tendance tracer '//ttext(iiq), "-", &
400                              (/ '', '', '', '', '', '', '', '', '' /))
[1813]401
[1828]402            o_dtr_con(iq-2) = ctrl_out((/ 5, 7, 7, 7, 10, 10, 11, 11, 11 /), &
403                              'd'//trim(tname(iq))//'_con', &
404                              'Tendance tracer '//ttext(iiq), "-", &
405                              (/ '', '', '', '', '', '', '', '', '' /))
[1813]406
[1828]407            o_dtr_lessi_impa(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
408                                     'd'//trim(tname(iq))//'_lessi_impa', &
409                                     'Tendance tracer '//ttext(iiq), "-", &
410                                     (/ '', '', '', '', '', '', '', '', '' /))
[1813]411
[1828]412            o_dtr_lessi_nucl(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
413                                     'd'//trim(tname(iq))//'_lessi_nucl', &
414                                     'Tendance tracer '//ttext(iiq), "-", &
415                                     (/ '', '', '', '', '', '', '', '', '' /))
[1813]416
[1828]417            o_dtr_insc(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
418                               'd'//trim(tname(iq))//'_insc', &
419                               'Tendance tracer '//ttext(iiq), "-", &
420                               (/ '', '', '', '', '', '', '', '', '' /))
[1813]421
[1828]422            o_dtr_bcscav(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
423                                 'd'//trim(tname(iq))//'_bcscav', &
424                                 'Tendance tracer '//ttext(iiq), "-", &
425                                 (/ '', '', '', '', '', '', '', '', '' /))
[1813]426
[1828]427            o_dtr_evapls(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
428                                 'd'//trim(tname(iq))//'_evapls', &
429                                 'Tendance tracer '//ttext(iiq), "-", &
430                                 (/ '', '', '', '', '', '', '', '', '' /))
[1813]431
[1828]432            o_dtr_ls(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
433                             'd'//trim(tname(iq))//'_ls', &
434                             'Tendance tracer '//ttext(iiq), "-", &
435                             (/ '', '', '', '', '', '', '', '', '' /))
[1813]436
[1828]437            o_dtr_trsp(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
438                               'd'//trim(tname(iq))//'_trsp', &
439                               'Tendance tracer '//ttext(iiq), "-", &
440                               (/ '', '', '', '', '', '', '', '', '' /))
[1813]441
[1828]442            o_dtr_sscav(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
443                                'd'//trim(tname(iq))//'_sscav', &
444                                'Tendance tracer '//ttext(iiq), "-", &
445                                (/ '', '', '', '', '', '', '', '', '' /))
[1813]446
[1828]447            o_dtr_sat(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
448                               'd'//trim(tname(iq))//'_sat', &
449                               'Tendance tracer '//ttext(iiq), "-", &
450                               (/ '', '', '', '', '', '', '', '', '' /))
[1813]451
[1828]452            o_dtr_uscav(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
453                                'd'//trim(tname(iq))//'_uscav', &
454                                'Tendance tracer '//ttext(iiq), "-", &
455                                 (/ '', '', '', '', '', '', '', '', '' /))
[1813]456
[1828]457            o_dtr_dry(iq-2) = ctrl_out((/ 7, 7, 7, 7, 10, 10, 11, 11, 11 /), &
458                              'cum'//'d'//trim(tname(iq))//'_dry', &
459                              'tracer tendency dry deposition'//ttext(iiq), "-", &
460                              (/ '', '', '', '', '', '', '', '', '' /))
[1813]461
[1828]462            o_trac_cum(iq-2) = ctrl_out((/ 3, 4, 10, 10, 10, 10, 11, 11, 11 /), &
463                               'cum'//tname(iiq),&
464                               'Cumulated tracer '//ttext(iiq), "-", &
465                               (/ '', '', '', '', '', '', '', '', '' /))
[1791]466            ENDDO
467      ENDIF
[1539]468
[1807]469    ENDIF ! clef_files
[1539]470
[1807]471    ENDDO !  iff
[1539]472
[1641]473    ! Updated write frequencies due to phys_out_filetimesteps.
474    ! Write frequencies are now in seconds. 
475    ecrit_mth = ecrit_files(1)
476    ecrit_day = ecrit_files(2)
477    ecrit_hf  = ecrit_files(3)
478    ecrit_ins = ecrit_files(4)
479    ecrit_LES = ecrit_files(5)
480    ecrit_ins = ecrit_files(6)
481
[1791]482    WRITE(lunout,*)'swaero_diag=',swaero_diag
483    WRITE(lunout,*)'Fin phys_output_mod.F90'
484  end SUBROUTINE phys_output_open
[907]485
486
487
[1562]488  SUBROUTINE convers_timesteps(str,dtime,timestep)
[1279]489
[1562]490    use ioipsl
491    USE phys_cal_mod
[1279]492
[1562]493    IMPLICIT NONE
[1279]494
[1791]495    CHARACTER(LEN=20)   :: str
496    CHARACTER(LEN=10)   :: type
497    INTEGER             :: ipos,il
[1562]498    real                :: ttt,xxx,timestep,dayseconde,dtime
499    parameter (dayseconde=86400.)
500    include "temps.h"
501    include "comconst.h"
[1575]502    include "iniprint.h"
[1279]503
[1791]504    ipos=scan(str,'0123456789.',.TRUE.)
[1562]505    ! 
506    il=len_trim(str)
[1791]507    WRITE(lunout,*)ipos,il
[1562]508    read(str(1:ipos),*) ttt
[1791]509    WRITE(lunout,*)ttt
[1562]510    type=str(ipos+1:il)
[1279]511
512
[1791]513    IF ( il == ipos ) then
[1562]514       type='day'
515    endif
[1279]516
[1791]517    IF ( type == 'day'.or.type == 'days'.or.type == 'jours'.or.type == 'jour' ) timestep = ttt * dayseconde
518    IF ( type == 'mounths'.or.type == 'mth'.or.type == 'mois' ) then
519       WRITE(lunout,*)'annee_ref,day_ref mon_len',annee_ref,day_ref,mth_len
[1562]520       timestep = ttt * dayseconde * mth_len
521    endif
[1791]522    IF ( type == 'hours'.or.type == 'hr'.or.type == 'heurs') timestep = ttt * dayseconde / 24.
523    IF ( type == 'mn'.or.type == 'minutes'  ) timestep = ttt * 60.
524    IF ( type == 's'.or.type == 'sec'.or.type == 'secondes'   ) timestep = ttt
525    IF ( type == 'TS' ) timestep = ttt * dtime
[1279]526
[1791]527    WRITE(lunout,*)'type =      ',type
528    WRITE(lunout,*)'nb j/h/m =  ',ttt
529    WRITE(lunout,*)'timestep(s)=',timestep
[1279]530
[1562]531  END SUBROUTINE convers_timesteps
[1279]532
[907]533END MODULE phys_output_mod
534
Note: See TracBrowser for help on using the repository browser.