source: LMDZ5/trunk/libf/phylmd/StratAer/traccoag_mod.F90 @ 2709

Last change on this file since 2709 was 2704, checked in by oboucher, 8 years ago

This revision concerns the StratAer? module and should not impact the rest of LMDz
Bug correction in interp_sulf_input.F90
Update of miecalc_aer.F90 and traccoag_mod.F90
Phytrac tracers are now dealt with in XIOS through the Fortran interface with minimal input in the xml
Making tracer groups in DefLists? for StratAer?

File size: 13.3 KB
Line 
1MODULE traccoag_mod
2!
3! This module calculates the concentration of aerosol particles in certain size bins
4! considering coagulation and sedimentation.
5!
6CONTAINS
7
8  SUBROUTINE traccoag(pdtphys, gmtime, debutphy, julien, &
9       presnivs, xlat, xlon, pphis, pphi, &
10       t_seri, pplay, paprs, sh, rh , &
11       tr_seri)
12
13    USE phys_local_var_mod, ONLY: mdw, sulf_convert, sulf_nucl, sulf_cond_evap, &
14        & sfluxaer, ocs_convert, R2SO4, DENSO4, f_r_wet, SO2_backgr_tend, OCS_backgr_tend, &
15        & OCS_lifetime, SO2_lifetime, surf_PM25_sulf
16
17    USE dimphy
18    USE infotrac
19    USE aerophys
20    USE geometry_mod, ONLY : cell_area
21    USE mod_grid_phy_lmdz
22    USE mod_phys_lmdz_mpi_data, ONLY :  is_mpi_root
23    USE mod_phys_lmdz_para, only: gather, scatter
24    USE phys_cal_mod
25    USE sulfate_aer_mod
26    USE phys_local_var_mod, ONLY: stratomask
27    USE YOMCST
28
29    IMPLICIT NONE
30
31! Input argument
32!---------------
33    REAL,INTENT(IN)    :: pdtphys    ! Pas d'integration pour la physique (seconde)
34    REAL,INTENT(IN)    :: gmtime     ! Heure courante
35    LOGICAL,INTENT(IN) :: debutphy   ! le flag de l'initialisation de la physique
36    INTEGER,INTENT(IN) :: julien     ! Jour julien
37
38    REAL,DIMENSION(klev),INTENT(IN)        :: presnivs! pressions approximat. des milieux couches (en PA)
39    REAL,DIMENSION(klon),INTENT(IN)        :: xlat    ! latitudes pour chaque point
40    REAL,DIMENSION(klon),INTENT(IN)        :: xlon    ! longitudes pour chaque point
41    REAL,DIMENSION(klon),INTENT(IN)        :: pphis   ! geopotentiel du sol
42    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pphi    ! geopotentiel de chaque couche
43
44    REAL,DIMENSION(klon,klev),INTENT(IN)   :: t_seri  ! Temperature
45    REAL,DIMENSION(klon,klev),INTENT(IN)   :: pplay   ! pression pour le mileu de chaque couche (en Pa)
46    REAL,DIMENSION(klon,klev+1),INTENT(IN) :: paprs   ! pression pour chaque inter-couche (en Pa)
47    REAL,DIMENSION(klon,klev),INTENT(IN)   :: sh      ! humidite specifique
48    REAL,DIMENSION(klon,klev),INTENT(IN)   :: rh      ! humidite relative   
49
50! Output argument
51!----------------
52    REAL,DIMENSION(klon,klev,nbtr),INTENT(INOUT)  :: tr_seri ! Concentration Traceur [U/KgA] 
53
54! Local variables
55!----------------
56! flag for sulfur emission scenario: (0) background aerosol ; (1) volcanic eruption ; (2) stratospheric aerosol injections (SAI)
57    INTEGER,PARAMETER  :: flag_sulf_emit=2
58!
59!--flag_sulf_emit=1 --example Pinatubo
60    INTEGER,PARAMETER :: year_emit_vol=1991          ! year of emission date
61    INTEGER,PARAMETER :: mth_emit_vol=6              ! month of emission date
62    INTEGER,PARAMETER :: day_emit_vol=15             ! day of emission date
63    REAL,PARAMETER    :: m_aer_emiss_vol=7.e9   ! emitted sulfur mass in kgS, e.g. 7Tg(S)=14Tg(SO2)
64    REAL,PARAMETER    :: altemiss_vol=17.e3     ! emission altitude in m
65    REAL,PARAMETER    :: sigma_alt_vol=1.e3     ! standard deviation of emission altitude in m
66    REAL,PARAMETER    :: xlat_vol=15.14         ! latitude of volcano in degree
67    REAL,PARAMETER    :: xlon_vol=120.35        ! longitude of volcano in degree
68
69!--flag_sulf_emit=2 --SAI
70    REAL,PARAMETER    :: m_aer_emiss_sai=1.e10  ! emitted sulfur mass in kgS, eg 1e9=1TgS, 1e10=10TgS
71    REAL,PARAMETER    :: altemiss_sai=17.e3     ! emission altitude in m
72    REAL,PARAMETER    :: sigma_alt_sai=1.e3     ! standard deviation of emission altitude in m
73    REAL,PARAMETER    :: xlat_sai=0.0           ! latitude of SAI in degree
74    REAL,PARAMETER    :: xlon_sai=120.35        ! longitude of SAI in degree
75
76!--other local variables
77    INTEGER                                :: it, k, i, ilon, ilev, itime, i_int
78    LOGICAL,DIMENSION(klon,klev)           :: is_strato           ! true = above tropopause, false = below
79    REAL,DIMENSION(klon,klev)              :: m_air_gridbox       ! mass of air in every grid box [kg]
80    REAL,DIMENSION(klon_glo,klev,nbtr)     :: tr_seri_glo         ! Concentration Traceur [U/KgA] 
81    REAL,DIMENSION(klev+1)                 :: altLMDz             ! altitude of layer interfaces in m
82    REAL,DIMENSION(klev)                   :: f_lay_emiss         ! fraction of emission for every vertical layer
83    REAL                                   :: f_lay_sum           ! sum of layer emission fractions
84    REAL                                   :: alt                 ! altitude for integral calculation
85    INTEGER,PARAMETER                      :: n_int_alt=10        ! number of subintervals for integration over Gaussian emission profile
86    REAL,DIMENSION(nbtr_bin)               :: r_bin               ! particle radius in size bin [m]
87    REAL,DIMENSION(nbtr_bin)               :: r_lower             ! particle radius at lower bin boundary [m]
88    REAL,DIMENSION(nbtr_bin)               :: r_upper             ! particle radius at upper bin boundary [m]
89    REAL,DIMENSION(nbtr_bin)               :: m_part_dry          ! mass of one dry particle in size bin [kg]
90    REAL                                   :: zrho                ! Density of air [kg/m3]
91    REAL                                   :: zdz                 ! thickness of atm. model layer in m
92    REAL,DIMENSION(klon,klev)              :: dens_aer            ! density of aerosol particles [kg/m3 aerosol] with default H2SO4 mass fraction
93    REAL                                   :: dlat, dlon          ! d latitude and d longitude of grid in degree
94
95    IF (is_mpi_root) THEN
96      PRINT *,'in traccoag: date from phys_cal_mod =',year_cur,'-',mth_cur,'-',day_cur,'-',hour
97    ENDIF
98
99    dlat=180./2./FLOAT(nbp_lat)   ! d latitude in degree
100    dlon=360./2./FLOAT(nbp_lon)   ! d longitude in degree
101
102    DO it=1, nbtr_bin
103      r_bin(it)=mdw(it)/2.
104    ENDDO
105
106!--set boundaries of size bins
107    DO it=1, nbtr_bin
108    IF (it.EQ.1) THEN
109      r_upper(it)=sqrt(r_bin(it+1)*r_bin(it))
110      r_lower(it)=r_bin(it)**2./r_upper(it)
111    ELSEIF (it.EQ.nbtr_bin) THEN
112      r_lower(it)=sqrt(r_bin(it)*r_bin(it-1))
113      r_upper(it)=r_bin(it)**2./r_lower(it)
114    ELSE
115      r_lower(it)=sqrt(r_bin(it)*r_bin(it-1))
116      r_upper(it)=sqrt(r_bin(it+1)*r_bin(it))
117    ENDIF
118    ENDDO
119
120    IF (debutphy .and. is_mpi_root) THEN
121      DO it=1, nbtr_bin
122        PRINT *,'radius bin', it, ':', r_bin(it), '(from',  r_lower(it), 'to', r_upper(it), ')'
123      ENDDO
124    ENDIF
125
126!--initialising logical is_strato from stratomask
127    is_strato(:,:)=.FALSE.
128    WHERE (stratomask.GT.0.5) is_strato=.TRUE.
129
130! STRACOMP (H2O, P, t_seri -> aerosol composition (R2SO4))
131! H2SO4 mass fraction in aerosol (%)
132    CALL stracomp(sh,t_seri,pplay)
133
134! aerosol density (gr/cm3)
135    CALL denh2sa(t_seri)
136
137! compute factor for converting dry to wet radius (for every grid box)
138    f_r_wet(:,:) = (dens_aer_dry/(DENSO4(:,:)*1000.)/(R2SO4(:,:)/100.))**(1./3.)
139
140!--calculate mass of air in every grid box
141    DO ilon=1, klon
142    DO ilev=1, klev
143      m_air_gridbox(ilon,ilev)=(paprs(ilon,ilev)-paprs(ilon,ilev+1)) / RG * cell_area(ilon)
144    ENDDO
145    ENDDO
146
147    IF (debutphy) THEN
148      CALL gather(tr_seri, tr_seri_glo)
149      IF (MAXVAL(tr_seri_glo).LT.1.e-30) THEN
150!--initialising tracer concentrations to zero
151        DO it=1, nbtr
152        tr_seri(:,:,it)=0.0
153        ENDDO
154      ENDIF
155    ENDIF
156
157!--sulfur emission, depending on chosen scenario (flag_sulf_emit)
158    SELECT CASE(flag_sulf_emit)
159
160    CASE(0) ! background aerosol
161      ! do nothing (no emission)
162
163    CASE(1) ! volcanic eruption
164      !--only emit on day of eruption
165      ! stretch emission over one day of Pinatubo eruption
166      IF (year_cur==year_emit_vol.AND.mth_cur==mth_emit_vol.AND.day_cur==day_emit_vol) THEN
167!
168        DO i=1,klon
169          !Pinatubo eruption at 15.14N, 120.35E
170          IF  ( xlat(i).GE.xlat_vol-dlat .AND. xlat(i).LT.xlat_vol+dlat .AND. &
171                xlon(i).GE.xlon_vol-dlon .AND. xlon(i).LT.xlon_vol+dlon ) THEN
172!         compute altLMDz
173            altLMDz(:)=0.0
174            DO k=1, klev
175              zrho=pplay(i,k)/t_seri(i,k)/RD            !air density in kg/m3
176              zdz=(paprs(i,k)-paprs(i,k+1))/zrho/RG     !thickness of layer in m
177              altLMDz(k+1)=altLMDz(k)+zdz
178            ENDDO
179            !compute distribution of emission to vertical model layers (based on Gaussian peak in altitude)
180            f_lay_sum=0.0
181            DO k=1, klev
182              f_lay_emiss(k)=0.0
183              DO i_int=1, n_int_alt
184                alt=altLMDz(k)+float(i_int)*(altLMDz(k+1)-altLMDz(k))/float(n_int_alt)
185                f_lay_emiss(k)=f_lay_emiss(k)+1./(sqrt(2.*RPI)*sigma_alt_vol)* &
186                &              exp(-0.5*((alt-altemiss_vol)/sigma_alt_vol)**2.)*   &
187                &              (altLMDz(k+1)-altLMDz(k))/float(n_int_alt)
188              ENDDO
189              f_lay_sum=f_lay_sum+f_lay_emiss(k)
190            ENDDO
191            !correct for step integration error
192            f_lay_emiss(:)=f_lay_emiss(:)/f_lay_sum
193            !emission as SO2 gas (with m(SO2)=64/32*m_aer_emiss)
194            !vertically distributed emission
195            DO k=1, klev
196              tr_seri(i,k,id_SO2_strat)=tr_seri(i,k,id_SO2_strat)+ &
197              & m_aer_emiss_vol*(mSO2mol/mSatom)/m_air_gridbox(i,k)*f_lay_emiss(k) &
198              & /(1.*86400./pdtphys) ! stretch emission over one day of Pinatubo eruption
199            ENDDO
200          ENDIF ! emission grid cell
201        ENDDO ! klon loop
202      ENDIF ! emission period
203
204    CASE(2) ! stratospheric aerosol injections (SAI)
205!
206      DO i=1,klon
207!       SAI standard scenario with continuous emission from 1 grid point at the equator
208!       SAI emission on single month
209!       IF  ((mth_cur==4 .AND. &
210!       SAI continuous emission o
211        IF  ( xlat(i).GE.xlat_sai-dlat .AND. xlat(i).LT.xlat_sai+dlat .AND. &
212          &   xlon(i).GE.xlon_sai-dlon .AND. xlon(i).LT.xlon_sai+dlon ) THEN
213!         compute altLMDz
214          altLMDz(:)=0.0
215          DO k=1, klev
216            zrho=pplay(i,k)/t_seri(i,k)/RD            !air density in kg/m3
217            zdz=(paprs(i,k)-paprs(i,k+1))/zrho/RG     !thickness of layer in m
218            altLMDz(k+1)=altLMDz(k)+zdz
219          ENDDO
220          !compute distribution of emission to vertical model layers (based on Gaussian peak in altitude)
221          f_lay_sum=0.0
222          DO k=1, klev
223            f_lay_emiss(k)=0.0
224            DO i_int=1, n_int_alt
225              alt=altLMDz(k)+float(i_int)*(altLMDz(k+1)-altLMDz(k))/float(n_int_alt)
226              f_lay_emiss(k)=f_lay_emiss(k)+1./(sqrt(2.*RPI)*sigma_alt_sai)* &
227              &              exp(-0.5*((alt-altemiss_sai)/sigma_alt_sai)**2.)*   &
228              &              (altLMDz(k+1)-altLMDz(k))/float(n_int_alt)
229            ENDDO
230            f_lay_sum=f_lay_sum+f_lay_emiss(k)
231          ENDDO
232          !correct for step integration error
233          f_lay_emiss(:)=f_lay_emiss(:)/f_lay_sum
234          !emission as SO2 gas (with m(SO2)=64/32*m_aer_emiss)
235          !vertically distributed emission
236          DO k=1, klev
237            tr_seri(i,k,id_SO2_strat)=tr_seri(i,k,id_SO2_strat)+ &
238            & m_aer_emiss_sai*(mSO2mol/mSatom)/m_air_gridbox(i,k)*f_lay_emiss(k) &
239            & /(360.*86400./pdtphys) ! stretch emission over whole year (360d)
240!            & /(60.*86400./pdtphys) ! stretch emission over 2 months (seasonal emission)
241!            & /7. ! distribute equally over 7 emission grid points
242          ENDDO
243!          !emission as monodisperse particles with 0.1um dry radius (BIN21)
244!          !vertically distributed emission
245!          DO k=1, klev
246!            tr_seri(i,k,id_BIN01_strat+20)=tr_seri(i,k,id_BIN01_strat+20)+ &
247!            & m_aer_emiss*(mH2SO4mol/mSatom)/m_part_dry(21)/m_air_gridbox(i,k)*f_lay_emiss(k) &
248!            & /(360.*86400./pdtphys) & ! stretch emission over whole year (360d)
249!            & /7. ! distribute equally over 7 emission grid points
250!          ENDDO
251        ENDIF ! emission grid cell
252      ENDDO ! klon loop
253
254    END SELECT ! emission scenario (flag_sulf_emit)
255
256!--read background concentrations of OCS and SO2 and lifetimes from input file
257!--update the variables defined in phys_local_var_mod
258    CALL interp_sulf_input(debutphy,pdtphys,paprs,tr_seri)
259
260!--convert OCS to SO2 in the stratosphere
261    CALL ocs_to_so2(pdtphys,tr_seri,t_seri,pplay,paprs,sh,is_strato)
262
263!--convert SO2 to H2SO4
264    CALL so2_to_h2so4(pdtphys,tr_seri,t_seri,pplay,paprs,sh,is_strato)
265
266!--common routine for nucleation and condensation/evaporation with adaptive timestep
267    CALL micphy_tstep(pdtphys,tr_seri,t_seri,pplay,paprs,rh,is_strato)
268
269!--call coagulation routine
270    CALL coagulate(pdtphys,mdw,tr_seri,t_seri,pplay,dens_aer,is_strato)
271
272!--call sedimentation routine
273    CALL aer_sedimnt(pdtphys, t_seri, pplay, paprs, tr_seri, dens_aer)
274
275!--compute mass concentration of PM2.5 sulfate particles (wet diameter and mass) at the surface for health studies
276    surf_PM25_sulf(:)=0.0
277    DO i=1,klon
278      DO it=1, nbtr_bin
279        IF (mdw(it) .LT. 2.5e-6) THEN
280          !surf_PM25_sulf(i)=surf_PM25_sulf(i)+tr_seri(i,1,it+nbtr_sulgas)*m_part(i,1,it) &
281          !assume that particles consist of ammonium sulfate at the surface (132g/mol) and are dry at T = 20 deg. C and 50 perc. humidity
282          surf_PM25_sulf(i)=surf_PM25_sulf(i)+tr_seri(i,1,it+nbtr_sulgas) &
283                           & *132./98.*dens_aer_dry*4./3.*RPI*(mdw(it)/2.)**3 &
284                           & *pplay(i,1)/t_seri(i,1)/RD*1e9
285        ENDIF
286      ENDDO
287    ENDDO
288
289!    CALL minmaxsimple(tr_seri(:,:,id_SO2_strat),0.0,0.0,'fin SO2')
290!    DO it=1, nbtr_bin
291!      CALL minmaxsimple(tr_seri(:,:,nbtr_sulgas+it),0.0,0.0,'fin bin ')
292!    ENDDO
293
294  END SUBROUTINE traccoag
295
296END MODULE traccoag_mod
Note: See TracBrowser for help on using the repository browser.