[3247] | 1 | subroutine prodhaze(ngrid,nlayer,nq,ptimestep, & |
---|
| 2 | pplev,pq,pdq,pdist_sol,mu0,declin,zdqfasthaze,zdqsfasthaze,gradflux, & |
---|
| 3 | flym_bot,flym_sol_bot,flym_ipm_bot,flym_sol,flym_ipm) |
---|
| 4 | |
---|
| 5 | use radinc_h, only : naerkind |
---|
| 6 | use comgeomfi_h |
---|
| 7 | use comcstfi_mod, only: pi, g |
---|
| 8 | use tracer_h, only: igcm_ch4_gas, igcm_prec_haze, noms, mmol |
---|
| 9 | use geometry_mod, only: longitude, latitude ! in radians |
---|
| 10 | use callkeys_mod, only: hazeconservch4, diurnal |
---|
| 11 | implicit none |
---|
| 12 | |
---|
| 13 | !================================================================== |
---|
| 14 | ! Purpose |
---|
| 15 | ! ------- |
---|
| 16 | ! Fast production of haze in the atmosphere, direct accumulation to surface |
---|
| 17 | ! |
---|
| 18 | ! Authors |
---|
| 19 | ! ------- |
---|
| 20 | ! Tanguy Bertrand and Francois Forget (2014) |
---|
| 21 | ! |
---|
| 22 | !================================================================== |
---|
| 23 | |
---|
| 24 | !----------------------------------------------------------------------- |
---|
| 25 | ! Arguments |
---|
| 26 | |
---|
| 27 | INTEGER ngrid, nlayer,nq |
---|
| 28 | LOGICAL firstcall |
---|
| 29 | SAVE firstcall |
---|
| 30 | DATA firstcall/.true./ |
---|
| 31 | REAL ptimestep |
---|
| 32 | REAL pplev(ngrid,nlayer+1) |
---|
| 33 | REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) |
---|
| 34 | REAL,INTENT(IN) :: pdq(ngrid,nlayer,nq) |
---|
| 35 | REAL,INTENT(IN) :: pdist_sol ! distance SUN-pluto in AU |
---|
| 36 | REAL,INTENT(IN) :: mu0(ngrid) ! cosinus of solar incident flux |
---|
| 37 | REAL,INTENT(IN) :: declin ! declination |
---|
| 38 | REAL,INTENT(OUT) :: zdqfasthaze(ngrid,nq) ! Final tendancy ch4 |
---|
| 39 | REAL,INTENT(OUT) :: zdqsfasthaze(ngrid) ! Final tendancy haze surf |
---|
| 40 | REAL, INTENT(OUT) :: gradflux(ngrid) ! gradient flux Lyman alpha ph.m-2.s-1 |
---|
| 41 | REAL, INTENT(OUT) :: flym_bot(ngrid) ! flux Lyman alpha bottom ph.m-2.s-1 |
---|
| 42 | REAL, INTENT(OUT) :: flym_sol_bot(ngrid) ! Solar flux Lyman alpha ph.m-2.s-1 reaching the surface |
---|
| 43 | REAL, INTENT(OUT) :: flym_ipm_bot(ngrid) ! IPM (Interplanetary) flux Lyman alpha ph.m-2.s-1 reaching the surface |
---|
| 44 | REAL, INTENT(OUT) :: flym_sol(ngrid) ! Incident Solar flux Lyman alpha ph.m-2.s-1 |
---|
| 45 | REAL, INTENT(OUT) :: flym_ipm(ngrid) ! Incident IPM (Interplanetary) flux Lyman alpha ph.m-2.s-1 |
---|
| 46 | !----------------------------------------------------------------------- |
---|
| 47 | ! Local variables |
---|
| 48 | |
---|
| 49 | INTEGER l,ig,iq |
---|
| 50 | REAL zq_ch4(ngrid) |
---|
| 51 | REAL zq_prec(ngrid) |
---|
| 52 | REAL tauch4 |
---|
| 53 | REAL sigch4 ! aborption cross section ch4 cm-2 mol-1 |
---|
| 54 | REAL flym_sol_earth ! initial flux Earth ph.m-2.s-1 |
---|
| 55 | REAL flym_sol_pluto ! initial Lyman alpha SOLAR flux Pluto ph.m-2.s-1 |
---|
| 56 | REAL mu_ipm(ngrid) ! local Mean incident flux for IPM Lyman alpha photons |
---|
| 57 | |
---|
| 58 | REAL kch4 ! fraction of Carbon from ch4 directly dissociated by prec_haze |
---|
| 59 | REAL ncratio ! ration Nitrogen/Carbon observed in tholins |
---|
| 60 | REAL avogadro ! avogadro constant |
---|
| 61 | |
---|
| 62 | CHARACTER(len=10) :: tname |
---|
| 63 | REAL pagg,pform |
---|
| 64 | REAL tconv |
---|
| 65 | REAL longit |
---|
| 66 | REAL latit |
---|
| 67 | REAL valmin |
---|
| 68 | REAL valmax |
---|
| 69 | REAL valmin_dl |
---|
| 70 | REAL puis |
---|
| 71 | REAL dist |
---|
| 72 | REAL zdqprec(ngrid) |
---|
| 73 | REAL zdqch4(ngrid) |
---|
| 74 | REAL zdqconv_prec(ngrid) |
---|
| 75 | !----------------------------------------------------------------------- |
---|
| 76 | |
---|
| 77 | !******** INPUT |
---|
| 78 | avogadro = 6.022141e23 |
---|
| 79 | sigch4 = 1.85e-17 ! aborption cross section ch4 cm-2 mol-1 |
---|
| 80 | !! Initial Solar flux Lyman alpha on Earth (1AU) |
---|
| 81 | flym_sol_earth=5.e15 ! ph.m-2.s-1 --> 5e+11 ph.cm-2.s-1 |
---|
| 82 | !! Parameter of conversion precurseur to haze : here accelaeration/deceleration of haze producition |
---|
| 83 | kch4=1. |
---|
| 84 | ncratio=0.5 ! boost for haze considering nitrogen contribution |
---|
| 85 | ! ratio n/c : =0.25 if there is 1N for 3C |
---|
| 86 | pagg=0.05 ! Minimum pressure required for aggregation |
---|
| 87 | pform=0.006 ! Minimum pressure required for haze formation |
---|
| 88 | tconv=1.e7 ! (s) time needed to convert gas precursors into aerosol |
---|
| 89 | !******** Incident Solar flux |
---|
| 90 | !! Initial Solar flux Lyman alpha on Pluto |
---|
| 91 | flym_sol_pluto=0.25*flym_sol_earth/(pdist_sol*pdist_sol) ! ph.m-2.s-1 |
---|
| 92 | ! Decreasing by 12% the initial IPM flux to account for Interplanetary H absorption: |
---|
| 93 | ! Fig 3 Gladstone et al 2014 : Lyalpha at Pluto |
---|
| 94 | flym_sol_pluto=flym_sol_pluto*(1.-0.12*pdist_sol/32.91) |
---|
| 95 | |
---|
| 96 | !******** Incident IPM (Interplanetary) flux Lyman alpha ph.m-2.s-1 |
---|
| 97 | ! fit Fig. 4 in Randall Gladstone et al. (2015) |
---|
| 98 | ! -- New version : Integration over semi sphere of Gladstone data. |
---|
| 99 | ! Fit of results |
---|
| 100 | |
---|
| 101 | if (diurnal) then |
---|
| 102 | ! 1) get longitude/latitude (rad) of anti-subsolar point (max de mu0 - 180) |
---|
| 103 | longit=longitude((MAXLOC(mu0,DIM=1,MASK=mu0.GT.0.9))) ! rad |
---|
| 104 | IF (longit.GE.0) THEN |
---|
| 105 | longit=longit-pi |
---|
| 106 | ELSE |
---|
| 107 | longit=longit+pi |
---|
| 108 | ENDIF |
---|
| 109 | latit=-declin ! rad |
---|
| 110 | |
---|
| 111 | ! 2) Define input parameter for the fit |
---|
| 112 | valmin=48.74e10 ! minimum value of ipm flux in ph/m2/s |
---|
| 113 | valmax=115.014e10 ! maximal value of ipm flux in ph/m2/s |
---|
| 114 | valmin_dl=74.5e10 ! daylight minimum value |
---|
| 115 | puis=3.5 |
---|
| 116 | |
---|
| 117 | ! 3) Loop for each location and fit |
---|
| 118 | DO ig=1,ngrid |
---|
| 119 | !!! Solar flux |
---|
| 120 | flym_sol(ig)=flym_sol_pluto*mu0(ig) |
---|
| 121 | |
---|
| 122 | !!! IPM flux |
---|
| 123 | ! calculation of cosinus of incident angle for IPM flux |
---|
| 124 | mu_ipm(ig) = max(mu0(ig), 0.5) |
---|
| 125 | IF (mu0(ig).LT.1.e-4) THEN |
---|
| 126 | dist=acos(sin(latitude(ig))*sin(latit)+cos(latitude(ig))* & |
---|
| 127 | cos(latit)*cos(longit-longitude(ig)))*180/pi |
---|
| 128 | IF (dist.gt.90) dist=90 |
---|
| 129 | flym_ipm(ig)=(valmin_dl-valmin)/(90.**puis)*(dist)**puis & |
---|
| 130 | +valmin |
---|
| 131 | ELSE |
---|
| 132 | flym_ipm(ig)= mu0(ig)*(valmax-valmin_dl)+valmin_dl |
---|
| 133 | |
---|
| 134 | !! IPM flux must be weighted by solar flux |
---|
| 135 | flym_ipm(ig)=flym_ipm(ig)*flym_sol_pluto/1.0156e12 ! Value in 2015 for flym_sol_pluto |
---|
| 136 | ENDIF |
---|
| 137 | ENDDO ! ngrid |
---|
| 138 | |
---|
| 139 | ELSE ! diurnal=false |
---|
| 140 | |
---|
| 141 | DO ig=1,ngrid |
---|
| 142 | flym_sol(ig)=flym_sol_pluto*mu0(ig) |
---|
| 143 | flym_ipm(ig)=flym_sol(ig)*1.15 |
---|
| 144 | mu_ipm(ig) = max(mu0(ig), 0.5) |
---|
| 145 | ENDDO |
---|
| 146 | |
---|
| 147 | ENDIF ! diurnal |
---|
| 148 | |
---|
| 149 | IF (firstcall) then |
---|
| 150 | write(*,*) 'Fast haze parameters:' |
---|
| 151 | write(*,*) 'kch4, ncratio = ',kch4,ncratio |
---|
| 152 | firstcall=.false. |
---|
| 153 | ENDIF |
---|
| 154 | |
---|
| 155 | ! note: mu0 = cos(solar zenith angle) |
---|
| 156 | ! max(mu0) = declin |
---|
| 157 | |
---|
| 158 | !************************** |
---|
| 159 | !******** LOOP ************ |
---|
| 160 | !************************** |
---|
| 161 | zdqfasthaze(:,:)=0. |
---|
| 162 | zdqsfasthaze(:)=0. |
---|
| 163 | zdqch4(:)=0. |
---|
| 164 | zdqprec(:)=0. |
---|
| 165 | DO ig=1,ngrid |
---|
| 166 | |
---|
| 167 | !! Update of tracer ch4 and prec_haze |
---|
| 168 | zq_ch4(ig)=pq(ig,1,igcm_ch4_gas)+pdq(ig,1,igcm_ch4_gas)*ptimestep |
---|
| 169 | zq_prec(ig)=pq(ig,1,igcm_prec_haze)+pdq(ig,1,igcm_prec_haze)*ptimestep |
---|
| 170 | !! Security as we loose CH4 molecules |
---|
| 171 | if (zq_ch4(ig).lt.0.) then |
---|
| 172 | zq_ch4(ig)=0. |
---|
| 173 | endif |
---|
| 174 | |
---|
| 175 | !! Calculation of CH4 optical depth in Lyman alpha for whole atm |
---|
| 176 | tauch4=sigch4*1.e-4*avogadro*(zq_ch4(ig)/(mmol(igcm_ch4_gas)*1.e-3))*(pplev(ig,1))/g |
---|
| 177 | |
---|
| 178 | !! Calculation of Flux reaching the surface |
---|
| 179 | if (mu0(ig).gt.1.e-6) then |
---|
| 180 | flym_sol_bot(ig)=flym_sol(ig)*exp(-tauch4/(mu0(ig))) |
---|
| 181 | endif |
---|
| 182 | flym_ipm_bot(ig)=flym_ipm(ig)*exp(-tauch4/mu_ipm(ig)) |
---|
| 183 | |
---|
| 184 | !! Total flux reaching the surface and flux absorbed by CH4 |
---|
| 185 | gradflux(ig)=flym_sol(ig)-flym_sol_bot(ig) + flym_ipm(ig)-flym_ipm_bot(ig) |
---|
| 186 | flym_bot(ig)=flym_ipm_bot(ig)+flym_sol_bot(ig) |
---|
| 187 | |
---|
| 188 | !! Tendancy on ch4 considering 1 photon destroys 1 ch4 by photolysis |
---|
| 189 | zdqch4(ig)=-mmol(igcm_ch4_gas)*1.e-3/avogadro*gradflux(ig)*g/(pplev(ig,1)) |
---|
| 190 | !! Tendency of precursors created |
---|
| 191 | zdqprec(ig)=-zdqch4(ig) |
---|
| 192 | |
---|
| 193 | !! update precurseur zq_prec |
---|
| 194 | zq_prec(ig)=zq_prec(ig)+zdqprec(ig)*ptimestep |
---|
| 195 | |
---|
| 196 | !! If pressure > pagg, aggregation is possible |
---|
| 197 | !! If pressure > pform, haze formation is possible |
---|
| 198 | if (pplev(ig,1).gt.pform) then |
---|
| 199 | ! New tendancy for prec_haze, turning into haze aerosols |
---|
| 200 | zdqconv_prec(ig) = -zq_prec(ig)*(1.-exp(-ptimestep/tconv))/ptimestep |
---|
| 201 | else |
---|
| 202 | zdqconv_prec(ig) = 0. |
---|
| 203 | endif |
---|
| 204 | |
---|
| 205 | !! Final tendencies |
---|
| 206 | DO iq=1,nq |
---|
| 207 | tname=noms(iq) |
---|
| 208 | if (tname(1:4).eq."haze") then |
---|
| 209 | ! Tendency haze in atmosphere kg/kg |
---|
| 210 | zdqfasthaze(ig,iq) = -zdqconv_prec(ig)*kch4*(1.+ncratio*14./12.) |
---|
| 211 | ! Tendancy for haze directly deposited at the surface kg/m2 |
---|
| 212 | zdqsfasthaze(ig) = zdqfasthaze(ig,iq)*pplev(ig,1)/g |
---|
| 213 | else if (noms(iq).eq."prec_haze") then |
---|
| 214 | zdqfasthaze(ig,igcm_prec_haze)= zdqprec(ig) + zdqconv_prec(ig) |
---|
| 215 | else if (noms(iq).eq."ch4_gas".and.(.not.hazeconservch4)) then |
---|
| 216 | ! Tendancy for haze directly deposited at the surface kg/m2 |
---|
| 217 | zdqfasthaze(ig,igcm_ch4_gas)= zdqch4(ig) |
---|
| 218 | endif |
---|
| 219 | ENDDO |
---|
| 220 | |
---|
| 221 | |
---|
| 222 | |
---|
| 223 | ENDDO ! ngrid |
---|
| 224 | |
---|
| 225 | end subroutine prodhaze |
---|
| 226 | |
---|