1 | subroutine phytrac_chimie ( |
---|
2 | $ debutphy, |
---|
3 | $ gmtime, |
---|
4 | $ nqmax, |
---|
5 | $ nlon, |
---|
6 | $ lat, |
---|
7 | $ lon, |
---|
8 | $ nlev, |
---|
9 | $ pdtphys, |
---|
10 | $ temp, |
---|
11 | $ pplay, |
---|
12 | $ trac, |
---|
13 | $ d_tr_chem, |
---|
14 | $ iter) |
---|
15 | |
---|
16 | use chemparam_mod |
---|
17 | use conc, only: mmean |
---|
18 | |
---|
19 | implicit none |
---|
20 | |
---|
21 | #include "clesphys.h" |
---|
22 | #include "YOMCST.h" |
---|
23 | |
---|
24 | !=================================================================== |
---|
25 | ! input |
---|
26 | !=================================================================== |
---|
27 | |
---|
28 | integer :: nlon, nlev ! number of gridpoints and levels |
---|
29 | integer :: nqmax ! number of tracers |
---|
30 | |
---|
31 | real :: gmtime ! day fraction |
---|
32 | real :: pdtphys ! phytrac_chimie timestep (s) |
---|
33 | real, dimension(nlon,nlev) :: temp ! temperature (k) |
---|
34 | real, dimension(nlon,nlev) :: pplay ! pressure (pa) |
---|
35 | real, dimension(nlon,nlev,nqmax) :: trac ! tracer mass mixing ratio |
---|
36 | |
---|
37 | logical :: debutphy ! first call flag |
---|
38 | |
---|
39 | !=================================================================== |
---|
40 | ! output |
---|
41 | !=================================================================== |
---|
42 | |
---|
43 | real, dimension(nlon,nlev,nqmax) :: d_tr_chem ! chemical tendency for each tracer |
---|
44 | integer, dimension(nlon,nlev) :: iter ! chemical iterations |
---|
45 | |
---|
46 | !=================================================================== |
---|
47 | ! local |
---|
48 | !=================================================================== |
---|
49 | |
---|
50 | real :: sza_local ! solar zenith angle (deg) |
---|
51 | real :: lon_sun |
---|
52 | |
---|
53 | integer :: i, iq |
---|
54 | integer :: ilon, ilev |
---|
55 | |
---|
56 | real lat(nlon), lat_local(nlon) |
---|
57 | real lon(nlon), lon_local(nlon) |
---|
58 | |
---|
59 | real, dimension(nlon,nlev) :: mrtwv, mrtsa ! total water and total sulfuric acid |
---|
60 | real, dimension(nlon,nlev) :: mrwv, mrsa ! gas-phase water and gas-phase sulfuric acid |
---|
61 | real, dimension(nlon,nlev) :: trac_sum |
---|
62 | real, dimension(nlon,nlev,nqmax) :: ztrac ! local tracer mixing ratio |
---|
63 | |
---|
64 | !=================================================================== |
---|
65 | ! first call : initialisations |
---|
66 | !=================================================================== |
---|
67 | |
---|
68 | if (debutphy) then |
---|
69 | |
---|
70 | !------------------------------------------------------------------- |
---|
71 | ! case of tracers re-initialisation with chemistry |
---|
72 | !------------------------------------------------------------------- |
---|
73 | if (reinit_trac .and. ok_chem) then |
---|
74 | |
---|
75 | print*, "Tracers are re-initialised" |
---|
76 | trac(:,:,:) = 1.0e-30 |
---|
77 | |
---|
78 | if ((i_ocs /= 0) .and. (i_co /= 0) .and. (i_hcl /= 0) |
---|
79 | $ .and. (i_so2 /= 0) .and. (i_h2o /= 0) .and. (i_n2/ = 0) |
---|
80 | $ .and. (i_co2 /= 0)) then |
---|
81 | |
---|
82 | trac(:,1:22,i_ocs) = 3.e-6 |
---|
83 | trac(:,1:22,i_co) = 25.e-6 |
---|
84 | trac(:,:,i_hcl) = 0.4e-6 |
---|
85 | trac(:,1:22,i_so2) = 10.e-6 |
---|
86 | trac(:,1:22,i_h2o) = 30.e-6 |
---|
87 | trac(:,:,i_n2) = 0.35e-1 |
---|
88 | |
---|
89 | ! ensure that sum of mixing ratios = 1 |
---|
90 | |
---|
91 | trac_sum(:,:) = 0. |
---|
92 | |
---|
93 | do iq = 1,nqmax - nmicro |
---|
94 | if (iq /= i_co2) then |
---|
95 | trac_sum(:,:) = trac_sum(:,:) + trac(:,:,iq) |
---|
96 | end if |
---|
97 | end do |
---|
98 | |
---|
99 | ! initialise co2 |
---|
100 | |
---|
101 | trac(:,:,i_co2) = 1. - trac_sum(:,:) |
---|
102 | |
---|
103 | else |
---|
104 | write(*,*) "at least one tracer is missing : stop" |
---|
105 | stop |
---|
106 | end if |
---|
107 | |
---|
108 | ! convert volume to mass mixing ratio |
---|
109 | |
---|
110 | do iq = 1,nqmax - nmicro |
---|
111 | trac(:,:,iq) = trac(:,:,iq)*m_tr(iq)/mmean(:,:) |
---|
112 | end do |
---|
113 | |
---|
114 | end if ! reinit_trac |
---|
115 | |
---|
116 | !------------------------------------------------------------------- |
---|
117 | ! case of detailed microphysics without chemistry |
---|
118 | !------------------------------------------------------------------- |
---|
119 | if (.not. ok_chem .and. ok_cloud .and. cl_scheme == 2) then |
---|
120 | |
---|
121 | ! convert mass to volume mixing ratio |
---|
122 | |
---|
123 | do iq = 1,nqmax - nmicro |
---|
124 | ztrac(:,:,iq) = trac(:,:,iq)*mmean(:,:)/m_tr(iq) |
---|
125 | end do |
---|
126 | |
---|
127 | ! initialise microphysics |
---|
128 | |
---|
129 | call vapors4muphy_ini(nlon,nlev,ztrac) |
---|
130 | |
---|
131 | ! convert volume to mass mixing ratio |
---|
132 | |
---|
133 | do iq = 1,nqmax - nmicro |
---|
134 | trac(:,:,iq) = ztrac(:,:,iq)*m_tr(iq)/mmean(:,:) |
---|
135 | end do |
---|
136 | |
---|
137 | end if |
---|
138 | |
---|
139 | end if ! debutphy |
---|
140 | |
---|
141 | !=================================================================== |
---|
142 | ! convert mass to volume mixing ratio : gas phase |
---|
143 | !=================================================================== |
---|
144 | |
---|
145 | do iq = 1,nqmax - nmicro |
---|
146 | ztrac(:,:,iq) = max(trac(:,:,iq)*mmean(:,:)/m_tr(iq), 1.e-30) |
---|
147 | end do |
---|
148 | |
---|
149 | !=================================================================== |
---|
150 | ! microphysics: simplified scheme (phd aurelien stolzenbach) |
---|
151 | !=================================================================== |
---|
152 | |
---|
153 | if (ok_cloud .and. cl_scheme == 1) then |
---|
154 | |
---|
155 | ! convert mass to volume mixing ratio : liquid phase |
---|
156 | |
---|
157 | ztrac(:,:,i_h2so4liq) = max(trac(:,:,i_h2so4liq) |
---|
158 | $ *mmean(:,:)/m_tr(i_h2so4liq), 1.e-30) |
---|
159 | ztrac(:,:,i_h2oliq) = max(trac(:,:,i_h2oliq) |
---|
160 | $ *mmean(:,:)/m_tr(i_h2oliq), 1.e-30) |
---|
161 | |
---|
162 | ! total water and sulfuric acid (gas + liquid) |
---|
163 | |
---|
164 | mrtwv(:,:) = ztrac(:,:,i_h2o) + ztrac(:,:,i_h2oliq) |
---|
165 | mrtsa(:,:) = ztrac(:,:,i_h2so4) + ztrac(:,:,i_h2so4liq) |
---|
166 | |
---|
167 | ! all water and sulfuric acid is put in the gas-phase |
---|
168 | |
---|
169 | mrwv(:,:) = mrtwv(:,:) |
---|
170 | mrsa(:,:) = mrtsa(:,:) |
---|
171 | |
---|
172 | ! call microphysics |
---|
173 | |
---|
174 | call new_cloud_venus(nlev, nlon, temp, pplay, |
---|
175 | $ mrtwv, mrtsa, mrwv, mrsa) |
---|
176 | |
---|
177 | ! update water vapour and sulfuric acid |
---|
178 | |
---|
179 | ztrac(:,:,i_h2o) = mrwv(:,:) |
---|
180 | ztrac(:,:,i_h2oliq) = mrtwv(:,:) - ztrac(:,:,i_h2o) |
---|
181 | |
---|
182 | ztrac(:,:,i_h2so4) = mrsa(:,:) |
---|
183 | ztrac(:,:,i_h2so4liq) = mrtsa(:,:) - ztrac(:,:,i_h2so4) |
---|
184 | |
---|
185 | end if ! simplified scheme |
---|
186 | |
---|
187 | !=================================================================== |
---|
188 | ! microphysics: detailed scheme (phd sabrina guilbon) |
---|
189 | !=================================================================== |
---|
190 | |
---|
191 | if (ok_cloud .and. cl_scheme == 2) then |
---|
192 | |
---|
193 | do ilon = 1,nlon |
---|
194 | do ilev = 1, nlev |
---|
195 | if (temp(ilon,ilev) < 500.) then |
---|
196 | call mad_muphy(pdtphys, ! timestep |
---|
197 | $ temp(ilon,ilev),pplay(ilon,ilev), ! temperature and pressure |
---|
198 | $ ztrac(ilon,ilev,i_h2o), |
---|
199 | $ ztrac(ilon,ilev,i_h2so4), |
---|
200 | $ ztrac(ilon,ilev,i_m0_aer), |
---|
201 | $ ztrac(ilon,ilev,i_m3_aer), |
---|
202 | $ ztrac(ilon,ilev,i_m0_mode1drop), |
---|
203 | $ ztrac(ilon,ilev,i_m0_mode1ccn), |
---|
204 | $ ztrac(ilon,ilev,i_m3_mode1sa), |
---|
205 | $ ztrac(ilon,ilev,i_m3_mode1w), |
---|
206 | $ ztrac(ilon,ilev,i_m3_mode1ccn), |
---|
207 | $ ztrac(ilon,ilev,i_m0_mode2drop), |
---|
208 | $ ztrac(ilon,ilev,i_m0_mode2ccn), |
---|
209 | $ ztrac(ilon,ilev,i_m3_mode2sa), |
---|
210 | $ ztrac(ilon,ilev,i_m3_mode2w), |
---|
211 | $ ztrac(ilon,ilev,i_m3_mode2ccn)) |
---|
212 | else |
---|
213 | ztrac(ilon,ilev,i_m0_aer) = 0. |
---|
214 | ztrac(ilon,ilev,i_m3_aer) = 0. |
---|
215 | ztrac(ilon,ilev,i_m0_mode1drop) = 0. |
---|
216 | ztrac(ilon,ilev,i_m0_mode1ccn) = 0. |
---|
217 | ztrac(ilon,ilev,i_m3_mode1sa) = 0. |
---|
218 | ztrac(ilon,ilev,i_m3_mode1w) = 0. |
---|
219 | ztrac(ilon,ilev,i_m3_mode1ccn) = 0. |
---|
220 | ztrac(ilon,ilev,i_m0_mode2drop) = 0. |
---|
221 | ztrac(ilon,ilev,i_m0_mode2ccn) = 0. |
---|
222 | ztrac(ilon,ilev,i_m3_mode2sa) = 0. |
---|
223 | ztrac(ilon,ilev,i_m3_mode2w) = 0. |
---|
224 | ztrac(ilon,ilev,i_m3_mode2ccn) = 0. |
---|
225 | end if |
---|
226 | end do |
---|
227 | end do |
---|
228 | |
---|
229 | end if ! detailed scheme |
---|
230 | |
---|
231 | !=================================================================== |
---|
232 | ! photochemistry |
---|
233 | !=================================================================== |
---|
234 | |
---|
235 | if (ok_chem) then |
---|
236 | |
---|
237 | lon_sun = (0.5 - gmtime)*2.*rpi |
---|
238 | lon_local(:) = lon(:)*rpi/180. |
---|
239 | lat_local(:) = lat(:)*rpi/180. |
---|
240 | |
---|
241 | do ilon = 1,nlon |
---|
242 | |
---|
243 | ! solar zenith angle |
---|
244 | |
---|
245 | sza_local = acos(cos(lat_local(ilon))*cos(lon_local(ilon)) |
---|
246 | $ *cos(lon_sun) + cos(lat_local(ilon)) |
---|
247 | $ *sin(lon_local(ilon))*sin(lon_sun))*180./rpi |
---|
248 | |
---|
249 | call photochemistry_venus(nlev, nlon, pdtphys, |
---|
250 | $ pplay(ilon,:)/100., |
---|
251 | $ temp(ilon,:), |
---|
252 | $ ztrac(ilon,:,:), |
---|
253 | $ mmean(ilon,:), |
---|
254 | $ sza_local, nqmax, iter(ilon,:)) |
---|
255 | |
---|
256 | end do |
---|
257 | |
---|
258 | end if ! ok_chem |
---|
259 | |
---|
260 | !=================================================================== |
---|
261 | ! compute tendencies |
---|
262 | !=================================================================== |
---|
263 | |
---|
264 | ! gas phase |
---|
265 | |
---|
266 | do iq = 1,nqmax - nmicro |
---|
267 | ztrac(:,:,iq) = ztrac(:,:,iq)*m_tr(iq)/mmean(:,:) |
---|
268 | d_tr_chem(:,:,iq) = (ztrac(:,:,iq) - trac(:,:,iq))/pdtphys |
---|
269 | end do |
---|
270 | |
---|
271 | ! liquid phase |
---|
272 | |
---|
273 | if (ok_cloud .and. cl_scheme == 1) then |
---|
274 | ztrac(:,:,i_h2so4liq) = ztrac(:,:,i_h2so4liq)*m_tr(i_h2so4liq) |
---|
275 | $ /mmean(:,:) |
---|
276 | ztrac(:,:,i_h2oliq) = ztrac(:,:,i_h2oliq)*m_tr(i_h2oliq) |
---|
277 | $ /mmean(:,:) |
---|
278 | d_tr_chem(:,:,i_h2so4liq) = (ztrac(:,:,i_h2so4liq) |
---|
279 | $ - trac(:,:,i_h2so4liq))/pdtphys |
---|
280 | d_tr_chem(:,:,i_h2oliq) = (ztrac(:,:,i_h2oliq) |
---|
281 | $ - trac(:,:,i_h2oliq))/pdtphys |
---|
282 | end if |
---|
283 | |
---|
284 | end subroutine phytrac_chimie |
---|