source: LMDZ6/trunk/libf/phylmd/phytrac_mod.F90 @ 4412

Last change on this file since 4412 was 4412, checked in by evignon, 16 months ago

travail de l'atelier nuages du 30/01/23: on renomme cldliq (ou radliq) en radocond
car le nom est tres trompeur + on ajoute des commentaires dans lscp_mod

  • 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:keywords set to Id
File size: 45.4 KB
Line 
1!$Id: phytrac_mod.F90 4412 2023-02-01 16:35:29Z evignon $
2MODULE phytrac_mod
3!=================================================================================
4! Interface between the LMDZ physical package and tracer computation.
5! Chemistry modules (INCA, Reprobus or the more specific traclmdz routine)
6! are called from phytrac.
7!
8!======================================================================
9! Auteur(s) FH
10! Objet: Moniteur general des tendances traceurs
11!
12! iflag_vdf_trac : Options for activating transport by vertical diffusion :
13!     1. notmal
14!     0. emission is injected in the first layer only, without diffusion
15!    -1  no emission & no diffusion
16! Modification 2013/07/22 : transformed into a module to pass tendencies to
17!     physics outputs. Additional keys for controling activation of sub processes.
18! Modification R. Pilon 10 octobre 2012 large scale scavenging incloud_scav + bc_scav
19! Modification R. Pilon 01 janvier 2012 transport+scavenging KE scheme : cvltr
20!=================================================================================
21
22!
23! Tracer tendencies, for outputs
24!-------------------------------
25  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_cl  ! Td couche limite/traceur
26  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_dec                            !RomP
27  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_cv  ! Td convection/traceur
28! RomP >>>
29  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_insc
30  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_bcscav
31  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_evapls
32  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_ls
33  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_trsp
34  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_sscav
35  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_sat
36  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_uscav
37  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: qPr,qDi ! concentration tra dans pluie,air descente insaturee
38  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: qPa,qMel
39  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: qTrdi,dtrcvMA ! conc traceur descente air insaturee et td convective MA
40! RomP <<<
41  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_th  ! Td thermique
42  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_lessi_impa ! Td du lessivage par impaction
43  REAL,DIMENSION(:,:,:),ALLOCATABLE,SAVE :: d_tr_lessi_nucl ! Td du lessivage par nucleation
44  REAL,DIMENSION(:,:),ALLOCATABLE,SAVE      :: qPrls      !jyg: concentration tra dans pluie LS a la surf.
45  REAL,DIMENSION(:,:),ALLOCATABLE,SAVE      :: d_tr_dry ! Td depot sec/traceur (1st layer),ALLOCATABLE,SAVE  jyg
46  REAL,DIMENSION(:,:),ALLOCATABLE,SAVE      :: flux_tr_dry ! depot sec/traceur (surface),ALLOCATABLE,SAVE    jyg
47
48!$OMP THREADPRIVATE(qPa,qMel,qTrdi,dtrcvMA,d_tr_th,d_tr_lessi_impa,d_tr_lessi_nucl)
49!$OMP THREADPRIVATE(d_tr_trsp,d_tr_sscav,d_tr_sat,d_tr_uscav,qPr,qDi)
50!$OMP THREADPRIVATE(d_tr_insc,d_tr_bcscav,d_tr_evapls,d_tr_ls,qPrls)
51!$OMP THREADPRIVATE(d_tr_cl,d_tr_dry,flux_tr_dry,d_tr_dec,d_tr_cv)
52
53
54CONTAINS
55
56  SUBROUTINE phytrac_init()
57    USE dimphy
58    USE infotrac_phy, ONLY: nbtr, type_trac
59    USE tracco2i_mod, ONLY: tracco2i_init
60    IMPLICIT NONE
61
62       ALLOCATE(d_tr_cl(klon,klev,nbtr),d_tr_dry(klon,nbtr))
63       ALLOCATE(flux_tr_dry(klon,nbtr),d_tr_dec(klon,klev,nbtr),d_tr_cv(klon,klev,nbtr))
64       ALLOCATE(d_tr_insc(klon,klev,nbtr),d_tr_bcscav(klon,klev,nbtr))
65       ALLOCATE(d_tr_evapls(klon,klev,nbtr),d_tr_ls(klon,klev,nbtr))
66       ALLOCATE(qPrls(klon,nbtr),d_tr_trsp(klon,klev,nbtr))
67       ALLOCATE(d_tr_sscav(klon,klev,nbtr),d_tr_sat(klon,klev,nbtr))
68       ALLOCATE(d_tr_uscav(klon,klev,nbtr),qPr(klon,klev,nbtr),qDi(klon,klev,nbtr))
69       ALLOCATE(qPa(klon,klev,nbtr),qMel(klon,klev,nbtr))
70       ALLOCATE(qTrdi(klon,klev,nbtr),dtrcvMA(klon,klev,nbtr))
71       ALLOCATE(d_tr_th(klon,klev,nbtr))
72       ALLOCATE(d_tr_lessi_impa(klon,klev,nbtr),d_tr_lessi_nucl(klon,klev,nbtr))
73
74
75
76    !===============================================================================
77    !    -- Do specific treatment according to chemestry model or local LMDZ tracers
78    !     
79    !===============================================================================
80    !   -- CO2 interactif --
81    IF(ANY(type_trac == ['co2i','inco'])) CALL tracco2i_init()
82
83       !   -- type_trac == 'co2i' ! PC
84       !   -- CO2 interactif --
85       !   -- source is updated with FF and BB emissions
86       !   -- and net fluxes from ocean and orchidee
87       !   -- sign convention : positive into the atmosphere
88
89  END SUBROUTINE phytrac_init
90
91  SUBROUTINE phytrac(                                 &
92       nstep,     julien,   gmtime,   debutphy,       &
93       lafin,     pdtphys,  u, v,     t_seri,         &
94       paprs,     pplay,    pmfu,     pmfd,           &
95       pen_u,     pde_u,    pen_d,    pde_d,          &
96       cdragh,    coefh,    fm_therm, entr_therm,     &
97       yu1,       yv1,      ftsol,    pctsrf,         &
98       ustar,     u10m,      v10m,                    &
99       wstar,     ale_bl,      ale_wake,              &
100       xlat,      xlon,                               &
101       frac_impa,frac_nucl,beta_fisrt,beta_v1,        &
102       presnivs,  pphis,    pphi,     albsol,         &
103       sh,        ch, rh,   cldfra,   rneb,           &
104       diafra,    cldliq,   itop_con, ibas_con,       &
105       pmflxr,    pmflxs,   prfl,     psfl,           &
106       da,        phi,      mp,       upwd,           &
107       phi2,      d1a,      dam,      sij, wght_cvfd, &   ! RomP +RL
108       wdtrainA,  wdtrainM, sigd,     clw, elij,      &   ! RomP
109       evap,      ep,       epmlmMm,  eplaMm,         &   ! RomP
110       dnwd,      aerosol_couple,     flxmass_w,      &
111       tau_aero,  piz_aero,  cg_aero, ccm,            &
112       rfname,                                        &
113       d_tr_dyn,                                      &   ! RomP
114       tr_seri, init_source)         
115    !
116    !======================================================================
117    ! Auteur(s) FH
118    ! Objet: Moniteur general des tendances traceurs
119    ! Modification R. Pilon 01 janvier 2012 transport+scavenging KE scheme : cvltr
120    ! Modification R. Pilon 10 octobre 2012 large scale scavenging incloud_scav + bc_scav
121    !======================================================================
122
123    USE ioipsl
124    USE phys_cal_mod, only : hour
125    USE dimphy
126    USE infotrac_phy, ONLY: nbtr, nqCO2, type_trac, conv_flg, pbl_flg
127    USE strings_mod,  ONLY: int2str
128    USE mod_grid_phy_lmdz
129    USE mod_phys_lmdz_para
130    USE iophy
131    USE traclmdz_mod
132    USE tracinca_mod
133    USE tracreprobus_mod
134    USE indice_sol_mod
135    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
136    USE print_control_mod, ONLY: lunout
137    USE aero_mod, ONLY : naero_grp
138
139    USE tracco2i_mod
140
141#ifdef CPP_StratAer
142    USE traccoag_mod
143    USE phys_local_var_mod, ONLY: mdw
144    USE phys_local_var_mod, ONLY: budg_dep_dry_ocs,   budg_dep_wet_ocs
145    USE phys_local_var_mod, ONLY: budg_dep_dry_so2,   budg_dep_wet_so2
146    USE phys_local_var_mod, ONLY: budg_dep_dry_h2so4, budg_dep_wet_h2so4
147    USE phys_local_var_mod, ONLY: budg_dep_dry_part,  budg_dep_wet_part
148    USE infotrac_phy, ONLY: nbtr_sulgas, id_OCS_strat, id_SO2_strat, id_H2SO4_strat
149    USE aerophys
150#endif
151
152    IMPLICIT NONE
153
154    INCLUDE "YOMCST.h"
155    INCLUDE "clesphys.h"
156    !==========================================================================
157    !                   -- ARGUMENT DESCRIPTION --
158    !==========================================================================
159
160    ! Input arguments
161    !----------------
162    !Configuration grille,temps:
163    INTEGER,INTENT(IN) :: nstep      ! Appel physique
164    INTEGER,INTENT(IN) :: julien     ! Jour julien
165    REAL,INTENT(IN)    :: gmtime     ! Heure courante
166    REAL,INTENT(IN)    :: pdtphys    ! Pas d'integration pour la physique (seconde)
167    LOGICAL,INTENT(IN) :: debutphy   ! le flag de l'initialisation de la physique
168    LOGICAL,INTENT(IN) :: lafin      ! le flag de la fin de la physique
169
170    REAL,DIMENSION(klon),INTENT(IN) :: xlat    ! latitudes pour chaque point
171    REAL,DIMENSION(klon),INTENT(IN) :: xlon    ! longitudes pour chaque point
172    !
173    !Physique:
174    !--------
175    REAL,DIMENSION(klon,klev),INTENT(IN)   :: t_seri  ! Temperature
176    REAL,DIMENSION(klon,klev),INTENT(IN)   :: u       ! variable not used
177    REAL,DIMENSION(klon,klev),INTENT(IN)   :: v       ! variable not used
178    REAL,DIMENSION(klon,klev),INTENT(IN)   :: sh      ! humidite specifique
179    REAL,DIMENSION(klon,klev),INTENT(IN)   :: rh      ! humidite relative
180    REAL,DIMENSION(klon,klev),INTENT(IN)   :: ch      ! eau liquide (+ glace si le traceur existe)
181    REAL,DIMENSION(klon,klev+1),INTENT(IN) :: paprs   ! pression pour chaque inter-couche (en Pa)
182    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pplay   ! pression pour le mileu de chaque couche (en Pa)
183    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pphi    ! geopotentiel
184    REAL,DIMENSION(klon),INTENT(IN)        :: pphis
185    REAL,DIMENSION(klev),INTENT(IN)        :: presnivs
186    REAL,DIMENSION(klon,klev),INTENT(IN)   :: cldliq  ! eau condensee totale
187    REAL,DIMENSION(klon,klev),INTENT(IN)   :: cldfra  ! fraction nuageuse (tous les nuages)
188    REAL,DIMENSION(klon,klev),INTENT(IN)   :: diafra  ! fraction nuageuse (convection ou stratus artificiels)
189    REAL,DIMENSION(klon,klev),INTENT(IN)   :: rneb    ! fraction nuageuse (grande echelle)
190    !
191    REAL                                   :: ql_incl ! contenu en eau liquide nuageuse dans le nuage ! ql_incl=oliq/rneb
192    REAL,DIMENSION(klon,klev),INTENT(IN)   :: beta_fisrt ! taux de conversion de l'eau cond (de fisrtilp)
193    REAL,DIMENSION(klon,klev),INTENT(out)  :: beta_v1    ! -- (originale version)
194
195    !
196    INTEGER,DIMENSION(klon),INTENT(IN)     :: itop_con
197    INTEGER,DIMENSION(klon),INTENT(IN)     :: ibas_con
198    REAL,DIMENSION(klon),INTENT(IN)        :: albsol  ! albedo surface
199    !
200    !Dynamique
201    !--------
202    REAL,DIMENSION(klon,klev,nbtr),INTENT(IN)    :: d_tr_dyn
203    !
204    !Convection:
205    !----------
206    REAL,DIMENSION(klon,klev),INTENT(IN) :: pmfu  ! flux de masse dans le panache montant
207    REAL,DIMENSION(klon,klev),INTENT(IN) :: pmfd  ! flux de masse dans le panache descendant
208    REAL,DIMENSION(klon,klev),INTENT(IN) :: pen_u ! flux entraine dans le panache montant
209    REAL,DIMENSION(klon,klev),INTENT(IN) :: pde_u ! flux detraine dans le panache montant
210    REAL,DIMENSION(klon,klev),INTENT(IN) :: pen_d ! flux entraine dans le panache descendant
211    REAL,DIMENSION(klon,klev),INTENT(IN) :: pde_d ! flux detraine dans le panache descendant
212
213    !...Tiedke     
214    REAL,DIMENSION(klon,klev+1),INTENT(IN)   :: pmflxr, pmflxs ! Flux precipitant de pluie, neige aux interfaces [convection]
215    REAL,DIMENSION(klon,klev+1),INTENT(IN)   :: prfl, psfl ! Flux precipitant de pluie, neige aux interfaces [large-scale]
216
217    LOGICAL,INTENT(IN)                       :: aerosol_couple
218    REAL,DIMENSION(klon,klev),INTENT(IN)     :: flxmass_w
219    REAL,DIMENSION(klon,klev,naero_grp,2),INTENT(IN) :: tau_aero
220    REAL,DIMENSION(klon,klev,naero_grp,2),INTENT(IN) :: piz_aero
221    REAL,DIMENSION(klon,klev,naero_grp,2),INTENT(IN) :: cg_aero
222    CHARACTER(len=4),DIMENSION(naero_grp),INTENT(IN) :: rfname
223    REAL,DIMENSION(klon,klev,2),INTENT(IN)   :: ccm
224    !... K.Emanuel
225    REAL,DIMENSION(klon,klev),INTENT(IN)     :: da
226    REAL,DIMENSION(klon,klev,klev),INTENT(IN):: phi
227    ! RomP >>>
228    REAL,DIMENSION(klon,klev),INTENT(IN)      :: d1a,dam
229    REAL,DIMENSION(klon,klev,klev),INTENT(IN) :: phi2
230    !
231    REAL,DIMENSION(klon,klev),INTENT(IN)      :: wdtrainA
232    REAL,DIMENSION(klon,klev),INTENT(IN)      :: wdtrainM
233    REAL,DIMENSION(klon),INTENT(IN)           :: sigd
234    ! ---- RomP flux entraine, detraine et precipitant kerry Emanuel
235    REAL,DIMENSION(klon,klev),INTENT(IN)      :: evap
236    REAL,DIMENSION(klon,klev),INTENT(IN)      :: ep
237    REAL,DIMENSION(klon,klev,klev),INTENT(IN) :: sij
238    REAL,DIMENSION(klon,klev),INTENT(IN)      :: wght_cvfd          !RL
239    REAL,DIMENSION(klon,klev,klev),INTENT(IN) :: elij
240    REAL,DIMENSION(klon,klev,klev),INTENT(IN) :: epmlmMm
241    REAL,DIMENSION(klon,klev),INTENT(IN)      :: eplaMm
242    REAL,DIMENSION(klon,klev),INTENT(IN)      :: clw
243    ! RomP <<<
244
245    !
246    REAL,DIMENSION(klon,klev),INTENT(IN)     :: mp
247    REAL,DIMENSION(klon,klev),INTENT(IN)     :: upwd      ! saturated updraft mass flux
248    REAL,DIMENSION(klon,klev),INTENT(IN)     :: dnwd      ! saturated downdraft mass flux
249    !
250    !Thermiques:
251    !----------
252    REAL,DIMENSION(klon,klev+1),INTENT(IN)   :: fm_therm
253    REAL,DIMENSION(klon,klev),INTENT(IN)     :: entr_therm
254    !
255    !Couche limite:
256    !--------------
257    !
258    REAL,DIMENSION(:),INTENT(IN)   :: cdragh          ! (klon) coeff drag pour T et Q
259    REAL,DIMENSION(:,:),INTENT(IN) :: coefh           ! (klon,klev) coeff melange CL (m**2/s)
260    REAL,DIMENSION(:),INTENT(IN)   :: ustar,u10m,v10m ! (klon) u* & vent a 10m (m/s)
261    REAL,DIMENSION(:),INTENT(IN)   :: wstar,ale_bl,ale_wake ! (klon) w* and Avail. Lifting Ener.
262    REAL,DIMENSION(:),INTENT(IN)   :: yu1             ! (klon) vents au premier niveau
263    REAL,DIMENSION(:),INTENT(IN)   :: yv1             ! (klon) vents au premier niveau
264    !
265    !Lessivage:
266    !----------
267    !
268    REAL, DIMENSION(:), ALLOCATABLE, SAVE :: ccntrAA
269    REAL, DIMENSION(:), ALLOCATABLE, SAVE :: ccntrENV
270    REAL, DIMENSION(:), ALLOCATABLE, SAVE :: coefcoli
271    LOGICAL, DIMENSION(:), ALLOCATABLE, SAVE :: flag_cvltr
272!$OMP THREADPRIVATE(ccntrAA,ccntrENV,coefcoli,flag_cvltr)
273    REAL, DIMENSION(klon,klev) :: ccntrAA_3d
274    REAL, DIMENSION(klon,klev) :: ccntrENV_3d
275    REAL, DIMENSION(klon,klev) :: coefcoli_3d
276    !
277    ! pour le ON-LINE
278    !
279    REAL,DIMENSION(klon,klev),INTENT(IN) :: frac_impa ! fraction d'aerosols non impactes
280    REAL,DIMENSION(klon,klev),INTENT(IN) :: frac_nucl ! fraction d'aerosols non nuclees
281
282    ! Arguments necessaires pour les sources et puits de traceur:
283    REAL,DIMENSION(klon,nbsrf),INTENT(IN) :: ftsol  ! Temperature du sol (surf)(Kelvin)
284    REAL,DIMENSION(klon,nbsrf),INTENT(IN) :: pctsrf ! Pourcentage de sol (nature du sol)
285
286#ifdef CPP_StratAer
287    REAL,DIMENSION(klon)           :: v_dep_dry !dry deposition velocity of stratospheric sulfate in m/s
288#endif
289    ! Output argument
290    !----------------
291    REAL,DIMENSION(klon,klev,nbtr),INTENT(INOUT) :: tr_seri ! Concentration Traceur [U/KgA]
292    REAL,DIMENSION(klon,klev)                    :: sourceBE
293    REAL,DIMENSION(klon,nbtr), INTENT(IN) :: init_source
294
295    !=======================================================================================
296    !                        -- LOCAL VARIABLES --
297    !=======================================================================================
298
299    INTEGER :: i, k, it
300    INTEGER :: nsplit
301
302    !Sources et Reservoirs de traceurs (ex:Radon):
303    !--------------------------------------------
304    !
305    REAL,DIMENSION(:,:),ALLOCATABLE,SAVE :: source  ! a voir lorsque le flux de surface est prescrit
306!$OMP THREADPRIVATE(source)
307
308    !
309    !Entrees/Sorties:
310    !---------------
311    INTEGER                   :: iiq, ierr
312    INTEGER                   :: nhori, nvert
313    REAL                      :: zsto, zout, zjulian
314    INTEGER,SAVE              :: nid_tra     ! pointe vers le fichier histrac.nc         
315!$OMP THREADPRIVATE(nid_tra)
316    REAL,DIMENSION(klon)      :: zx_tmp_fi2d ! variable temporaire grille physique
317    INTEGER                   :: itau_w      ! pas de temps ecriture = nstep + itau_phy
318    LOGICAL,PARAMETER         :: ok_sync=.TRUE.
319    !
320    ! Nature du traceur
321    !------------------
322    LOGICAL,DIMENSION(:),ALLOCATABLE,SAVE :: aerosol  ! aerosol(it) = true  => aerosol => lessivage
323!$OMP THREADPRIVATE(aerosol)                        ! aerosol(it) = false => gaz
324    REAL,DIMENSION(klon,klev)             :: delp     ! epaisseur de couche (Pa)
325    !
326    ! Tendances de traceurs (Td) et flux de traceurs:
327    !------------------------
328    REAL,DIMENSION(klon,klev)      :: d_tr     ! Td dans l'atmosphere
329    REAL,DIMENSION(klon,klev)      :: Mint
330    REAL,DIMENSION(klon,klev,nbtr) :: zmfd1a
331    REAL,DIMENSION(klon,klev,nbtr) :: zmfdam
332    REAL,DIMENSION(klon,klev,nbtr) :: zmfphi2
333    ! Physique
334    !----------
335    REAL,DIMENSION(klon,klev,nbtr) :: flestottr ! flux de lessivage dans chaque couche
336    REAL,DIMENSION(klon,klev)      :: zmasse    ! densite atmospherique Kg/m2
337    REAL,DIMENSION(klon,klev)      :: ztra_th
338    !PhH
339    REAL,DIMENSION(klon,klev)      :: zrho
340    REAL,DIMENSION(klon,klev)      :: zdz
341    REAL                           :: evaplsc,dx,beta ! variable pour lessivage Genthon
342    REAL,DIMENSION(klon)           :: his_dh          ! ---
343    ! in-cloud scav variables
344    REAL           :: ql_incloud_ref     ! ref value of in-cloud condensed water content
345
346    !Controles:
347    !---------
348    INTEGER,SAVE :: iflag_vdf_trac,iflag_con_trac,iflag_the_trac
349    INTEGER,SAVE  :: iflag_con_trac_omp, iflag_vdf_trac_omp,iflag_the_trac_omp
350!$OMP THREADPRIVATE(iflag_vdf_trac,iflag_con_trac,iflag_the_trac)
351
352    LOGICAL,SAVE :: lessivage
353!$OMP THREADPRIVATE(lessivage)
354
355    !RomP >>>
356    INTEGER,SAVE  :: iflag_lscav_omp,iflag_lscav
357    REAL, SAVE ::   ccntrAA_in,ccntrAA_omp
358    REAL, SAVE ::   ccntrENV_in,ccntrENV_omp
359    REAL, SAVE ::   coefcoli_in,coefcoli_omp
360
361    LOGICAL,SAVE  :: convscav_omp,convscav
362!$OMP THREADPRIVATE(iflag_lscav)
363!$OMP THREADPRIVATE(ccntrAA_in,ccntrENV_in,coefcoli_in)
364!$OMP THREADPRIVATE(convscav)
365    !RomP <<<
366    !######################################################################
367    !                    -- INITIALIZATION --
368    !######################################################################
369
370    DO k=1,klev
371       DO i=1,klon
372          sourceBE(i,k)=0.
373          Mint(i,k)=0.
374          zrho(i,k)=0.
375          zdz(i,k)=0.
376       END DO
377    END DO
378
379    DO it=1, nbtr
380       DO k=1,klev
381          DO i=1,klon
382             d_tr_insc(i,k,it)=0.
383             d_tr_bcscav(i,k,it)=0.
384             d_tr_evapls(i,k,it)=0.
385             d_tr_ls(i,k,it)=0.
386             d_tr_cv(i,k,it)=0.
387             d_tr_cl(i,k,it)=0.
388             d_tr_trsp(i,k,it)=0.
389             d_tr_sscav(i,k,it)=0.
390             d_tr_sat(i,k,it)=0.
391             d_tr_uscav(i,k,it)=0.
392             d_tr_lessi_impa(i,k,it)=0.
393             d_tr_lessi_nucl(i,k,it)=0.
394             qDi(i,k,it)=0.
395             qPr(i,k,it)=0.
396             qPa(i,k,it)=0.
397             qMel(i,k,it)=0.
398             qTrdi(i,k,it)=0.
399             dtrcvMA(i,k,it)=0.
400             zmfd1a(i,k,it)=0.
401             zmfdam(i,k,it)=0.
402             zmfphi2(i,k,it)=0.
403          END DO
404       END DO
405    END DO
406
407    DO it=1, nbtr
408       DO i=1,klon
409          d_tr_dry(i,it)=0.
410          flux_tr_dry(i,it)=0.
411       END DO
412    END DO
413
414    DO k = 1, klev
415       DO i = 1, klon
416          delp(i,k) = paprs(i,k)-paprs(i,k+1)
417       END DO
418    END DO
419
420    IF (debutphy) THEN
421       !!jyg
422!$OMP BARRIER
423       ecrit_tra=86400. ! frequence de stokage en dur
424       ! obsolete car remplace par des ecritures dans phys_output_write
425       !RomP >>>
426       !
427       !Config Key  = convscav
428       !Config Desc = Convective scavenging switch: 0=off, 1=on.
429       !Config Def  = .FALSE.
430       !Config Help =
431       !
432!$OMP MASTER
433       convscav_omp=.FALSE.
434       call getin('convscav', convscav_omp)
435       iflag_vdf_trac_omp=1
436       call getin('iflag_vdf_trac', iflag_vdf_trac_omp)
437       iflag_con_trac_omp=1
438       call getin('iflag_con_trac', iflag_con_trac_omp)
439       iflag_the_trac_omp=1
440       call getin('iflag_the_trac', iflag_the_trac_omp)
441!$OMP END MASTER
442!$OMP BARRIER
443       convscav=convscav_omp
444       iflag_vdf_trac=iflag_vdf_trac_omp
445       iflag_con_trac=iflag_con_trac_omp
446       iflag_the_trac=iflag_the_trac_omp
447       write(lunout,*) 'phytrac passage dans routine conv avec lessivage', convscav
448       !
449       !Config Key  = iflag_lscav
450       !Config Desc = Large scale scavenging parametrization: 0=none, 1=old(Genthon92),
451       !              2=1+PHeinrich, 3=Reddy_Boucher2004, 4=3+RPilon.
452       !Config Def  = 1
453       !Config Help =
454       !
455!$OMP MASTER
456       iflag_lscav_omp=1
457       call getin('iflag_lscav', iflag_lscav_omp)
458       ccntrAA_omp=1
459       ccntrENV_omp=1.
460       coefcoli_omp=0.001
461       call getin('ccntrAA', ccntrAA_omp)
462       call getin('ccntrENV', ccntrENV_omp)
463       call getin('coefcoli', coefcoli_omp)
464!$OMP END MASTER
465!$OMP BARRIER
466       iflag_lscav=iflag_lscav_omp
467       ccntrAA_in=ccntrAA_omp
468       ccntrENV_in=ccntrENV_omp
469       coefcoli_in=coefcoli_omp
470       !
471       SELECT CASE(iflag_lscav)
472       CASE(0)
473          WRITE(lunout,*)  'Large scale scavenging: none'
474       CASE(1)
475          WRITE(lunout,*)  'Large scale scavenging: C. Genthon, Tellus(1992), 44B, 371-389'
476       CASE(2)
477          WRITE(lunout,*)  'Large scale scavenging: C. Genthon, modified P. Heinrich'
478       CASE(3)
479          WRITE(lunout,*)  'Large scale scavenging: M. Shekkar Reddy and O. Boucher, JGR(2004), 109, D14202'
480       CASE(4)
481          WRITE(lunout,*)  'Large scale scavenging: Reddy and Boucher, modified R. Pilon'
482       END SELECT
483       !RomP <<<
484       WRITE(*,*) 'FIRST TIME IN PHYTRAC : pdtphys(sec) = ',pdtphys,'ecrit_tra (sec) = ',ecrit_tra
485       ALLOCATE( source(klon,nbtr), stat=ierr)
486       IF (ierr /= 0) CALL abort_physic('phytrac', 'pb in allocation 1',1)
487
488       ALLOCATE( aerosol(nbtr), stat=ierr)
489       IF (ierr /= 0) CALL abort_physic('phytrac', 'pb in allocation 2',1)
490
491
492       ! Initialize module for specific tracers
493       IF(type_trac == 'inca') THEN
494          source(:,:)=init_source(:,:)
495          CALL tracinca_init(aerosol,lessivage)
496       ELSE IF(type_trac == 'repr') THEN
497          source(:,:)=0.
498       ELSE IF(type_trac == 'co2i') THEN
499          source(:,:)=0.
500          lessivage  = .FALSE.
501          aerosol(:) = .FALSE.
502          pbl_flg(:) = 1
503          iflag_the_trac= 1
504          iflag_vdf_trac= 1
505          iflag_con_trac= 1
506       ELSE IF(type_trac == 'inco') THEN
507          source(:,1:nqCO2) = 0.                          ! from CO2i ModThL
508          source(:,nqCO2+1:nbtr)=init_source(:,:)         ! from INCA ModThL
509          aerosol(1:nqCO2) = .FALSE.                      ! from CO2i ModThL
510          CALL tracinca_init(aerosol(nqCO2+1:nbtr),lessivage)     ! from INCA ModThL
511          pbl_flg(1:nqCO2) = 1                            ! From CO2i ModThL
512          iflag_the_trac = 1                              ! From CO2i
513          iflag_vdf_trac = 1                              ! From CO2i
514          iflag_con_trac = 1                              ! From CO2i
515#ifdef CPP_StratAer
516       ELSE IF(type_trac == 'coag') THEN
517          source(:,:)=0.
518          DO it= 1, nbtr_sulgas
519            aerosol(it)=.FALSE.
520            IF (it==id_H2SO4_strat) aerosol(it)=.TRUE.
521          ENDDO
522          DO it= nbtr_sulgas+1, nbtr
523            aerosol(it)=.TRUE.
524          ENDDO
525#endif
526       ELSE IF(type_trac == 'lmdz') THEN
527          CALL traclmdz_init(pctsrf,xlat,xlon,ftsol,tr_seri,t_seri,pplay,sh,pdtphys,aerosol,lessivage)
528       END IF
529
530       !
531       !--initialising coefficients for scavenging in the case of NP
532       !
533       ALLOCATE(flag_cvltr(nbtr))
534       IF (iflag_con.EQ.3) THEN
535          !
536          ALLOCATE(ccntrAA(nbtr))
537          ALLOCATE(ccntrENV(nbtr))
538          ALLOCATE(coefcoli(nbtr))
539          !
540          DO it=1, nbtr
541             IF(type_trac == 'repr') THEN
542                 flag_cvltr(it)=.FALSE.
543             ELSE IF(type_trac == 'inca') THEN
544!                IF ((it.EQ.id_Rn222) .OR. ((it.GE.id_SO2) .AND. (it.LE.id_NH3)) ) THEN
545!                   !--gas-phase species
546!                   flag_cvltr(it)=.FALSE.
547!
548!                ELSEIF ( (it.GE.id_CIDUSTM) .AND. (it.LE.id_AIN) ) THEN
549!                   !--insoluble aerosol species
550!                   flag_cvltr(it)=.TRUE.
551!                   ccntrAA(it)=0.7
552!                   ccntrENV(it)=0.7
553!                   coefcoli(it)=0.001
554!                ELSEIF ( (it.EQ.id_Pb210) .OR. ((it.GE.id_CSSSM) .AND. (it.LE.id_SSN))) THEN
555!                   !--soluble aerosol species
556!                   flag_cvltr(it)=.TRUE.
557!                   ccntrAA(it)=0.9
558!                   ccntrENV(it)=0.9
559!                   coefcoli(it)=0.001
560!                ELSE
561!                   WRITE(lunout,*) 'pb it=', it
562!                   CALL abort_physic('phytrac','pb it scavenging',1)
563!                ENDIF
564                !--test OB
565                !--for now we do not scavenge in cvltr
566                flag_cvltr(it)=.FALSE.
567             ELSE IF(type_trac == 'co2i') THEN
568                !--co2 tracers are not scavenged
569                flag_cvltr(it)=.FALSE.
570             ELSE IF(type_trac == 'inco') THEN     ! Add ThL
571                flag_cvltr(it)=.FALSE.
572#ifdef CPP_StratAer
573             ELSE IF(type_trac == 'coag') THEN
574                IF (convscav.and.aerosol(it)) THEN
575                   flag_cvltr(it)=.TRUE.
576                   ccntrAA(it) =ccntrAA_in   
577                   ccntrENV(it)=ccntrENV_in
578                   coefcoli(it)=coefcoli_in
579                ELSE
580                   flag_cvltr(it)=.FALSE.
581                ENDIF
582#endif
583             ELSE IF(type_trac == 'lmdz') THEN
584                IF (convscav.and.aerosol(it)) THEN
585                   flag_cvltr(it)=.TRUE.
586                   ccntrAA(it) =ccntrAA_in    !--a modifier par JYG a lire depuis fichier
587                   ccntrENV(it)=ccntrENV_in
588                   coefcoli(it)=coefcoli_in
589                ELSE
590                   flag_cvltr(it)=.FALSE.
591                ENDIF
592             END IF
593          ENDDO
594          !
595       ELSE ! iflag_con .ne. 3
596          flag_cvltr(:) = .FALSE.
597       ENDIF
598       !
599       ! print out all tracer flags
600       !
601       WRITE(lunout,*) 'print out all tracer flags'
602       WRITE(lunout,*) 'type_trac      =', type_trac
603       WRITE(lunout,*) 'config_inca    =', config_inca
604       WRITE(lunout,*) 'iflag_con_trac =', iflag_con_trac
605       WRITE(lunout,*) 'iflag_con      =', iflag_con
606       WRITE(lunout,*) 'convscav       =', convscav
607       WRITE(lunout,*) 'iflag_lscav    =', iflag_lscav
608       WRITE(lunout,*) 'aerosol        =', aerosol
609       WRITE(lunout,*) 'iflag_the_trac =', iflag_the_trac
610       WRITE(lunout,*) 'iflag_thermals =', iflag_thermals
611       WRITE(lunout,*) 'iflag_vdf_trac =', iflag_vdf_trac
612       WRITE(lunout,*) 'pbl_flg        =', pbl_flg
613       WRITE(lunout,*) 'lessivage      =', lessivage
614       write(lunout,*)  'flag_cvltr    = ', flag_cvltr
615
616       IF (lessivage .AND. ANY(type_trac == ['inca','inco'])) &
617          CALL abort_physic('phytrac', 'lessivage=T config_inca=inca impossible',1)
618       !
619    ENDIF ! of IF (debutphy)
620    !############################################ END INITIALIZATION #######
621
622    DO k=1,klev
623       DO i=1,klon
624          zmasse(i,k)=(paprs(i,k)-paprs(i,k+1))/rg
625       END DO
626    END DO
627    !
628    IF (id_be .GT. 0) THEN
629       DO k=1,klev
630          DO i=1,klon
631             sourceBE(i,k)=srcbe(i,k)       !RomP  -> pour sortie histrac
632          END DO
633       END DO
634    ENDIF
635
636    !===============================================================================
637    !    -- Do specific treatment according to chemestry model or local LMDZ tracers
638    !     
639    !===============================================================================
640    IF(type_trac == 'inca') THEN
641       !    -- CHIMIE INCA  config_inca = aero or chem --
642       ! Appel fait en fin de phytrac pour avoir les emissions modifiees par
643       ! la couche limite et la convection avant le calcul de la chimie
644
645    ELSE IF(type_trac == 'repr') THEN
646       !   -- CHIMIE REPROBUS --
647       CALL tracreprobus(pdtphys, gmtime, debutphy, julien, &
648            presnivs, xlat, xlon, pphis, pphi, &
649            t_seri, pplay, paprs, sh , &
650            tr_seri)
651
652    ELSE IF(type_trac == 'co2i') THEN
653       !   -- CO2 interactif --
654       !   -- source is updated with FF and BB emissions
655       !   -- and net fluxes from ocean and orchidee
656       !   -- sign convention : positive into the atmosphere
657
658       CALL tracco2i(pdtphys, debutphy, &
659            xlat, xlon, pphis, pphi, &
660            t_seri, pplay, paprs, tr_seri, source)
661    ELSE IF(type_trac == 'inco') THEN      ! Add ThL
662       CALL tracco2i(pdtphys, debutphy, &
663            xlat, xlon, pphis, pphi, &
664            t_seri, pplay, paprs, tr_seri, source)
665
666#ifdef CPP_StratAer
667    ELSE IF(type_trac == 'coag') THEN
668       !   --STRATOSPHERIC AER IN THE STRAT --
669       CALL traccoag(pdtphys, gmtime, debutphy, julien, &
670            presnivs, xlat, xlon, pphis, pphi, &
671            t_seri, pplay, paprs, sh, rh , &
672            tr_seri)
673#endif
674    ELSE IF(type_trac == 'lmdz') THEN
675       !    -- Traitement des traceurs avec traclmdz
676       CALL traclmdz(nstep, julien, gmtime, pdtphys, t_seri, paprs, pplay, &
677            cdragh,  coefh, yu1, yv1, ftsol, pctsrf, xlat, xlon,iflag_vdf_trac>=0,sh, &
678            rh, pphi, ustar, wstar, ale_bl, ale_wake,  u10m, v10m, &
679            tr_seri, source, d_tr_cl,d_tr_dec, zmasse)               !RomP
680    END IF
681    !======================================================================
682    !       -- Calcul de l'effet de la convection --
683    !======================================================================
684
685    IF (iflag_con_trac==1) THEN
686
687       DO it=1, nbtr
688          IF ( conv_flg(it) == 0 ) CYCLE
689          IF (iflag_con.LT.2) THEN
690             !--pas de transport convectif
691             d_tr_cv(:,:,it)=0.
692
693          ELSE IF (iflag_con.EQ.2) THEN
694             !--ancien transport convectif de Tiedtke
695
696             CALL nflxtr(pdtphys, pmfu, pmfd, pen_u, pde_u, pen_d, pde_d, &
697                  pplay, paprs, tr_seri(:,:,it), d_tr_cv(:,:,it))
698          ELSE   
699             !--nouveau transport convectif de Emanuel
700
701             IF (flag_cvltr(it)) THEN
702                !--nouveau transport convectif de Emanuel avec lessivage convectif
703                !
704                !
705                ccntrAA_3d(:,:) =ccntrAA(it)
706                ccntrENV_3d(:,:)=ccntrENV(it)
707                coefcoli_3d(:,:)=coefcoli(it)
708
709                !--beware this interface is a bit weird because it is called for each tracer
710                !--with the full array tr_seri even if only item it is processed
711
712                CALL cvltr_scav(pdtphys, da, phi,phi2,d1a,dam, mp,ep,         &
713                     sigd,sij,wght_cvfd,clw,elij,epmlmMm,eplaMm,              &     
714                     pmflxr,pmflxs,evap,t_seri,wdtrainA,wdtrainM,             &   
715                     paprs,it,tr_seri,upwd,dnwd,itop_con,ibas_con,            &
716                     ccntrAA_3d,ccntrENV_3d,coefcoli_3d,                      &
717                     d_tr_cv,d_tr_trsp,d_tr_sscav,d_tr_sat,d_tr_uscav,qDi,qPr,&
718                     qPa,qMel,qTrdi,dtrcvMA,Mint,                             &
719                     zmfd1a,zmfphi2,zmfdam)
720
721
722             ELSE  !---flag_cvltr(it).EQ.FALSE
723                !--nouveau transport convectif de Emanuel mais pas de lessivage convectif
724
725                !--beware this interface is a bit weird because it is called for each tracer
726                !--with the full array tr_seri even if only item it is processed
727                !
728                CALL cvltr_noscav(it,pdtphys, da, phi,mp,wght_cvfd,paprs,pplay, &  !jyg
729                     tr_seri,upwd,dnwd,d_tr_cv)                                !jyg
730
731             ENDIF
732
733          ENDIF !--iflag
734
735          !--on ajoute les tendances
736
737          DO k = 1, klev
738             DO i = 1, klon       
739                tr_seri(i,k,it) = tr_seri(i,k,it) + d_tr_cv(i,k,it)
740             END DO
741          END DO
742
743          CALL minmaxqfi(tr_seri(:,:,it),0.,1.e33,'convection it = '//TRIM(int2str(it)))
744
745       END DO ! nbtr
746
747#ifdef CPP_StratAer
748       IF (type_trac=='coag') THEN
749         ! initialize wet deposition flux of sulfur
750         budg_dep_wet_ocs(:)=0.0
751         budg_dep_wet_so2(:)=0.0
752         budg_dep_wet_h2so4(:)=0.0
753         budg_dep_wet_part(:)=0.0
754         ! compute wet deposition flux of sulfur (sum over gases and particles)
755         ! and convert to kg(S)/m2/s
756         DO i = 1, klon
757         DO k = 1, klev
758         DO it = 1, nbtr
759         !do not include SO2 because most of it comes trom the troposphere
760           IF (it==id_OCS_strat) THEN
761             budg_dep_wet_ocs(i)=budg_dep_wet_ocs(i)+d_tr_cv(i,k,it)*(mSatom/mOCSmol) &
762                            & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
763           ELSEIF (it==id_SO2_strat) THEN
764             budg_dep_wet_so2(i)=budg_dep_wet_so2(i)+d_tr_cv(i,k,it)*(mSatom/mSO2mol) &
765                            & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
766           ELSEIF (it==id_H2SO4_strat) THEN
767             budg_dep_wet_h2so4(i)=budg_dep_wet_h2so4(i)+d_tr_cv(i,k,it)*(mSatom/mH2SO4mol) &
768                            & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
769           ELSEIF (it.GT.nbtr_sulgas) THEN
770             budg_dep_wet_part(i)=budg_dep_wet_part(i)+d_tr_cv(i,k,it)*(mSatom/mH2SO4mol)  &
771                            & *dens_aer_dry*4./3.*RPI*(mdw(it-nbtr_sulgas)/2.)**3 &
772                            & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
773           ENDIF
774         ENDDO
775         ENDDO
776         ENDDO
777       ENDIF
778#endif
779
780    ENDIF ! convection
781
782    !======================================================================
783    !    -- Calcul de l'effet des thermiques --
784    !======================================================================
785
786    DO it=1,nbtr
787       DO k=1,klev
788          DO i=1,klon
789             d_tr_th(i,k,it)=0.
790             tr_seri(i,k,it)=MAX(tr_seri(i,k,it),0.)
791! the next safeguard causes some problem for stratospheric aerosol tracers (particle number)
792! and there is little justification for it so it is commented out (4 December 2017) by OB
793! if reinstated please keep the ifndef CPP_StratAer
794!#ifndef CPP_StratAer
795!             tr_seri(i,k,it)=MIN(tr_seri(i,k,it),1.e10)
796!#endif
797          END DO
798       END DO
799    END DO
800
801    IF (iflag_thermals.GT.0.AND.iflag_the_trac>0) THEN   
802
803       DO it=1, nbtr
804
805          CALL thermcell_dq(klon,klev,1,pdtphys,fm_therm,entr_therm, &
806               zmasse,tr_seri(1:klon,1:klev,it),        &
807               d_tr_th(1:klon,1:klev,it),ztra_th,0 )
808
809          DO k=1,klev
810             DO i=1,klon
811                d_tr_th(i,k,it)=pdtphys*d_tr_th(i,k,it)
812                tr_seri(i,k,it)=MAX(tr_seri(i,k,it)+d_tr_th(i,k,it),0.)
813             END DO
814          END DO
815
816       END DO ! it
817
818    ENDIF ! Thermiques
819
820    !======================================================================
821    !     -- Calcul de l'effet de la couche limite --
822    !======================================================================
823
824    IF (iflag_vdf_trac==1) THEN
825
826       !  Injection during BL mixing
827       !
828#ifdef CPP_StratAer
829       IF (type_trac=='coag') THEN
830
831         ! initialize dry deposition flux of sulfur
832         budg_dep_dry_ocs(:)=0.0
833         budg_dep_dry_so2(:)=0.0
834         budg_dep_dry_h2so4(:)=0.0
835         budg_dep_dry_part(:)=0.0
836
837         ! compute dry deposition velocity as function of surface type (numbers
838         ! from IPSL note 23, 2002)
839         v_dep_dry(:) =  pctsrf(:,is_ter) * 2.5e-3 &
840                     & + pctsrf(:,is_oce) * 0.5e-3 &
841                     & + pctsrf(:,is_lic) * 2.5e-3 &
842                     & + pctsrf(:,is_sic) * 2.5e-3
843
844         ! compute surface dry deposition flux
845         zrho(:,1)=pplay(:,1)/t_seri(:,1)/RD
846
847         DO it=1, nbtr
848          source(:,it) = - v_dep_dry(:) * tr_seri(:,1,it) * zrho(:,1)
849         ENDDO
850
851       ENDIF
852#endif
853
854       DO it=1, nbtr
855          !
856          IF( pbl_flg(it) /= 0 ) THEN
857             !
858             CALL cltrac(pdtphys, coefh,t_seri,       &
859                  tr_seri(:,:,it), source(:,it),      &
860                  paprs, pplay, delp,                 &
861                  d_tr_cl(:,:,it),d_tr_dry(:,it),flux_tr_dry(:,it))
862             !
863             tr_seri(:,:,it)=tr_seri(:,:,it)+d_tr_cl(:,:,it)
864             !
865#ifdef CPP_StratAer
866             IF (type_trac=='coag') THEN
867               ! compute dry deposition flux of sulfur (sum over gases and particles)
868               IF (it==id_OCS_strat) THEN
869                 budg_dep_dry_ocs(:)=budg_dep_dry_ocs(:)-source(:,it)*(mSatom/mOCSmol)
870               ELSEIF (it==id_SO2_strat) THEN
871                 budg_dep_dry_so2(:)=budg_dep_dry_so2(:)-source(:,it)*(mSatom/mSO2mol)
872               ELSEIF (it==id_H2SO4_strat) THEN
873                 budg_dep_dry_h2so4(:)=budg_dep_dry_h2so4(:)-source(:,it)*(mSatom/mH2SO4mol)
874               ELSEIF (it.GT.nbtr_sulgas) THEN
875                 budg_dep_dry_part(:)=budg_dep_dry_part(:)-source(:,it)*(mSatom/mH2SO4mol)*dens_aer_dry &
876                                & *4./3.*RPI*(mdw(it-nbtr_sulgas)/2.)**3
877               ENDIF
878             ENDIF
879#endif
880             !
881          ENDIF
882          !
883       ENDDO
884       !
885    ELSE IF (iflag_vdf_trac==0) THEN
886       !
887       !   Injection of source in the first model layer
888       !
889       DO it=1,nbtr
890          d_tr_cl(:,1,it)=source(:,it)*RG/delp(:,1)*pdtphys
891          tr_seri(:,1,it)=tr_seri(:,1,it)+d_tr_cl(:,1,it)
892       ENDDO
893       d_tr_cl(:,2:klev,1:nbtr)=0.
894       !
895    ELSE IF (iflag_vdf_trac==-1) THEN
896       !
897       ! Nothing happens
898       d_tr_cl=0.
899       !
900    ELSE
901       !
902       CALL abort_physic('iflag_vdf_trac', 'cas non prevu',1)
903       !
904    ENDIF ! couche limite
905
906    !======================================================================
907    !   Calcul de l'effet de la precipitation grande echelle
908    !   POUR INCA le lessivage est fait directement dans INCA
909    !======================================================================
910
911    IF (lessivage) THEN
912
913       ql_incloud_ref = 10.e-4
914       ql_incloud_ref =  5.e-4
915
916
917       ! calcul du contenu en eau liquide au sein du nuage
918       ql_incl = ql_incloud_ref
919       ! choix du lessivage
920       !
921       IF (iflag_lscav .EQ. 3 .OR. iflag_lscav .EQ. 4) THEN
922          ! ********  Olivier Boucher version (3) possibly with modified ql_incl (4)
923          !
924          DO it = 1, nbtr
925
926             IF (aerosol(it)) THEN
927             !  incloud scavenging and removal by large scale rain ! orig : ql_incl was replaced by 0.5e-3 kg/kg
928             ! the value 0.5e-3 kg/kg is from Giorgi and Chameides (1986), JGR
929             ! Liu (2001) proposed to use 1.5e-3 kg/kg
930
931!jyg<
932!!             CALL lsc_scav(pdtphys,it,iflag_lscav,ql_incl,prfl,psfl,rneb,beta_fisrt,  &
933             CALL lsc_scav(pdtphys,it,iflag_lscav,aerosol,ql_incl,prfl,psfl,rneb,beta_fisrt,  &
934!>jyg
935                  beta_v1,pplay,paprs,t_seri,tr_seri,d_tr_insc,   &
936                  d_tr_bcscav,d_tr_evapls,qPrls)
937
938             !large scale scavenging tendency
939             DO k = 1, klev
940                DO i = 1, klon
941                   d_tr_ls(i,k,it)=d_tr_insc(i,k,it)+d_tr_bcscav(i,k,it)+d_tr_evapls(i,k,it)
942                   tr_seri(i,k,it)=tr_seri(i,k,it)+d_tr_ls(i,k,it)
943                ENDDO
944             ENDDO
945             CALL minmaxqfi(tr_seri(:,:,it),0.,1.e33,'lsc scav it = '//TRIM(int2str(it)))
946             ENDIF
947
948          END DO  !tr
949
950#ifdef CPP_StratAer
951         IF (type_trac=='coag') THEN
952           ! compute wet deposition flux of sulfur (sum over gases and
953           ! particles) and convert to kg(S)/m2/s
954           ! adding contribution of d_tr_ls to d_tr_cv (above)
955           DO i = 1, klon
956           DO k = 1, klev
957           DO it = 1, nbtr
958             IF (it==id_OCS_strat) THEN
959               budg_dep_wet_ocs(i)=budg_dep_wet_ocs(i)+d_tr_ls(i,k,it)*(mSatom/mOCSmol) &
960                              & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
961             ELSEIF (it==id_SO2_strat) THEN
962               budg_dep_wet_so2(i)=budg_dep_wet_so2(i)+d_tr_ls(i,k,it)*(mSatom/mSO2mol) &
963                              & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
964             ELSEIF (it==id_H2SO4_strat) THEN
965               budg_dep_wet_h2so4(i)=budg_dep_wet_h2so4(i)+d_tr_ls(i,k,it)*(mSatom/mH2SO4mol) &
966                              & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
967             ELSEIF (it.GT.nbtr_sulgas) THEN
968               budg_dep_wet_part(i)=budg_dep_wet_part(i)+d_tr_ls(i,k,it)*(mSatom/mH2SO4mol)  &
969                              & *dens_aer_dry*4./3.*RPI*(mdw(it-nbtr_sulgas)/2.)**3 &
970                              & *(paprs(i,k)-paprs(i,k+1))/RG/pdtphys
971             ENDIF
972           ENDDO
973           ENDDO
974           ENDDO
975         ENDIF
976#endif
977
978       ELSE IF (iflag_lscav .EQ. 2) THEN ! frac_impa, frac_nucl
979          ! *********   modified  old version
980
981          d_tr_lessi_nucl(:,:,:) = 0.
982          d_tr_lessi_impa(:,:,:) = 0.
983          flestottr(:,:,:) = 0.
984          ! Tendance des aerosols nuclees et impactes
985          DO it = 1, nbtr
986             IF (aerosol(it)) THEN
987                his_dh(:)=0.
988                DO k = 1, klev
989                   DO i = 1, klon
990                      !PhH
991                      zrho(i,k)=pplay(i,k)/t_seri(i,k)/RD
992                      zdz(i,k)=(paprs(i,k)-paprs(i,k+1))/zrho(i,k)/RG
993                      !
994                   ENDDO
995                ENDDO
996
997                DO k=klev-1, 1, -1
998                   DO i=1, klon
999                      !             d_tr_ls(i,k,it)=tr_seri(i,k,it)*(frac_impa(i,k)*frac_nucl(i,k)-1.)
1000                      dx=d_tr_ls(i,k,it)
1001                      his_dh(i)=his_dh(i)-dx*zrho(i,k)*zdz(i,k)/pdtphys !  kg/m2/s
1002                      evaplsc = prfl(i,k) - prfl(i,k+1) + psfl(i,k) - psfl(i,k+1)
1003                      ! Evaporation Partielle -> Liberation Partielle 0.5*evap
1004                      IF ( evaplsc .LT.0..and.abs(prfl(i,k+1)+psfl(i,k+1)).gt.1.e-10) THEN
1005                         evaplsc = (-evaplsc)/(prfl(i,k+1)+psfl(i,k+1))
1006                         ! evaplsc est donc positif, his_dh(i) est positif
1007                         !-------------- 
1008                         d_tr_evapls(i,k,it)=0.5*evaplsc*(d_tr_lessi_nucl(i,k+1,it) &
1009                              +d_tr_lessi_impa(i,k+1,it))
1010                         !-------------   d_tr_evapls(i,k,it)=-0.5*evaplsc*(d_tr_lsc(i,k+1,it))
1011                         beta=0.5*evaplsc
1012                         if ((prfl(i,k)+psfl(i,k)).lt.1.e-10) THEN
1013                            beta=1.0*evaplsc
1014                         endif
1015                         dx=beta*his_dh(i)/zrho(i,k)/zdz(i,k)*pdtphys
1016                         his_dh(i)=(1.-beta)*his_dh(i)   ! tracer from
1017                         d_tr_evapls(i,k,it)=dx
1018                      ENDIF
1019                      d_tr_ls(i,k,it)=tr_seri(i,k,it)*(frac_impa(i,k)*frac_nucl(i,k)-1.) &
1020                           +d_tr_evapls(i,k,it)
1021
1022                      !--------------
1023                      d_tr_lessi_nucl(i,k,it) = d_tr_lessi_nucl(i,k,it) +    &
1024                           ( 1 - frac_nucl(i,k) )*tr_seri(i,k,it)
1025                      d_tr_lessi_impa(i,k,it) = d_tr_lessi_impa(i,k,it) +    &
1026                           ( 1 - frac_impa(i,k) )*tr_seri(i,k,it)
1027                      !
1028                      ! Flux lessivage total
1029                      flestottr(i,k,it) = flestottr(i,k,it) -   &
1030                           ( d_tr_lessi_nucl(i,k,it)   +        &
1031                           d_tr_lessi_impa(i,k,it) ) *          &
1032                           ( paprs(i,k)-paprs(i,k+1) ) /        &
1033                           (RG * pdtphys)
1034                      !! Mise a jour des traceurs due a l'impaction,nucleation
1035                      !                 tr_seri(i,k,it)=tr_seri(i,k,it)*frac_impa(i,k)*frac_nucl(i,k)
1036                      !!  calcul de la tendance liee au lessivage stratiforme
1037                      !                 d_tr_ls(i,k,it)=tr_seri(i,k,it)*&
1038                      !                                (1.-1./(frac_impa(i,k)*frac_nucl(i,k)))
1039                      !--------------
1040                   ENDDO
1041                ENDDO
1042             ENDIF
1043          ENDDO
1044          ! *********   end modified old version
1045
1046       ELSE IF (iflag_lscav .EQ. 1) THEN ! frac_impa, frac_nucl
1047          ! *********    old version
1048
1049          d_tr_lessi_nucl(:,:,:) = 0.
1050          d_tr_lessi_impa(:,:,:) = 0.
1051          flestottr(:,:,:) = 0.
1052          !=========================
1053          ! LESSIVAGE LARGE SCALE :
1054          !=========================
1055
1056          ! Tendance des aerosols nuclees et impactes
1057          ! -----------------------------------------
1058          DO it = 1, nbtr
1059             IF (aerosol(it)) THEN
1060                DO k = 1, klev
1061                   DO i = 1, klon
1062                      d_tr_lessi_nucl(i,k,it) = d_tr_lessi_nucl(i,k,it) +    &
1063                           ( 1 - frac_nucl(i,k) )*tr_seri(i,k,it)
1064                      d_tr_lessi_impa(i,k,it) = d_tr_lessi_impa(i,k,it) +    &
1065                           ( 1 - frac_impa(i,k) )*tr_seri(i,k,it)
1066
1067                      !
1068                      ! Flux lessivage total
1069                      ! ------------------------------------------------------------
1070                      flestottr(i,k,it) = flestottr(i,k,it) -   &
1071                           ( d_tr_lessi_nucl(i,k,it)   +        &
1072                           d_tr_lessi_impa(i,k,it) ) *          &
1073                           ( paprs(i,k)-paprs(i,k+1) ) /        &
1074                           (RG * pdtphys)
1075                      !
1076                      ! Mise a jour des traceurs due a l'impaction,nucleation
1077                      ! ----------------------------------------------------------------------
1078                      tr_seri(i,k,it)=tr_seri(i,k,it)*frac_impa(i,k)*frac_nucl(i,k)
1079                   ENDDO
1080                ENDDO
1081             ENDIF
1082          ENDDO
1083
1084          ! *********   end old version
1085       ENDIF  !  iflag_lscav . EQ. 1, 2, 3 or 4
1086       !
1087    ENDIF !  lessivage
1088
1089
1090    !    -- CHIMIE INCA  config_inca = aero or chem --
1091    IF (ANY(type_trac == ['inca','inco'])) THEN  ! ModThL
1092
1093       CALL tracinca(&
1094            nstep,    julien,   gmtime,         lafin,     &
1095            pdtphys,  t_seri,   paprs,          pplay,     &
1096            pmfu,     upwd,     ftsol,  pctsrf, pphis,     &
1097            pphi,     albsol,   sh,    ch,     rh,        &
1098            cldfra,   rneb,     diafra,         cldliq,    &
1099            itop_con, ibas_con, pmflxr,         pmflxs,    &
1100            prfl,     psfl,     aerosol_couple, flxmass_w, &
1101            tau_aero, piz_aero, cg_aero,        ccm,       &
1102            rfname,                                        &
1103            tr_seri(:,:,1+nqCO2:nbtr),  source(:,1+nqCO2:nbtr))  ! ModThL 
1104    ENDIF
1105
1106  END SUBROUTINE phytrac
1107
1108END MODULE
Note: See TracBrowser for help on using the repository browser.