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