MODULE traccoag_mod ! ! This module calculates the concentration of aerosol particles in certain size bins ! considering coagulation and sedimentation. ! CONTAINS SUBROUTINE traccoag(pdtphys, gmtime, debutphy, julien, & presnivs, xlat, xlon, pphis, pphi, & t_seri, pplay, paprs, sh, rh , & tr_seri) USE phys_local_var_mod, ONLY: mdw, sulf_convert, sulf_nucl, sulf_cond_evap, & & sfluxaer, ocs_convert, R2SO4, DENSO4, f_r_wet, SO2_backgr_tend, OCS_backgr_tend, & & OCS_lifetime, SO2_lifetime, surf_PM25_sulf USE dimphy USE infotrac USE aerophys USE geometry_mod, ONLY : cell_area USE mod_grid_phy_lmdz USE mod_phys_lmdz_mpi_data, ONLY : is_mpi_root USE mod_phys_lmdz_para, only: gather, scatter USE phys_cal_mod USE sulfate_aer_mod USE phys_local_var_mod, ONLY: stratomask USE sulfate_aer_mod, ONLY : stracomp, denh2sa USE YOMCST IMPLICIT NONE ! Input argument !--------------- REAL,INTENT(IN) :: pdtphys ! Pas d'integration pour la physique (seconde) REAL,INTENT(IN) :: gmtime ! Heure courante LOGICAL,INTENT(IN) :: debutphy ! le flag de l'initialisation de la physique INTEGER,INTENT(IN) :: julien ! Jour julien REAL,DIMENSION(klev),INTENT(IN) :: presnivs! pressions approximat. des milieux couches (en PA) REAL,DIMENSION(klon),INTENT(IN) :: xlat ! latitudes pour chaque point REAL,DIMENSION(klon),INTENT(IN) :: xlon ! longitudes pour chaque point REAL,DIMENSION(klon),INTENT(IN) :: pphis ! geopotentiel du sol REAL,DIMENSION(klon,klev),INTENT(IN) :: pphi ! geopotentiel de chaque couche REAL,DIMENSION(klon,klev),INTENT(IN) :: t_seri ! Temperature REAL,DIMENSION(klon,klev),INTENT(IN) :: pplay ! pression pour le mileu de chaque couche (en Pa) REAL,DIMENSION(klon,klev+1),INTENT(IN) :: paprs ! pression pour chaque inter-couche (en Pa) REAL,DIMENSION(klon,klev),INTENT(IN) :: sh ! humidite specifique REAL,DIMENSION(klon,klev),INTENT(IN) :: rh ! humidite relative ! Output argument !---------------- REAL,DIMENSION(klon,klev,nbtr),INTENT(INOUT) :: tr_seri ! Concentration Traceur [U/KgA] ! Local variables !---------------- ! flag for sulfur emission scenario: (0) background aerosol ; (1) volcanic eruption ; (2) stratospheric aerosol injections (SAI) INTEGER,PARAMETER :: flag_sulf_emit=2 ! !--flag_sulf_emit=1 --example Pinatubo INTEGER,PARAMETER :: year_emit_vol=1991 ! year of emission date INTEGER,PARAMETER :: mth_emit_vol=6 ! month of emission date INTEGER,PARAMETER :: day_emit_vol=15 ! day of emission date REAL,PARAMETER :: m_aer_emiss_vol=7.e9 ! emitted sulfur mass in kgS, e.g. 7Tg(S)=14Tg(SO2) REAL,PARAMETER :: altemiss_vol=17.e3 ! emission altitude in m REAL,PARAMETER :: sigma_alt_vol=1.e3 ! standard deviation of emission altitude in m REAL,PARAMETER :: xlat_vol=15.14 ! latitude of volcano in degree REAL,PARAMETER :: xlon_vol=120.35 ! longitude of volcano in degree !--flag_sulf_emit=2 --SAI REAL,PARAMETER :: m_aer_emiss_sai=1.e10 ! emitted sulfur mass in kgS, eg 1e9=1TgS, 1e10=10TgS REAL,PARAMETER :: altemiss_sai=17.e3 ! emission altitude in m REAL,PARAMETER :: sigma_alt_sai=1.e3 ! standard deviation of emission altitude in m REAL,PARAMETER :: xlat_sai=0.0 ! latitude of SAI in degree REAL,PARAMETER :: xlon_sai=120.35 ! longitude of SAI in degree !--be careful - this needs to be changed with resolution - here for 96x95 - seems to work for 48x36 as well REAL,PARAMETER :: dlat=0.9474 ! d latitude in degree REAL,PARAMETER :: dlon=1.875 ! d longitude in degree !--other local variables INTEGER :: it, k, i, ilon, ilev, itime, i_int LOGICAL,DIMENSION(klon,klev) :: is_strato ! true = above tropopause, false = below REAL,DIMENSION(klon,klev) :: m_air_gridbox ! mass of air in every grid box [kg] REAL,DIMENSION(klon_glo,klev,nbtr) :: tr_seri_glo ! Concentration Traceur [U/KgA] REAL,DIMENSION(klev+1) :: altLMDz ! altitude of layer interfaces in m REAL,DIMENSION(klev) :: f_lay_emiss ! fraction of emission for every vertical layer REAL :: f_lay_sum ! sum of layer emission fractions REAL :: alt ! altitude for integral calculation INTEGER,PARAMETER :: n_int_alt=10 ! number of subintervals for integration over Gaussian emission profile REAL,DIMENSION(nbtr_bin) :: r_bin ! particle radius in size bin [m] REAL,DIMENSION(nbtr_bin) :: r_lower ! particle radius at lower bin boundary [m] REAL,DIMENSION(nbtr_bin) :: r_upper ! particle radius at upper bin boundary [m] REAL,DIMENSION(nbtr_bin) :: m_part_dry ! mass of one dry particle in size bin [kg] REAL :: zrho ! Density of air [kg/m3] REAL :: zdz ! thickness of atm. model layer in m REAL,DIMENSION(klon,klev) :: dens_aer ! density of aerosol particles [kg/m3 aerosol] with default H2SO4 mass fraction IF (is_mpi_root) THEN PRINT *,'in traccoag: date from phys_cal_mod =',year_cur,'-',mth_cur,'-',day_cur,'-',hour ENDIF DO it=1, nbtr_bin r_bin(it)=mdw(it)/2. ENDDO !--set boundaries of size bins DO it=1, nbtr_bin IF (it.EQ.1) THEN r_upper(it)=sqrt(r_bin(it+1)*r_bin(it)) r_lower(it)=r_bin(it)**2./r_upper(it) ELSEIF (it.EQ.nbtr_bin) THEN r_lower(it)=sqrt(r_bin(it)*r_bin(it-1)) r_upper(it)=r_bin(it)**2./r_lower(it) ELSE r_lower(it)=sqrt(r_bin(it)*r_bin(it-1)) r_upper(it)=sqrt(r_bin(it+1)*r_bin(it)) ENDIF ENDDO IF (debutphy .and. is_mpi_root) THEN DO it=1, nbtr_bin PRINT *,'radius bin', it, ':', r_bin(it), '(from', r_lower(it), 'to', r_upper(it), ')' ENDDO ENDIF !--initialising logical is_strato from stratomask is_strato(:,:)=.FALSE. WHERE (stratomask.GT.0.5) is_strato=.TRUE. ! STRACOMP (H2O, P, t_seri -> aerosol composition (R2SO4)) ! H2SO4 mass fraction in aerosol (%) CALL stracomp(sh,t_seri,pplay) ! aerosol density (gr/cm3) CALL denh2sa(t_seri) ! compute factor for converting dry to wet radius (for every grid box) f_r_wet(:,:) = (dens_aer_dry/(DENSO4(:,:)*1000.)/(R2SO4(:,:)/100.))**(1./3.) !--calculate mass of air in every grid box DO ilon=1, klon DO ilev=1, klev m_air_gridbox(ilon,ilev)=(paprs(ilon,ilev)-paprs(ilon,ilev+1)) / RG * cell_area(ilon) ENDDO ENDDO IF (debutphy) THEN CALL gather(tr_seri, tr_seri_glo) IF (MAXVAL(tr_seri_glo).LT.1.e-30) THEN !--initialising tracer concentrations to zero DO it=1, nbtr tr_seri(:,:,it)=0.0 ENDDO ENDIF ENDIF !--sulfur emission, depending on chosen scenario (flag_sulf_emit) SELECT CASE(flag_sulf_emit) CASE(0) ! background aerosol ! do nothing (no emission) CASE(1) ! volcanic eruption !--only emit on day of eruption ! stretch emission over one day of Pinatubo eruption IF (year_cur==year_emit_vol.AND.mth_cur==mth_emit_vol.AND.day_cur==day_emit_vol) THEN ! DO i=1,klon !Pinatubo eruption at 15.14N, 120.35E IF ( xlat(i).GT.xlat_vol-dlat .AND. xlat(i).LT.xlat_vol+dlat .AND. & xlon(i).GT.xlon_vol-dlon .AND. xlon(i).LT.xlon_vol+dlon ) THEN ! compute altLMDz altLMDz(:)=0.0 DO k=1, klev zrho=pplay(i,k)/t_seri(i,k)/RD !air density in kg/m3 zdz=(paprs(i,k)-paprs(i,k+1))/zrho/RG !thickness of layer in m altLMDz(k+1)=altLMDz(k)+zdz ENDDO !compute distribution of emission to vertical model layers (based on Gaussian peak in altitude) f_lay_sum=0.0 DO k=1, klev f_lay_emiss(k)=0.0 DO i_int=1, n_int_alt alt=altLMDz(k)+float(i_int)*(altLMDz(k+1)-altLMDz(k))/float(n_int_alt) f_lay_emiss(k)=f_lay_emiss(k)+1./(sqrt(2.*RPI)*sigma_alt_vol)* & & exp(-0.5*((alt-altemiss_vol)/sigma_alt_vol)**2.)* & & (altLMDz(k+1)-altLMDz(k))/float(n_int_alt) ENDDO f_lay_sum=f_lay_sum+f_lay_emiss(k) ENDDO !correct for step integration error f_lay_emiss(:)=f_lay_emiss(:)/f_lay_sum !emission as SO2 gas (with m(SO2)=64/32*m_aer_emiss) !vertically distributed emission DO k=1, klev tr_seri(i,k,id_SO2_strat)=tr_seri(i,k,id_SO2_strat)+ & & m_aer_emiss_vol*(mSO2mol/mSatom)/m_air_gridbox(i,k)*f_lay_emiss(k) & & /(1.*86400./pdtphys) ! stretch emission over one day of Pinatubo eruption ENDDO ENDIF ! emission grid cell ENDDO ! klon loop ENDIF ! emission period CASE(2) ! stratospheric aerosol injections (SAI) ! DO i=1,klon ! SAI standard scenario with continuous emission from 1 grid point at the equator ! SAI emission on single month ! IF ((mth_cur==4 .AND. & ! SAI continuous emission o IF ( xlat(i).GT.xlat_sai-dlat .AND. xlat(i).LT.xlat_sai+dlat .AND. & & xlon(i).GT.xlon_sai-dlon .AND. xlon(i).LT.xlon_sai+dlon ) THEN ! compute altLMDz altLMDz(:)=0.0 DO k=1, klev zrho=pplay(i,k)/t_seri(i,k)/RD !air density in kg/m3 zdz=(paprs(i,k)-paprs(i,k+1))/zrho/RG !thickness of layer in m altLMDz(k+1)=altLMDz(k)+zdz ENDDO !compute distribution of emission to vertical model layers (based on Gaussian peak in altitude) f_lay_sum=0.0 DO k=1, klev f_lay_emiss(k)=0.0 DO i_int=1, n_int_alt alt=altLMDz(k)+float(i_int)*(altLMDz(k+1)-altLMDz(k))/float(n_int_alt) f_lay_emiss(k)=f_lay_emiss(k)+1./(sqrt(2.*RPI)*sigma_alt_sai)* & & exp(-0.5*((alt-altemiss_sai)/sigma_alt_sai)**2.)* & & (altLMDz(k+1)-altLMDz(k))/float(n_int_alt) ENDDO f_lay_sum=f_lay_sum+f_lay_emiss(k) ENDDO !correct for step integration error f_lay_emiss(:)=f_lay_emiss(:)/f_lay_sum !emission as SO2 gas (with m(SO2)=64/32*m_aer_emiss) !vertically distributed emission DO k=1, klev tr_seri(i,k,id_SO2_strat)=tr_seri(i,k,id_SO2_strat)+ & & m_aer_emiss_sai*(mSO2mol/mSatom)/m_air_gridbox(i,k)*f_lay_emiss(k) & & /(360.*86400./pdtphys) ! stretch emission over whole year (360d) ! & /(60.*86400./pdtphys) ! stretch emission over 2 months (seasonal emission) ! & /7. ! distribute equally over 7 emission grid points ENDDO ! !emission as monodisperse particles with 0.1um dry radius (BIN21) ! !vertically distributed emission ! DO k=1, klev ! tr_seri(i,k,id_BIN01_strat+20)=tr_seri(i,k,id_BIN01_strat+20)+ & ! & m_aer_emiss*(mH2SO4mol/mSatom)/m_part_dry(21)/m_air_gridbox(i,k)*f_lay_emiss(k) & ! & /(360.*86400./pdtphys) & ! stretch emission over whole year (360d) ! & /7. ! distribute equally over 7 emission grid points ! ENDDO ENDIF ! emission grid cell ENDDO ! klon loop END SELECT ! emission scenario (flag_sulf_emit) !--read background concentrations of OCS and SO2 and lifetimes from input file !--update the variables defined in phys_local_var_mod CALL interp_sulf_input(debutphy,pdtphys,paprs,tr_seri) !--convert OCS to SO2 in the stratosphere CALL ocs_to_so2(pdtphys,tr_seri,t_seri,pplay,paprs,sh,is_strato) !--convert SO2 to H2SO4 CALL so2_to_h2so4(pdtphys,tr_seri,t_seri,pplay,paprs,sh,is_strato) !--common routine for nucleation and condensation/evaporation with adaptive timestep CALL micphy_tstep(pdtphys,tr_seri,t_seri,pplay,paprs,rh,is_strato) !--call coagulation routine CALL coagulate(pdtphys,mdw,tr_seri,t_seri,pplay,dens_aer,is_strato) !--call sedimentation routine CALL aer_sedimnt(pdtphys, t_seri, pplay, paprs, tr_seri, dens_aer) !--compute mass concentration of PM2.5 sulfate particles (wet diameter and mass) at the surface for health studies surf_PM25_sulf(:)=0.0 DO i=1,klon DO it=1, nbtr_bin IF (mdw(it) .LT. 2.5e-6) THEN !surf_PM25_sulf(i)=surf_PM25_sulf(i)+tr_seri(i,1,it+nbtr_sulgas)*m_part(i,1,it) & !assume that particles consist of ammonium sulfate at the surface (132g/mol) and are dry at T = 20 deg. C and 50 perc. humidity surf_PM25_sulf(i)=surf_PM25_sulf(i)+tr_seri(i,1,it+nbtr_sulgas) & & *132./98.*dens_aer_dry*4./3.*RPI*(mdw(it)/2.)**3 & & *pplay(i,1)/t_seri(i,1)/RD*1e9 ENDIF ENDDO ENDDO ! CALL minmaxsimple(tr_seri(:,:,id_SO2_strat),0.0,0.0,'fin SO2') ! DO it=1, nbtr_bin ! CALL minmaxsimple(tr_seri(:,:,nbtr_sulgas+it),0.0,0.0,'fin bin ') ! ENDDO END SUBROUTINE traccoag END MODULE traccoag_mod