Line data Source code
1 : !-----------------------------------------------------------------------
2 : !
3 : ! Manages the adjustment of ClOy and BrOy family components in response
4 : ! to conservation issues resulting from advection.
5 : !
6 : ! Created by: Francis Vitt
7 : ! Date: 21 May 2008
8 : ! Modified by Stacy Walters
9 : ! Date: 13 August 2008
10 : !-----------------------------------------------------------------------
11 :
12 : module clybry_fam
13 :
14 : use shr_kind_mod, only : r8 => shr_kind_r8
15 : use ppgrid, only : pcols, pver
16 : use chem_mods, only : gas_pcnst, adv_mass
17 : use constituents, only : pcnst
18 : use short_lived_species,only: set_short_lived_species,get_short_lived_species
19 :
20 : implicit none
21 :
22 : save
23 :
24 : private
25 : public :: clybry_fam_set
26 : public :: clybry_fam_adj
27 : public :: clybry_fam_init
28 :
29 : integer :: id_cly,id_bry
30 :
31 : integer :: id_cl,id_clo,id_hocl,id_cl2,id_cl2o2,id_oclo,id_hcl,id_clono2
32 : integer :: id_br,id_bro,id_hbr,id_brono2,id_brcl,id_hobr
33 :
34 : logical :: has_clybry
35 :
36 : contains
37 :
38 : !------------------------------------------
39 : !------------------------------------------
40 0 : subroutine clybry_fam_init
41 :
42 : use mo_chem_utls, only : get_spc_ndx
43 : implicit none
44 :
45 : integer :: ids(16)
46 :
47 0 : id_cly = get_spc_ndx('CLY')
48 0 : id_bry = get_spc_ndx('BRY')
49 :
50 0 : id_cl = get_spc_ndx('CL')
51 0 : id_clo = get_spc_ndx('CLO')
52 0 : id_hocl = get_spc_ndx('HOCL')
53 0 : id_cl2 = get_spc_ndx('CL2')
54 0 : id_cl2o2 = get_spc_ndx('CL2O2')
55 0 : id_oclo = get_spc_ndx('OCLO')
56 0 : id_hcl = get_spc_ndx('HCL')
57 0 : id_clono2 = get_spc_ndx('CLONO2')
58 :
59 0 : id_br = get_spc_ndx('BR')
60 0 : id_bro = get_spc_ndx('BRO')
61 0 : id_hbr = get_spc_ndx('HBR')
62 0 : id_brono2 = get_spc_ndx('BRONO2')
63 0 : id_brcl = get_spc_ndx('BRCL')
64 0 : id_hobr = get_spc_ndx('HOBR')
65 :
66 : ids = (/ id_cly,id_bry, &
67 : id_cl,id_clo,id_hocl,id_cl2,id_cl2o2,id_oclo,id_hcl,id_clono2, &
68 0 : id_br,id_bro,id_hbr,id_brono2,id_brcl,id_hobr /)
69 :
70 0 : has_clybry = all( ids(:) > 0 )
71 :
72 0 : endsubroutine clybry_fam_init
73 :
74 : !--------------------------------------------------------------
75 : ! set the ClOy and BrOy mass mixing ratios
76 : ! - this is call before advection
77 : !--------------------------------------------------------------
78 1489176 : subroutine clybry_fam_set( ncol, lchnk, map2chm, q, pbuf )
79 :
80 : use time_manager, only : get_nstep
81 : use physics_buffer, only : physics_buffer_desc
82 :
83 : implicit none
84 :
85 : !--------------------------------------------------------------
86 : ! ... dummy arguments
87 : !--------------------------------------------------------------
88 : integer, intent(in) :: ncol, lchnk
89 : integer, intent(in) :: map2chm(pcnst)
90 : real(r8), intent(inout) :: q(pcols,pver,pcnst)
91 : type(physics_buffer_desc), pointer :: pbuf(:)
92 :
93 2978352 : real(r8) :: wrk(ncol,pver,2)
94 : real(r8) :: mmr(pcols,pver,gas_pcnst)
95 : integer :: n, m
96 :
97 1489176 : if (.not. has_clybry) return
98 :
99 0 : do n = 1,pcnst
100 0 : m = map2chm(n)
101 0 : if( m > 0 ) then
102 0 : mmr(:ncol,:,m) = q(:ncol,:, n)
103 : endif
104 : enddo
105 0 : call get_short_lived_species( mmr, lchnk, ncol, pbuf )
106 :
107 : !--------------------------------------------------------------
108 : ! ... form updated chlorine, bromine atom mass mixing ratios
109 : !--------------------------------------------------------------
110 0 : wrk(:,:,1) = cloy( mmr, pcols, ncol )
111 : wrk(:,:,2) = broy( mmr, pcols, ncol )
112 :
113 : mmr(:ncol,:,id_cly) = wrk(:,:,1)
114 : mmr(:ncol,:,id_bry) = wrk(:,:,2)
115 :
116 : call set_short_lived_species( mmr, lchnk, ncol, pbuf )
117 : do n = 1,pcnst
118 : m = map2chm(n)
119 : if( m > 0 ) then
120 : q(:ncol,:, n) = mmr(:ncol,:,m)
121 : endif
122 : enddo
123 :
124 1489176 : end subroutine clybry_fam_set
125 :
126 : !--------------------------------------------------------------
127 : ! adjust the ClOy and BrOy individual family members
128 : ! - this is call after advection
129 : !--------------------------------------------------------------
130 1495368 : subroutine clybry_fam_adj( ncol, lchnk, map2chm, q, pbuf )
131 :
132 1489176 : use time_manager, only : is_first_step
133 : use physics_buffer, only : physics_buffer_desc
134 :
135 : implicit none
136 :
137 : !--------------------------------------------------------------
138 : ! ... dummy arguments
139 : !--------------------------------------------------------------
140 : integer, intent(in) :: ncol, lchnk
141 : integer, intent(in) :: map2chm(pcnst)
142 : real(r8), intent(inout) :: q(pcols,pver,pcnst)
143 : type(physics_buffer_desc), pointer :: pbuf(:)
144 :
145 : !--------------------------------------------------------------
146 : ! ... local variables
147 : !--------------------------------------------------------------
148 1495368 : real(r8) :: factor(ncol,pver)
149 2990736 : real(r8) :: wrk(ncol,pver)
150 : real(r8) :: mmr(pcols,pver,gas_pcnst)
151 :
152 : integer :: n, m
153 :
154 1495368 : if (.not. has_clybry) return
155 :
156 : !--------------------------------------------------------------
157 : ! ... CLY,BRY are not adjusted until the end of the first timestep
158 : !--------------------------------------------------------------
159 0 : if (is_first_step()) return
160 :
161 0 : do n = 1,pcnst
162 0 : m = map2chm(n)
163 0 : if( m > 0 ) then
164 0 : mmr(:ncol,:,m) = q(:ncol,:, n)
165 : endif
166 : enddo
167 0 : call get_short_lived_species( mmr, lchnk, ncol, pbuf )
168 :
169 : !--------------------------------------------------------------
170 : ! ... form updated chlorine atom mass mixing ratio
171 : !--------------------------------------------------------------
172 0 : wrk(:,:) = cloy( mmr, pcols, ncol )
173 :
174 : factor(:ncol,:) = mmr(:ncol,:,id_cly) / wrk(:ncol,:)
175 : !--------------------------------------------------------------
176 : ! ... adjust "group" members
177 : !--------------------------------------------------------------
178 : mmr(:ncol,:,id_cl) = factor(:ncol,:)*mmr(:ncol,:,id_cl)
179 : mmr(:ncol,:,id_clo) = factor(:ncol,:)*mmr(:ncol,:,id_clo)
180 : mmr(:ncol,:,id_hocl) = factor(:ncol,:)*mmr(:ncol,:,id_hocl)
181 : mmr(:ncol,:,id_cl2) = factor(:ncol,:)*mmr(:ncol,:,id_cl2)
182 : mmr(:ncol,:,id_cl2o2) = factor(:ncol,:)*mmr(:ncol,:,id_cl2o2)
183 : mmr(:ncol,:,id_oclo) = factor(:ncol,:)*mmr(:ncol,:,id_oclo)
184 : mmr(:ncol,:,id_hcl) = factor(:ncol,:)*mmr(:ncol,:,id_hcl)
185 : mmr(:ncol,:,id_clono2) = factor(:ncol,:)*mmr(:ncol,:,id_clono2)
186 :
187 : !--------------------------------------------------------------
188 : ! ... form updated bromine atom mass mixing ratio
189 : !--------------------------------------------------------------
190 : wrk(:,:) = broy( mmr, pcols, ncol )
191 :
192 : factor(:ncol,:) = mmr(:ncol,:,id_bry) / wrk(:ncol,:)
193 : !--------------------------------------------------------------
194 : ! ... adjust "group" members
195 : !--------------------------------------------------------------
196 : mmr(:ncol,:,id_br) = factor(:ncol,:)*mmr(:ncol,:,id_br)
197 : mmr(:ncol,:,id_bro) = factor(:ncol,:)*mmr(:ncol,:,id_bro)
198 : mmr(:ncol,:,id_hbr) = factor(:ncol,:)*mmr(:ncol,:,id_hbr)
199 : mmr(:ncol,:,id_brono2) = factor(:ncol,:)*mmr(:ncol,:,id_brono2)
200 : mmr(:ncol,:,id_brcl) = factor(:ncol,:)*mmr(:ncol,:,id_brcl)
201 : mmr(:ncol,:,id_hobr) = factor(:ncol,:)*mmr(:ncol,:,id_hobr)
202 :
203 : call set_short_lived_species( mmr, lchnk, ncol, pbuf )
204 : do n = 1,pcnst
205 : m = map2chm(n)
206 : if( m > 0 ) then
207 : q(:ncol,:, n) = mmr(:ncol,:,m)
208 : endif
209 : enddo
210 :
211 1495368 : end subroutine clybry_fam_adj
212 :
213 : !--------------------------------------------------------------
214 : ! private methods
215 : !--------------------------------------------------------------
216 :
217 : !--------------------------------------------------------------
218 : ! compute the mass mixing retio of ClOy
219 : !--------------------------------------------------------------
220 1495368 : function cloy( q, pcols, ncol )
221 :
222 : !--------------------------------------------------------------
223 : ! ... dummy arguments
224 : !--------------------------------------------------------------
225 : integer, intent(in) :: pcols
226 : integer, intent(in) :: ncol
227 : real(r8), intent(in) :: q(pcols,pver,gas_pcnst)
228 :
229 : !--------------------------------------------------------------
230 : ! ... function declaration
231 : !--------------------------------------------------------------
232 : real(r8) :: cloy(ncol,pver)
233 :
234 : !--------------------------------------------------------------
235 : ! ... local variables
236 : !--------------------------------------------------------------
237 0 : real(r8) :: wrk(ncol)
238 : integer :: k
239 :
240 : do k = 1,pver
241 0 : wrk(:) = q(:ncol,k,id_cl) /adv_mass(id_cl) &
242 : + q(:ncol,k,id_clo) /adv_mass(id_clo) &
243 : + q(:ncol,k,id_hocl) /adv_mass(id_hocl) &
244 : + 2._r8*( q(:ncol,k,id_cl2) /adv_mass(id_cl2) &
245 : + q(:ncol,k,id_cl2o2)/adv_mass(id_cl2o2) ) &
246 : + q(:ncol,k,id_oclo) /adv_mass(id_oclo) &
247 : + q(:ncol,k,id_hcl) /adv_mass(id_hcl) &
248 0 : + q(:ncol,k,id_clono2) /adv_mass(id_clono2)
249 : cloy(:,k) = adv_mass(id_cl) * wrk(:)
250 : end do
251 :
252 : end function cloy
253 :
254 : !--------------------------------------------------------------
255 : ! compute the mass mixing retio of BrOy
256 : !--------------------------------------------------------------
257 : function broy( q, pcols, ncol )
258 :
259 : !--------------------------------------------------------------
260 : ! ... dummy arguments
261 : !--------------------------------------------------------------
262 : integer, intent(in) :: pcols
263 : integer, intent(in) :: ncol
264 : real(r8), intent(in) :: q(pcols,pver,gas_pcnst)
265 :
266 : !--------------------------------------------------------------
267 : ! ... function declaration
268 : !--------------------------------------------------------------
269 : real(r8) :: broy(ncol,pver)
270 :
271 : !--------------------------------------------------------------
272 : ! ... local variables
273 : !--------------------------------------------------------------
274 : real(r8) :: wrk(ncol)
275 : integer :: k
276 :
277 : do k = 1,pver
278 : wrk(:) = q(:ncol,k,id_br) /adv_mass(id_br) &
279 : + q(:ncol,k,id_bro) /adv_mass(id_bro) &
280 : + q(:ncol,k,id_hbr) /adv_mass(id_hbr) &
281 : + q(:ncol,k,id_brono2)/adv_mass(id_brono2) &
282 : + q(:ncol,k,id_brcl) /adv_mass(id_brcl) &
283 : + q(:ncol,k,id_hobr) /adv_mass(id_hobr)
284 : broy(:,k) = adv_mass(id_br) * wrk(:)
285 : end do
286 :
287 : end function broy
288 :
289 : end module clybry_fam
|