Line data Source code
1 : module mo_chm_diags
2 :
3 : use shr_kind_mod, only : r8 => shr_kind_r8
4 : use chem_mods, only : gas_pcnst
5 : use mo_tracname, only : solsym
6 : use chem_mods, only : rxntot, nfs, gas_pcnst, indexm, adv_mass
7 : use ppgrid, only : pver
8 : use mo_constants, only : rgrav, rearth
9 : use mo_chem_utls, only : get_rxt_ndx, get_spc_ndx
10 : use cam_history, only : fieldname_len
11 : use mo_jeuv, only : neuv
12 : use gas_wetdep_opts,only : gas_wetdep_method
13 : use mo_drydep, only : has_drydep
14 :
15 : implicit none
16 : private
17 :
18 : public :: chm_diags_inti
19 : public :: chm_diags
20 : public :: het_diags
21 :
22 : integer :: id_n,id_no,id_no2,id_no3,id_n2o5,id_hno3,id_ho2no2,id_clono2,id_brono2
23 : integer :: id_isopfdn, id_isopfdnc, id_terpfdn !these are dinitrates
24 : integer :: id_cl,id_clo,id_hocl,id_cl2,id_cl2o2,id_oclo,id_hcl,id_brcl
25 : integer :: id_ccl4,id_cfc11,id_cfc113,id_ch3ccl3,id_cfc12,id_ch3cl,id_hcfc22,id_cf3br,id_cf2clbr
26 : integer :: id_cfc114,id_cfc115,id_hcfc141b,id_hcfc142b,id_h1202,id_h2402,id_ch2br2,id_chbr3
27 : integer :: id_hf,id_f,id_cof2,id_cofcl,id_ch3br
28 : integer :: id_br,id_bro,id_hbr,id_hobr,id_ch4,id_h2o,id_h2
29 : integer :: id_o,id_o2,id_h, id_h2o2, id_n2o
30 : integer :: id_co2,id_o3,id_oh,id_ho2,id_so4_a1,id_so4_a2,id_so4_a3
31 : integer :: id_num_a2,id_num_a3,id_dst_a3,id_ncl_a3
32 : integer :: id_ndep,id_nhdep
33 :
34 : integer, parameter :: NJEUV = neuv
35 : integer :: rid_jeuv(NJEUV), rid_jno_i, rid_jno
36 :
37 : logical :: has_jeuvs, has_jno_i, has_jno
38 :
39 : integer :: nox_species(3), noy_species(56)
40 : integer :: clox_species(6), cloy_species(9), tcly_species(21)
41 : integer :: brox_species(4), broy_species(6), tbry_species(13)
42 : integer :: foy_species(4), tfy_species(16)
43 : integer :: hox_species(4)
44 : integer :: toth_species(3)
45 : integer :: sox_species(3)
46 : integer :: nhx_species(2)
47 : integer :: aer_species(gas_pcnst)
48 :
49 : character(len=fieldname_len) :: dtchem_name(gas_pcnst)
50 : character(len=fieldname_len) :: depvel_name(gas_pcnst)
51 : character(len=fieldname_len) :: depflx_name(gas_pcnst)
52 : character(len=fieldname_len) :: wetdep_name(gas_pcnst)
53 : character(len=fieldname_len) :: wtrate_name(gas_pcnst)
54 :
55 : real(r8), parameter :: N_molwgt = 14.00674_r8
56 : real(r8), parameter :: S_molwgt = 32.066_r8
57 :
58 : contains
59 :
60 0 : subroutine chm_diags_inti
61 : !--------------------------------------------------------------------
62 : ! ... initialize utility routine
63 : !--------------------------------------------------------------------
64 :
65 : use cam_history, only : addfld, add_default, horiz_only
66 : use constituents, only : cnst_get_ind, cnst_longname
67 : use phys_control, only : phys_getopts
68 : use species_sums_diags, only : species_sums_init
69 :
70 : integer :: j, k, m, n
71 : character(len=16) :: jname, spc_name, attr
72 : character(len=2) :: jchar
73 : character(len=2) :: unit_basename ! Units 'kg' or '1'
74 :
75 : integer :: id_pan, id_onit, id_mpan, id_isopno3, id_onitr, id_nh4no3
76 : integer :: id_so2, id_so4, id_h2so4
77 : integer :: id_nh3, id_nh4
78 : integer :: id_honitr
79 : integer :: id_alknit
80 : integer :: id_isopnita
81 : integer :: id_isopnitb
82 : integer :: id_isopnooh
83 : integer :: id_nc4ch2oh
84 : integer :: id_nc4cho
85 : integer :: id_noa
86 : integer :: id_nterpooh
87 : integer :: id_pbznit
88 : integer :: id_terpnit
89 : integer :: id_dst01, id_dst02, id_dst03, id_dst04, id_sslt01, id_sslt02, id_sslt03, id_sslt04
90 : integer :: id_soa, id_oc1, id_oc2, id_cb1, id_cb2
91 : integer :: id_soam,id_soai,id_soat,id_soab,id_soax
92 : integer :: id_bry, id_cly
93 : integer :: id_isopn2b, id_isopn3b, id_isopn1d, id_isopn4d, id_isopnbno3
94 : integer :: id_isopfnp, id_isopnoohb, id_isopnoohd, id_inheb, id_inhed
95 : integer :: id_no3ch2cho, id_macrn, id_mvkn, id_isopfnc, id_terpns
96 : integer :: id_terpnt, id_terpnt1, id_terpns1, id_terpnpt, id_terpnps
97 : integer :: id_terpnpt1, id_terpnps1, id_sqtn, id_terphfn
98 : integer :: id_terpapan, id_terpa2pan, id_terpa3pan
99 :
100 : logical :: history_aerosol ! Output the MAM aerosol tendencies
101 : logical :: history_chemistry
102 : logical :: history_cesm_forcing
103 : logical :: history_scwaccm_forcing
104 : logical :: history_chemspecies_srf ! output the chemistry constituents species in the surface layer
105 : logical :: history_dust
106 : integer :: bulkaero_species(20)
107 :
108 : !-----------------------------------------------------------------------
109 :
110 : call phys_getopts( history_aerosol_out = history_aerosol, &
111 : history_chemistry_out = history_chemistry, &
112 : history_chemspecies_srf_out = history_chemspecies_srf, &
113 : history_cesm_forcing_out = history_cesm_forcing, &
114 : history_scwaccm_forcing_out = history_scwaccm_forcing, &
115 0 : history_dust_out = history_dust )
116 :
117 : id_bry = get_spc_ndx( 'BRY' )
118 : id_cly = get_spc_ndx( 'CLY' )
119 :
120 0 : id_n = get_spc_ndx( 'N' )
121 0 : id_no = get_spc_ndx( 'NO' )
122 0 : id_no2 = get_spc_ndx( 'NO2' )
123 0 : id_no3 = get_spc_ndx( 'NO3' )
124 0 : id_n2o5 = get_spc_ndx( 'N2O5' )
125 0 : id_n2o = get_spc_ndx( 'N2O' )
126 0 : id_hno3 = get_spc_ndx( 'HNO3' )
127 0 : id_ho2no2 = get_spc_ndx( 'HO2NO2' )
128 0 : id_clono2 = get_spc_ndx( 'CLONO2' )
129 0 : id_brono2 = get_spc_ndx( 'BRONO2' )
130 0 : id_cl = get_spc_ndx( 'CL' )
131 0 : id_clo = get_spc_ndx( 'CLO' )
132 0 : id_hocl = get_spc_ndx( 'HOCL' )
133 0 : id_cl2 = get_spc_ndx( 'CL2' )
134 0 : id_cl2o2 = get_spc_ndx( 'CL2O2' )
135 0 : id_oclo = get_spc_ndx( 'OCLO' )
136 0 : id_hcl = get_spc_ndx( 'HCL' )
137 0 : id_brcl = get_spc_ndx( 'BRCL' )
138 :
139 0 : id_co2 = get_spc_ndx( 'CO2' )
140 0 : id_o3 = get_spc_ndx( 'O3' )
141 0 : id_oh = get_spc_ndx( 'OH' )
142 0 : id_ho2 = get_spc_ndx( 'HO2' )
143 0 : id_h2o2 = get_spc_ndx( 'H2O2' )
144 0 : id_so4_a1 = get_spc_ndx( 'so4_a1' )
145 0 : id_so4_a2 = get_spc_ndx( 'so4_a2' )
146 0 : id_so4_a3 = get_spc_ndx( 'so4_a3' )
147 0 : id_num_a2 = get_spc_ndx( 'num_a2' )
148 0 : id_num_a3 = get_spc_ndx( 'num_a3' )
149 0 : id_dst_a3 = get_spc_ndx( 'dst_a3' )
150 0 : id_ncl_a3 = get_spc_ndx( 'ncl_a3' )
151 :
152 0 : id_f = get_spc_ndx( 'F' )
153 0 : id_hf = get_spc_ndx( 'HF' )
154 0 : id_cofcl = get_spc_ndx( 'COFCL' )
155 0 : id_cof2 = get_spc_ndx( 'COF2' )
156 :
157 0 : id_ccl4 = get_spc_ndx( 'CCL4' )
158 0 : id_cfc11 = get_spc_ndx( 'CFC11' )
159 :
160 0 : id_cfc113 = get_spc_ndx( 'CFC113' )
161 0 : id_cfc114 = get_spc_ndx( 'CFC114' )
162 0 : id_cfc115 = get_spc_ndx( 'CFC115' )
163 :
164 0 : id_ch3ccl3 = get_spc_ndx( 'CH3CCL3' )
165 0 : id_cfc12 = get_spc_ndx( 'CFC12' )
166 0 : id_ch3cl = get_spc_ndx( 'CH3CL' )
167 :
168 0 : id_hcfc22 = get_spc_ndx( 'HCFC22' )
169 0 : id_hcfc141b= get_spc_ndx( 'HCFC141B' )
170 0 : id_hcfc142b= get_spc_ndx( 'HCFC142B' )
171 :
172 0 : id_cf2clbr = get_spc_ndx( 'CF2CLBR' )
173 0 : id_cf3br = get_spc_ndx( 'CF3BR' )
174 0 : id_ch3br = get_spc_ndx( 'CH3BR' )
175 0 : id_h1202 = get_spc_ndx( 'H1202' )
176 0 : id_h2402 = get_spc_ndx( 'H2402' )
177 0 : id_ch2br2 = get_spc_ndx( 'CH2BR2' )
178 0 : id_chbr3 = get_spc_ndx( 'CHBR3' )
179 :
180 0 : id_br = get_spc_ndx( 'BR' )
181 0 : id_bro = get_spc_ndx( 'BRO' )
182 0 : id_hbr = get_spc_ndx( 'HBR' )
183 0 : id_hobr = get_spc_ndx( 'HOBR' )
184 0 : id_ch4 = get_spc_ndx( 'CH4' )
185 0 : id_h2o = get_spc_ndx( 'H2O' )
186 0 : id_h2 = get_spc_ndx( 'H2' )
187 0 : id_o = get_spc_ndx( 'O' )
188 0 : id_o2 = get_spc_ndx( 'O2' )
189 0 : id_h = get_spc_ndx( 'H' )
190 :
191 0 : id_pan = get_spc_ndx( 'PAN' )
192 0 : id_onit = get_spc_ndx( 'ONIT' )
193 0 : id_mpan = get_spc_ndx( 'MPAN' )
194 0 : id_isopno3 = get_spc_ndx( 'ISOPNO3' )
195 0 : id_onitr = get_spc_ndx( 'ONITR' )
196 0 : id_nh4no3 = get_spc_ndx( 'NH4NO3' )
197 :
198 0 : id_honitr = get_spc_ndx( 'HONITR' )
199 0 : id_alknit = get_spc_ndx( 'ALKNIT' )
200 0 : id_isopnita = get_spc_ndx( 'ISOPNITA' )
201 0 : id_isopnitb = get_spc_ndx( 'ISOPNITB' )
202 0 : id_isopnooh = get_spc_ndx( 'ISOPNOOH' )
203 0 : id_nc4ch2oh = get_spc_ndx( 'NC4CH2OH' )
204 0 : id_nc4cho = get_spc_ndx( 'NC4CHO' )
205 0 : id_noa = get_spc_ndx( 'NOA' )
206 0 : id_nterpooh = get_spc_ndx( 'NTERPOOH' )
207 0 : id_pbznit = get_spc_ndx( 'PBZNIT' )
208 0 : id_terpnit = get_spc_ndx( 'TERPNIT' )
209 0 : id_ndep = get_spc_ndx( 'NDEP' )
210 0 : id_nhdep = get_spc_ndx( 'NHDEP' )
211 :
212 0 : id_so2 = get_spc_ndx( 'SO2' )
213 0 : id_so4 = get_spc_ndx( 'SO4' )
214 0 : id_h2so4 = get_spc_ndx( 'H2SO4' )
215 :
216 0 : id_nh3 = get_spc_ndx( 'NH3' )
217 0 : id_nh4 = get_spc_ndx( 'NH4' )
218 0 : id_nh4no3 = get_spc_ndx( 'NH4NO3' )
219 :
220 0 : id_isopn2b = get_spc_ndx( 'ISOPN2B' )
221 0 : id_isopn3b = get_spc_ndx( 'ISOPN3B' )
222 0 : id_isopn1d = get_spc_ndx( 'ISOPN1D' )
223 0 : id_isopn4d = get_spc_ndx( 'ISOPN4D' )
224 0 : id_isopnbno3 = get_spc_ndx('ISOPNBNO3' )
225 0 : id_isopfdn = get_spc_ndx( 'ISOPFDN' )
226 0 : id_isopfdnc = get_spc_ndx( 'ISOPFDNC')
227 0 : id_terpfdn = get_spc_ndx( 'TERPFDN' )
228 0 : id_isopfnp = get_spc_ndx( 'ISOPFNP' )
229 0 : id_isopnoohb = get_spc_ndx( 'ISOPNOOHB' )
230 0 : id_isopnoohd = get_spc_ndx( 'ISOPNOOHD' )
231 0 : id_inheb = get_spc_ndx( 'INHEB' )
232 0 : id_inhed = get_spc_ndx( 'INHED' )
233 0 : id_no3ch2cho = get_spc_ndx( 'NO3CH2CHO' )
234 0 : id_macrn = get_spc_ndx( 'MACRN' )
235 0 : id_mvkn = get_spc_ndx( 'MVKN' )
236 0 : id_isopfnc = get_spc_ndx( 'ISOPFNC' )
237 0 : id_terpns = get_spc_ndx( 'TERPNS' )
238 0 : id_terpnt = get_spc_ndx( 'TERPNT' )
239 0 : id_terpnt1 = get_spc_ndx( 'TERPNT1' )
240 0 : id_terpns1 = get_spc_ndx( 'TERPNS1' )
241 0 : id_terpnpt = get_spc_ndx( 'TERPNPT' )
242 0 : id_terpnps = get_spc_ndx( 'TERPNPS' )
243 0 : id_terpnpt1 = get_spc_ndx( 'TERPNPT1' )
244 0 : id_terpnps1 = get_spc_ndx( 'TERPNPS1' )
245 0 : id_sqtn = get_spc_ndx( 'SQTN' )
246 0 : id_terphfn = get_spc_ndx( 'TERPHFN' )
247 0 : id_terpapan = get_spc_ndx( 'TERPAPAN' )
248 0 : id_terpa2pan = get_spc_ndx( 'TERPA2PAN' )
249 0 : id_terpa3pan = get_spc_ndx( 'TERPA3PAN' )
250 :
251 0 : id_dst01 = get_spc_ndx( 'DST01' )
252 0 : id_dst02 = get_spc_ndx( 'DST02' )
253 0 : id_dst03 = get_spc_ndx( 'DST03' )
254 0 : id_dst04 = get_spc_ndx( 'DST04' )
255 0 : id_sslt01 = get_spc_ndx( 'SSLT01' )
256 0 : id_sslt02 = get_spc_ndx( 'SSLT02' )
257 0 : id_sslt03 = get_spc_ndx( 'SSLT03' )
258 0 : id_sslt04 = get_spc_ndx( 'SSLT04' )
259 0 : id_soa = get_spc_ndx( 'SOA' )
260 0 : id_so4 = get_spc_ndx( 'SO4' )
261 0 : id_oc1 = get_spc_ndx( 'OC1' )
262 0 : id_oc2 = get_spc_ndx( 'OC2' )
263 0 : id_cb1 = get_spc_ndx( 'CB1' )
264 0 : id_cb2 = get_spc_ndx( 'CB2' )
265 :
266 0 : rid_jno = get_rxt_ndx( 'jno' )
267 0 : rid_jno_i = get_rxt_ndx( 'jno_i' )
268 :
269 0 : id_soam = get_spc_ndx( 'SOAM' )
270 0 : id_soai = get_spc_ndx( 'SOAI' )
271 0 : id_soat = get_spc_ndx( 'SOAT' )
272 0 : id_soab = get_spc_ndx( 'SOAB' )
273 0 : id_soax = get_spc_ndx( 'SOAX' )
274 :
275 :
276 : !... NOY species
277 0 : nox_species = (/ id_n, id_no, id_no2 /)
278 : noy_species = (/ id_n, id_no, id_no2, id_no3, id_n2o5, id_hno3, id_ho2no2, id_clono2, &
279 : id_brono2, id_pan, id_onit, id_mpan, id_isopno3, id_onitr, id_nh4no3, &
280 : id_honitr, id_alknit, id_isopnita, id_isopnitb, id_isopnooh, id_nc4ch2oh, &
281 : id_nc4cho, id_noa, id_nterpooh, id_pbznit, id_terpnit, &
282 : id_isopn2b, id_isopn3b, id_isopn1d, id_isopn4d, id_isopnbno3, &
283 : id_isopfdn, id_isopfdnc, id_terpfdn, &
284 : id_isopfnp, id_isopnoohb, id_isopnoohd, id_inheb, id_inhed, &
285 : id_no3ch2cho, id_macrn, id_mvkn, id_isopfnc, id_terpns, &
286 : id_terpnt, id_terpnt1, id_terpns1, id_terpnpt, id_terpnps, &
287 : id_terpnpt1, id_terpnps1, id_sqtn, id_terphfn, &
288 0 : id_terpapan, id_terpa2pan, id_terpa3pan /)
289 : !... HOX species
290 0 : hox_species = (/ id_h, id_oh, id_ho2, id_h2o2 /)
291 :
292 : !... CLOY species
293 0 : clox_species = (/ id_cl, id_clo, id_hocl, id_cl2, id_cl2o2, id_oclo /)
294 0 : cloy_species = (/ id_cl, id_clo, id_hocl, id_cl2, id_cl2o2, id_oclo, id_hcl, id_clono2, id_brcl /)
295 : tcly_species = (/ id_cl, id_clo, id_hocl, id_cl2, id_cl2o2, id_oclo, id_hcl, id_clono2, id_brcl, &
296 : id_ccl4, id_cfc11, id_cfc113, id_cfc114, id_cfc115, id_ch3ccl3, id_cfc12, id_ch3cl, &
297 0 : id_hcfc22, id_hcfc141b, id_hcfc142b, id_cf2clbr /)
298 :
299 : !... FOY species
300 0 : foy_species = (/ id_F, id_hf, id_cofcl, id_cof2 /)
301 : tfy_species = (/ id_f, id_hf, id_cofcl, id_cof2, id_cfc11, id_cfc12, id_cfc113, id_cfc114, id_cfc115, &
302 0 : id_hcfc22, id_hcfc141b, id_hcfc142b, id_cf2clbr, id_cf3br, id_h1202, id_h2402 /)
303 :
304 : !... BROY species
305 0 : brox_species = (/ id_br, id_bro, id_brcl, id_hobr /)
306 0 : broy_species = (/ id_br, id_bro, id_hbr, id_brono2, id_brcl, id_hobr /)
307 : tbry_species = (/ id_br, id_bro, id_hbr, id_brono2, id_brcl, id_hobr, id_cf2clbr, id_cf3br, id_ch3br, id_h1202, &
308 0 : id_h2402, id_ch2br2, id_chbr3 /)
309 :
310 0 : sox_species = (/ id_so2, id_so4, id_h2so4 /)
311 0 : nhx_species = (/ id_nh3, id_nh4 /)
312 0 : bulkaero_species(:) = -1
313 : bulkaero_species(1:20) = (/ id_dst01, id_dst02, id_dst03, id_dst04, &
314 : id_sslt01, id_sslt02, id_sslt03, id_sslt04, &
315 : id_soa, id_so4, id_oc1, id_oc2, id_cb1, id_cb2, id_nh4no3, &
316 0 : id_soam,id_soai,id_soat,id_soab,id_soax /)
317 :
318 : aer_species(:) = -1
319 0 : n = 1
320 : do m = 1,gas_pcnst
321 : k=0
322 : if ( any(bulkaero_species(:)==m) ) k=1
323 : if ( k==0 ) k = index(trim(solsym(m)), '_a')
324 : if ( k==0 ) k = index(trim(solsym(m)), '_c')
325 : if ( k>0 ) then ! must be aerosol species
326 : aer_species(n) = m
327 : n = n+1
328 : endif
329 : enddo
330 :
331 0 : toth_species = (/ id_ch4, id_h2o, id_h2 /)
332 :
333 0 : call addfld( 'NOX', (/ 'lev' /), 'A', 'mol/mol', 'nox (N+NO+NO2)' )
334 : call addfld( 'NOY', (/ 'lev' /), 'A', 'mol/mol', &
335 0 : 'noy = total nitrogen (N+NO+NO2+NO3+2N2O5+HNO3+HO2NO2+ORGNOY+NH4NO3)' )
336 0 : call addfld( 'NOY_SRF', horiz_only, 'A', 'mol/mol', 'surface noy volume mixing ratio' )
337 0 : call addfld( 'HOX', (/ 'lev' /), 'A', 'mol/mol', 'HOx (H+OH+HO2+2H2O2)' )
338 :
339 0 : call addfld( 'BROX', (/ 'lev' /), 'A', 'mol/mol', 'brox (Br+BrO+BRCl+HOBr)' )
340 0 : call addfld( 'BROY', (/ 'lev' /), 'A', 'mol/mol', 'total inorganic bromine (Br+BrO+HOBr+BrONO2+HBr+BrCl)' )
341 0 : call addfld( 'TBRY', (/ 'lev' /), 'A', 'mol/mol', 'total Br (ORG+INORG) volume mixing ratio' )
342 :
343 0 : call addfld( 'CLOX', (/ 'lev' /), 'A', 'mol/mol', 'clox (Cl+CLO+HOCl+2Cl2+2Cl2O2+OClO)' )
344 0 : call addfld( 'CLOY', (/ 'lev' /), 'A', 'mol/mol', 'total inorganic chlorine (Cl+ClO+2Cl2+2Cl2O2+OClO+HOCl+ClONO2+HCl+BrCl)' )
345 0 : call addfld( 'TCLY', (/ 'lev' /), 'A', 'mol/mol', 'total Cl (ORG+INORG) volume mixing ratio' )
346 :
347 0 : call addfld( 'FOY', (/ 'lev' /), 'A', 'mol/mol', 'total inorganic fluorine (F+HF+COFCL+2COF2)' )
348 0 : call addfld( 'TFY', (/ 'lev' /), 'A', 'mol/mol', 'total F (ORG+INORG) volume mixing ratio' )
349 :
350 0 : call addfld( 'TOTH', (/ 'lev' /), 'A', 'mol/mol', 'total H2 volume mixing ratio' )
351 :
352 0 : call addfld( 'NOY_mmr', (/ 'lev' /), 'A', 'kg/kg', 'NOy mass mixing ratio' )
353 0 : call addfld( 'SOX_mmr', (/ 'lev' /), 'A', 'kg/kg', 'SOx mass mixing ratio' )
354 0 : call addfld( 'NHX_mmr', (/ 'lev' /), 'A', 'kg/kg', 'NHx mass mixing ratio' )
355 :
356 0 : do j = 1,NJEUV
357 0 : write( jchar, '(I2)' ) j
358 0 : jname = 'jeuv_'//trim(adjustl(jchar))
359 0 : rid_jeuv(j) = get_rxt_ndx( trim(jname) )
360 : enddo
361 :
362 0 : has_jeuvs = all( rid_jeuv(:) > 0 )
363 0 : has_jno_i = rid_jno_i>0
364 0 : has_jno = rid_jno>0
365 :
366 0 : if ( has_jeuvs ) then
367 0 : call addfld( 'PION_EUV', (/ 'lev' /), 'I', '/cm^3/s', 'total euv ionization rate' )
368 0 : call addfld( 'PEUV1', (/ 'lev' /), 'I', '/cm^3/s', '(j1+j2+j3)*o' )
369 0 : call addfld( 'PEUV1e', (/ 'lev' /), 'I', '/cm^3/s', '(j14+j15+j16)*o' )
370 0 : call addfld( 'PEUV2', (/ 'lev' /), 'I', '/cm^3/s', 'j4*n' )
371 0 : call addfld( 'PEUV3', (/ 'lev' /), 'I', '/cm^3/s', '(j5+j7+j8+j9)*o2' )
372 0 : call addfld( 'PEUV3e', (/ 'lev' /), 'I', '/cm^3/s', '(j17+j19+j20+j21)*o2' )
373 0 : call addfld( 'PEUV4', (/ 'lev' /), 'I', '/cm^3/s', '(j10+j11)*n2' )
374 0 : call addfld( 'PEUV4e', (/ 'lev' /), 'I', '/cm^3/s', '(j22+j23)*n2' )
375 0 : call addfld( 'PEUVN2D', (/ 'lev' /), 'I', '/cm^3/s', '(j11+j13)*n2' )
376 0 : call addfld( 'PEUVN2De', (/ 'lev' /), 'I', '/cm^3/s', '(j23+j25)*n2' )
377 : endif
378 0 : if ( has_jno ) then
379 0 : call addfld( 'PJNO', (/ 'lev' /), 'I', '/cm^3/s', 'jno*no' )
380 : endif
381 0 : if ( has_jno_i ) then
382 0 : call addfld( 'PJNO_I', (/ 'lev' /), 'I', '/cm^3/s', 'jno_i*no' )
383 : endif
384 : !
385 : ! CCMI
386 : !
387 : call addfld( 'DO3CHM_TRP', horiz_only, 'A', 'kg/s', 'integrated net tendency from chem in troposphere', &
388 0 : flag_xyfill=.True. )
389 : call addfld( 'DO3CHM_LMS', horiz_only, 'A', 'kg/s', 'integrated net tendency from chem in lowermost stratosphere', &
390 0 : flag_xyfill=.True. )
391 : !
392 : do m = 1,gas_pcnst
393 :
394 : spc_name = trim(solsym(m))
395 :
396 : call cnst_get_ind(spc_name, n, abort=.false. )
397 : if ( n > 0 ) then
398 : attr = cnst_longname(n)
399 : elseif ( trim(spc_name) == 'H2O' ) then
400 : attr = 'water vapor'
401 : else
402 : attr = spc_name
403 : endif
404 :
405 : dtchem_name(m) = 'D'//trim(spc_name)//'CHM'
406 : call addfld( dtchem_name(m), (/ 'lev' /), 'A', 'kg/s', 'net tendency from chem' )
407 :
408 : depvel_name(m) = 'DV_'//trim(spc_name)
409 : depflx_name(m) = 'DF_'//trim(spc_name)
410 : if (has_drydep(spc_name)) then
411 : call addfld( depvel_name(m), horiz_only, 'A', 'cm/s', 'deposition velocity ' )
412 : call addfld( depflx_name(m), horiz_only, 'A', 'kg/m2/s', 'dry deposition flux ' )
413 : if (history_chemistry) then
414 : call add_default( depflx_name(m), 1, ' ' )
415 : endif
416 : endif
417 :
418 : if (gas_wetdep_method=='MOZ') then
419 : wetdep_name(m) = 'WD_'//trim(spc_name)
420 : wtrate_name(m) = 'WDR_'//trim(spc_name)
421 :
422 : call addfld( wetdep_name(m), horiz_only, 'A', 'kg/s', spc_name//' wet deposition' )
423 : call addfld( wtrate_name(m), (/ 'lev' /), 'A', '/s', spc_name//' wet deposition rate' )
424 : endif
425 :
426 : if (spc_name(1:3) == 'num') then
427 : unit_basename = ' 1'
428 : else
429 : unit_basename = 'kg'
430 : endif
431 :
432 : if ( any( aer_species == m ) ) then
433 : call addfld( spc_name, (/ 'lev' /), 'A', unit_basename//'/kg ', trim(attr)//' concentration')
434 : call addfld( trim(spc_name)//'_SRF', horiz_only, 'A', unit_basename//'/kg', trim(attr)//" in bottom layer")
435 : else
436 : call addfld( spc_name, (/ 'lev' /), 'A', 'mol/mol', trim(attr)//' concentration')
437 : call addfld( trim(spc_name)//'_SRF', horiz_only, 'A', 'mol/mol', trim(attr)//" in bottom layer")
438 : endif
439 :
440 : if ((m /= id_cly) .and. (m /= id_bry)) then
441 : if (history_aerosol.or.history_chemistry) then
442 : call add_default( spc_name, 1, ' ' )
443 : endif
444 : if (history_chemspecies_srf) then
445 : call add_default( trim(spc_name)//'_SRF', 1, ' ' )
446 : endif
447 : endif
448 :
449 : if ( history_cesm_forcing ) then
450 : if (m==id_o3) call add_default( spc_name, 1, ' ')
451 : if (m==id_oh) call add_default( spc_name, 1, ' ')
452 : if (m==id_no3) call add_default( spc_name, 1, ' ')
453 : if (m==id_ho2) call add_default( spc_name, 1, ' ')
454 :
455 : if (m==id_o3) call add_default( spc_name, 8, ' ')
456 : if (m==id_so4_a1) call add_default( spc_name, 8, ' ')
457 : if (m==id_so4_a2) call add_default( spc_name, 8, ' ')
458 : if (m==id_so4_a3) call add_default( spc_name, 8, ' ')
459 :
460 : if (m==id_num_a2) call add_default( spc_name, 8, ' ')
461 : if (m==id_num_a3) call add_default( spc_name, 8, ' ')
462 : if (m==id_dst_a3) call add_default( spc_name, 8, ' ')
463 : if (m==id_ncl_a3) call add_default( spc_name, 8, ' ')
464 :
465 : endif
466 : if ( history_scwaccm_forcing ) then
467 : if (m==id_co2) call add_default( spc_name, 8, ' ')
468 : if (m==id_h) call add_default( spc_name, 8, ' ')
469 : if (m==id_no) call add_default( spc_name, 8, ' ')
470 : if (m==id_o) call add_default( spc_name, 8, ' ')
471 : if (m==id_o2) call add_default( spc_name, 8, ' ')
472 : if (m==id_o3) call add_default( spc_name, 8, ' ')
473 : if (m==id_h2o) call add_default( spc_name, 1, ' ')
474 : if (m==id_ch4 ) call add_default( spc_name, 1, ' ')
475 : if (m==id_n2o ) call add_default( spc_name, 1, ' ')
476 : if (m==id_cfc11 ) call add_default( spc_name, 1, ' ')
477 : if (m==id_cfc12 ) call add_default( spc_name, 1, ' ')
478 : endif
479 :
480 : if (history_dust .and. (index(spc_name,'dst_') > 0)) call add_default( spc_name, 1, ' ')
481 :
482 : enddo
483 :
484 0 : call addfld( 'MASS', (/ 'lev' /), 'A', 'kg', 'mass of grid box' )
485 0 : call addfld( 'AREA', horiz_only, 'A', 'm2', 'area of grid box' )
486 :
487 0 : call addfld( 'dry_deposition_NOy_as_N', horiz_only, 'I', 'kg/m2/s', 'NOy dry deposition flux ' )
488 0 : call addfld( 'DF_SOX', horiz_only, 'I', 'kg/m2/s', 'SOx dry deposition flux ' )
489 0 : call addfld( 'dry_deposition_NHx_as_N', horiz_only, 'I', 'kg/m2/s', 'NHx dry deposition flux ' )
490 0 : if (gas_wetdep_method=='NEU') then
491 0 : call addfld( 'wet_deposition_NOy_as_N', horiz_only, 'A', 'kg/m2/s', 'NOy wet deposition' )
492 0 : call addfld( 'wet_deposition_NHx_as_N', horiz_only, 'A', 'kg/m2/s', 'NHx wet deposition' )
493 0 : elseif (gas_wetdep_method=='MOZ') then
494 0 : call addfld( 'wet_deposition_NOy_as_N', horiz_only, 'A', 'kg/s', 'NOy wet deposition' )
495 0 : call addfld( 'WD_SOX', horiz_only, 'A', 'kg/s', 'SOx wet deposition' )
496 0 : call addfld( 'wet_deposition_NHx_as_N', horiz_only, 'A', 'kg/s', 'NHx wet deposition' )
497 : endif
498 0 : if ( history_cesm_forcing ) then
499 0 : call add_default('dry_deposition_NOy_as_N', 1, ' ')
500 0 : call add_default('dry_deposition_NHx_as_N', 1, ' ')
501 0 : call add_default('wet_deposition_NOy_as_N', 1, ' ')
502 0 : call add_default('wet_deposition_NHx_as_N', 1, ' ')
503 : endif
504 :
505 0 : call species_sums_init()
506 :
507 0 : end subroutine chm_diags_inti
508 :
509 0 : subroutine chm_diags( lchnk, ncol, vmr, mmr, rxt_rates, invariants, depvel, depflx, mmr_tend, pdel, pmid, ltrop, &
510 0 : wetdepflx, nhx_nitrogen_flx, noy_nitrogen_flx )
511 : !--------------------------------------------------------------------
512 : ! ... utility routine to output chemistry diagnostic variables
513 : !--------------------------------------------------------------------
514 :
515 0 : use cam_history, only : outfld
516 : use phys_grid, only : get_area_all_p
517 : use species_sums_diags, only : species_sums_output
518 : !
519 : ! CCMI
520 : !
521 : use cam_history_support, only : fillvalue
522 :
523 : !--------------------------------------------------------------------
524 : ! ... dummy arguments
525 : !--------------------------------------------------------------------
526 : integer, intent(in) :: lchnk
527 : integer, intent(in) :: ncol
528 : real(r8), intent(in) :: vmr(ncol,pver,gas_pcnst)
529 : real(r8), intent(in) :: mmr(ncol,pver,gas_pcnst)
530 : real(r8), intent(in) :: rxt_rates(ncol,pver,rxntot)
531 : real(r8), intent(in) :: invariants(ncol,pver,max(1,nfs))
532 : real(r8), intent(in) :: depvel(ncol, gas_pcnst)
533 : real(r8), intent(in) :: depflx(ncol, gas_pcnst)
534 : real(r8), intent(in) :: mmr_tend(ncol,pver,gas_pcnst)
535 : real(r8), intent(in) :: pdel(ncol,pver)
536 : real(r8), intent(in) :: pmid(ncol,pver)
537 : integer, intent(in) :: ltrop(ncol)
538 : real(r8), intent(in) :: wetdepflx(ncol, gas_pcnst)
539 : real(r8), intent(out) :: nhx_nitrogen_flx(ncol) ! kgN/m2/sec
540 : real(r8), intent(out) :: noy_nitrogen_flx(ncol) ! kgN/m2/sec
541 :
542 : !--------------------------------------------------------------------
543 : ! ... local variables
544 : !--------------------------------------------------------------------
545 : integer :: i, k, m
546 0 : real(r8) :: wrk(ncol,pver)
547 : ! real(r8) :: tmp(ncol,pver)
548 : ! real(r8) :: m(ncol,pver)
549 0 : real(r8) :: un2(ncol)
550 :
551 0 : real(r8), dimension(ncol,pver) :: vmr_nox, vmr_noy, vmr_clox, vmr_cloy, vmr_tcly, vmr_brox, vmr_broy, vmr_toth
552 0 : real(r8), dimension(ncol,pver) :: vmr_tbry, vmr_foy, vmr_tfy
553 0 : real(r8), dimension(ncol,pver) :: mmr_noy, mmr_sox, mmr_nhx, net_chem
554 0 : real(r8), dimension(ncol) :: df_noy, df_sox, df_nhx, do3chm_trp, do3chm_lms
555 0 : real(r8), dimension(ncol) :: wd_noy, wd_nhx
556 0 : real(r8), dimension(ncol,pver) :: vmr_hox
557 :
558 0 : real(r8) :: area(ncol), mass(ncol,pver)
559 : real(r8) :: wgt
560 :
561 : !--------------------------------------------------------------------
562 : ! ... "diagnostic" groups
563 : !--------------------------------------------------------------------
564 0 : vmr_nox(:ncol,:) = 0._r8
565 0 : vmr_noy(:ncol,:) = 0._r8
566 0 : vmr_hox(:ncol,:) = 0._r8
567 0 : vmr_clox(:ncol,:) = 0._r8
568 0 : vmr_cloy(:ncol,:) = 0._r8
569 0 : vmr_tcly(:ncol,:) = 0._r8
570 0 : vmr_brox(:ncol,:) = 0._r8
571 0 : vmr_broy(:ncol,:) = 0._r8
572 0 : vmr_tbry(:ncol,:) = 0._r8
573 0 : vmr_foy(:ncol,:) = 0._r8
574 0 : vmr_tfy(:ncol,:) = 0._r8
575 0 : vmr_toth(:ncol,:) = 0._r8
576 0 : mmr_noy(:ncol,:) = 0._r8
577 0 : mmr_sox(:ncol,:) = 0._r8
578 0 : mmr_nhx(:ncol,:) = 0._r8
579 0 : df_noy(:ncol) = 0._r8
580 0 : df_sox(:ncol) = 0._r8
581 0 : df_nhx(:ncol) = 0._r8
582 :
583 0 : wd_noy(:ncol) = 0._r8
584 0 : wd_nhx(:ncol) = 0._r8
585 :
586 0 : call get_area_all_p(lchnk, ncol, area)
587 0 : area = area * rearth**2
588 :
589 0 : do k = 1,pver
590 0 : mass(:ncol,k) = pdel(:ncol,k) * area(:ncol) * rgrav
591 : enddo
592 :
593 0 : call outfld( 'AREA', area(:ncol), ncol, lchnk )
594 0 : call outfld( 'MASS', mass(:ncol,:), ncol, lchnk )
595 :
596 : do m = 1,gas_pcnst
597 :
598 : !...FOY (counting Fluorines, not chlorines or bromines)
599 : if ( m == id_cfc12 .or. m == id_hcfc22 .or. m == id_cf2clbr .or. m == id_h1202 .or. m == id_hcfc142b &
600 : .or. m == id_cof2 ) then
601 : wgt = 2._r8
602 : elseif ( m == id_cfc113 .or. m == id_cf3br ) then
603 : wgt = 3._r8
604 : elseif ( m == id_cfc114 .or. m == id_h2402 ) then
605 : wgt = 4._r8
606 : elseif ( m == id_cfc115 ) then
607 : wgt = 5._r8
608 : else
609 : wgt = 1._r8
610 : endif
611 : if ( any( foy_species == m ) ) then
612 : vmr_foy(:ncol,:) = vmr_foy(:ncol,:) + wgt * vmr(:ncol,:,m)
613 : endif
614 : if ( any( tfy_species == m ) ) then
615 : vmr_tfy(:ncol,:) = vmr_tfy(:ncol,:) + wgt * vmr(:ncol,:,m)
616 : endif
617 :
618 : !... counting chlorine and bromines, etc... (and total H2 species)
619 : if ( m == id_ch4 .or. m == id_n2o5 .or. m == id_cfc12 .or. m == id_cl2 .or. m == id_cl2o2 .or. m==id_h2o2 ) then
620 : wgt = 2._r8
621 : elseif (m == id_cfc114 .or. m == id_hcfc141b .or. m == id_h1202 .or. m == id_h2402 .or. m == id_ch2br2 ) then
622 : wgt = 2._r8
623 : elseif (m == id_isopfdn .or. m == id_isopfdnc .or. m == id_terpfdn ) then
624 : wgt = 2._r8
625 : elseif ( m == id_cfc11 .or. m == id_cfc113 .or. m == id_ch3ccl3 .or. m == id_chbr3 ) then
626 : wgt = 3._r8
627 : elseif ( m == id_ccl4 ) then
628 : wgt = 4._r8
629 : else
630 : wgt = 1._r8
631 : endif
632 : !...NOY
633 : if ( any( nox_species == m ) ) then
634 : vmr_nox(:ncol,:) = vmr_nox(:ncol,:) + wgt * vmr(:ncol,:,m)
635 : endif
636 : if ( any( noy_species == m ) ) then
637 : vmr_noy(:ncol,:) = vmr_noy(:ncol,:) + wgt * vmr(:ncol,:,m)
638 : endif
639 : !...NOY, SOX, NHX
640 : if ( any( noy_species == m ) ) then
641 : mmr_noy(:ncol,:) = mmr_noy(:ncol,:) + wgt * mmr(:ncol,:,m)
642 : endif
643 : if ( any( sox_species == m ) ) then
644 : mmr_sox(:ncol,:) = mmr_sox(:ncol,:) + wgt * mmr(:ncol,:,m)
645 : endif
646 : if ( any( nhx_species == m ) ) then
647 : mmr_nhx(:ncol,:) = mmr_nhx(:ncol,:) + wgt * mmr(:ncol,:,m)
648 : endif
649 : !...CLOY
650 : if ( any( clox_species == m ) ) then
651 : vmr_clox(:ncol,:) = vmr_clox(:ncol,:) + wgt * vmr(:ncol,:,m)
652 : endif
653 : if ( any( cloy_species == m ) ) then
654 : vmr_cloy(:ncol,:) = vmr_cloy(:ncol,:) + wgt * vmr(:ncol,:,m)
655 : endif
656 : if ( any( tcly_species == m ) ) then
657 : vmr_tcly(:ncol,:) = vmr_tcly(:ncol,:) + wgt * vmr(:ncol,:,m)
658 : endif
659 : !...BROY
660 : if ( any( brox_species == m ) ) then
661 : vmr_brox(:ncol,:) = vmr_brox(:ncol,:) + wgt * vmr(:ncol,:,m)
662 : endif
663 : if ( any( broy_species == m ) ) then
664 : vmr_broy(:ncol,:) = vmr_broy(:ncol,:) + wgt * vmr(:ncol,:,m)
665 : endif
666 : if ( any( tbry_species == m ) ) then
667 : vmr_tbry(:ncol,:) = vmr_tbry(:ncol,:) + wgt * vmr(:ncol,:,m)
668 : endif
669 : !...HOY
670 : if ( any ( toth_species == m ) ) then
671 : vmr_toth(:ncol,:) = vmr_toth(:ncol,:) + wgt * vmr(:ncol,:,m)
672 : endif
673 : !...HOx
674 : if ( any( hox_species == m ) ) then
675 : vmr_hox(:ncol,:) = vmr_hox(:ncol,:) + wgt * vmr(:ncol,:,m)
676 : endif
677 :
678 : if ( any( aer_species == m ) ) then
679 : call outfld( solsym(m), mmr(:ncol,:,m), ncol ,lchnk )
680 : call outfld( trim(solsym(m))//'_SRF', mmr(:ncol,pver,m), ncol ,lchnk )
681 : else
682 : call outfld( solsym(m), vmr(:ncol,:,m), ncol ,lchnk )
683 : call outfld( trim(solsym(m))//'_SRF', vmr(:ncol,pver,m), ncol ,lchnk )
684 : endif
685 :
686 : if (has_drydep(solsym(m))) then
687 : call outfld( depvel_name(m), depvel(:ncol,m), ncol ,lchnk )
688 : call outfld( depflx_name(m), depflx(:ncol,m), ncol ,lchnk )
689 : endif
690 :
691 : if ( any( noy_species == m ) ) then
692 : df_noy(:ncol) = df_noy(:ncol) + wgt * depflx(:ncol,m)*N_molwgt/adv_mass(m)
693 : endif
694 : if ( any( sox_species == m ) ) then
695 : df_sox(:ncol) = df_sox(:ncol) + wgt * depflx(:ncol,m)*S_molwgt/adv_mass(m)
696 : endif
697 : if ( any( nhx_species == m ) ) then
698 : df_nhx(:ncol) = df_nhx(:ncol) + wgt * depflx(:ncol,m)*N_molwgt/adv_mass(m)
699 : endif
700 :
701 : if ( any( noy_species == m ) ) then
702 : wd_noy(:ncol) = wd_noy(:ncol) + wgt * wetdepflx(:ncol,m)*N_molwgt/adv_mass(m)
703 : endif
704 : if ( any( nhx_species == m ) ) then
705 : wd_nhx(:ncol) = wd_nhx(:ncol) + wgt * wetdepflx(:ncol,m)*N_molwgt/adv_mass(m)
706 : endif
707 : !
708 : ! add contribution from non-conservation tracers
709 : !
710 : if ( id_ndep == m ) then
711 : wd_noy(:ncol) = wd_noy(:ncol) + wgt * wetdepflx(:ncol,m)*N_molwgt/adv_mass(m)
712 : end if
713 : if ( id_nhdep == m ) then
714 : wd_nhx(:ncol) = wd_nhx(:ncol) + wgt * wetdepflx(:ncol,m)*N_molwgt/adv_mass(m)
715 : end if
716 :
717 : do k=1,pver
718 : do i=1,ncol
719 : net_chem(i,k) = mmr_tend(i,k,m) * mass(i,k)
720 : end do
721 : end do
722 : call outfld( dtchem_name(m), net_chem(:ncol,:), ncol, lchnk )
723 : !
724 : ! CCMI
725 : !
726 : if ( trim(dtchem_name(m)) == 'DO3CHM' ) then
727 : do3chm_trp(:) = 0._r8
728 : do i=1,ncol
729 : do k=ltrop(i),pver
730 : do3chm_trp(i) = do3chm_trp(i) + net_chem(i,k)
731 : end do
732 : end do
733 : where ( do3chm_trp == 0._r8 )
734 : do3chm_trp = fillvalue
735 : end where
736 : call outfld('DO3CHM_TRP',do3chm_trp(:ncol), ncol, lchnk )
737 : do3chm_lms(:) = 0._r8
738 : do i=1,ncol
739 : do k=1,pver
740 : if ( pmid(i,k) > 100.e2_r8 .and. k < ltrop(i) ) then
741 : do3chm_lms(i) = do3chm_lms(i) + net_chem(i,k)
742 : end if
743 : end do
744 : end do
745 : where ( do3chm_lms == 0._r8 )
746 : do3chm_lms = fillvalue
747 : end where
748 : call outfld('DO3CHM_LMS',do3chm_lms(:ncol), ncol, lchnk )
749 : end if
750 : !
751 : enddo
752 :
753 :
754 0 : call outfld( 'NOX', vmr_nox (:ncol,:), ncol, lchnk )
755 0 : call outfld( 'NOY', vmr_noy (:ncol,:), ncol, lchnk )
756 0 : call outfld( 'HOX', vmr_hox (:ncol,:), ncol, lchnk )
757 0 : call outfld( 'NOY_SRF', vmr_noy(:ncol,pver), ncol, lchnk )
758 0 : call outfld( 'CLOX', vmr_clox (:ncol,:), ncol, lchnk )
759 0 : call outfld( 'CLOY', vmr_cloy (:ncol,:), ncol, lchnk )
760 0 : call outfld( 'BROX', vmr_brox (:ncol,:), ncol, lchnk )
761 0 : call outfld( 'BROY', vmr_broy (:ncol,:), ncol, lchnk )
762 0 : call outfld( 'TCLY', vmr_tcly (:ncol,:), ncol, lchnk )
763 0 : call outfld( 'TBRY', vmr_tbry (:ncol,:), ncol, lchnk )
764 0 : call outfld( 'FOY', vmr_foy (:ncol,:), ncol, lchnk )
765 0 : call outfld( 'TFY', vmr_tfy (:ncol,:), ncol, lchnk )
766 0 : call outfld( 'TOTH', vmr_toth (:ncol,:), ncol, lchnk )
767 :
768 0 : call outfld( 'NOY_mmr', mmr_noy(:ncol,:), ncol ,lchnk )
769 0 : call outfld( 'SOX_mmr', mmr_sox(:ncol,:), ncol ,lchnk )
770 0 : call outfld( 'NHX_mmr', mmr_nhx(:ncol,:), ncol ,lchnk )
771 0 : call outfld( 'dry_deposition_NOy_as_N', df_noy(:ncol), ncol ,lchnk )
772 0 : call outfld( 'DF_SOX', df_sox(:ncol), ncol ,lchnk )
773 0 : call outfld( 'dry_deposition_NHx_as_N', df_nhx(:ncol), ncol ,lchnk )
774 0 : if (gas_wetdep_method=='NEU') then
775 0 : wd_noy(:ncol) = -wd_noy(:ncol) ! downward is possitive
776 0 : wd_nhx(:ncol) = -wd_nhx(:ncol)
777 0 : call outfld( 'wet_deposition_NOy_as_N', wd_noy(:ncol), ncol, lchnk )
778 0 : call outfld( 'wet_deposition_NHx_as_N', wd_nhx(:ncol), ncol, lchnk )
779 : end if
780 :
781 0 : nhx_nitrogen_flx = df_nhx + wd_nhx
782 0 : noy_nitrogen_flx = df_noy + wd_noy
783 :
784 : !--------------------------------------------------------------------
785 : ! ... euv ion production
786 : !--------------------------------------------------------------------
787 :
788 0 : jeuvs: if ( has_jeuvs ) then
789 : do k = 1,pver
790 0 : un2(:) = 1._r8 - (vmr(:,k,id_o) + vmr(:,k,id_o2) + vmr(:,k,id_h))
791 : wrk(:,k) = vmr(:,k,id_o)*(rxt_rates(:,k,rid_jeuv(1)) + rxt_rates(:,k,rid_jeuv(2)) &
792 : + rxt_rates(:,k,rid_jeuv(3)) + rxt_rates(:,k,rid_jeuv(14)) &
793 : + rxt_rates(:,k,rid_jeuv(15)) + rxt_rates(:,k,rid_jeuv(16))) &
794 : + vmr(:,k,id_n)*rxt_rates(:,k,rid_jeuv(4)) &
795 : + vmr(:,k,id_o2)*(rxt_rates(:,k,rid_jeuv(5)) + rxt_rates(:,k,rid_jeuv(7)) &
796 : + rxt_rates(:,k,rid_jeuv(8)) + rxt_rates(:,k,rid_jeuv(9)) &
797 : + rxt_rates(:,k,rid_jeuv(17)) + rxt_rates(:,k,rid_jeuv(19)) &
798 : + rxt_rates(:,k,rid_jeuv(20)) + rxt_rates(:,k,rid_jeuv(21))) &
799 : + un2(:)*(rxt_rates(:,k,rid_jeuv(6)) + rxt_rates(:,k,rid_jeuv(10)) &
800 : + rxt_rates(:,k,rid_jeuv(11)) + rxt_rates(:,k,rid_jeuv(18)) &
801 : + rxt_rates(:,k,rid_jeuv(22)) + rxt_rates(:,k,rid_jeuv(23)))
802 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
803 : end do
804 : call outfld( 'PION_EUV', wrk, ncol, lchnk )
805 :
806 : do k = 1,pver
807 : wrk(:,k) = vmr(:,k,id_o)*(rxt_rates(:,k,rid_jeuv(1)) + rxt_rates(:,k,rid_jeuv(2)) &
808 : + rxt_rates(:,k,rid_jeuv(3)))
809 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
810 : end do
811 : call outfld( 'PEUV1', wrk, ncol, lchnk )
812 : do k = 1,pver
813 : wrk(:,k) = vmr(:,k,id_o)*(rxt_rates(:,k,rid_jeuv(14)) + rxt_rates(:,k,rid_jeuv(15)) &
814 : + rxt_rates(:,k,rid_jeuv(16)))
815 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
816 : end do
817 : call outfld( 'PEUV1e', wrk, ncol, lchnk )
818 : do k = 1,pver
819 : wrk(:,k) = vmr(:,k,id_n)*rxt_rates(:,k,rid_jeuv(4))
820 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
821 : end do
822 : call outfld( 'PEUV2', wrk, ncol, lchnk )
823 : do k = 1,pver
824 : wrk(:,k) = vmr(:,k,id_o2)*(rxt_rates(:,k,rid_jeuv(5)) + rxt_rates(:,k,rid_jeuv(7)) &
825 : + rxt_rates(:,k,rid_jeuv(8)) + rxt_rates(:,k,rid_jeuv(9)))
826 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
827 : end do
828 : call outfld( 'PEUV3', wrk, ncol, lchnk )
829 : do k = 1,pver
830 : wrk(:,k) = vmr(:,k,id_o2)*(rxt_rates(:,k,rid_jeuv(17)) + rxt_rates(:,k,rid_jeuv(19)) &
831 : + rxt_rates(:,k,rid_jeuv(20)) + rxt_rates(:,k,rid_jeuv(21)))
832 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
833 : end do
834 : call outfld( 'PEUV3e', wrk, ncol, lchnk )
835 : do k = 1,pver
836 : un2(:) = 1._r8 - (vmr(:,k,id_o) + vmr(:,k,id_o2) + vmr(:,k,id_h))
837 : wrk(:,k) = un2(:)*(rxt_rates(:,k,rid_jeuv(6)) + rxt_rates(:,k,rid_jeuv(10)) + rxt_rates(:,k,rid_jeuv(11)))
838 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
839 : end do
840 : call outfld( 'PEUV4', wrk, ncol, lchnk )
841 : do k = 1,pver
842 : un2(:) = 1._r8 - (vmr(:,k,id_o) + vmr(:,k,id_o2) + vmr(:,k,id_h))
843 : wrk(:,k) = un2(:)*(rxt_rates(:,k,rid_jeuv(18)) + rxt_rates(:,k,rid_jeuv(22)) + rxt_rates(:,k,rid_jeuv(23)))
844 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
845 : end do
846 : call outfld( 'PEUV4e', wrk, ncol, lchnk )
847 : do k = 1,pver
848 : un2(:) = 1._r8 - (vmr(:,k,id_o) + vmr(:,k,id_o2) + vmr(:,k,id_h))
849 : wrk(:,k) = un2(:)*(rxt_rates(:,k,rid_jeuv(11)) + rxt_rates(:,k,rid_jeuv(13)))
850 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
851 : end do
852 : call outfld( 'PEUVN2D', wrk, ncol, lchnk )
853 : do k = 1,pver
854 : un2(:) = 1._r8 - (vmr(:,k,id_o) + vmr(:,k,id_o2) + vmr(:,k,id_h))
855 : wrk(:,k) = un2(:)*(rxt_rates(:,k,rid_jeuv(23)) + rxt_rates(:,k,rid_jeuv(25)))
856 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
857 : end do
858 : call outfld( 'PEUVN2De', wrk, ncol, lchnk )
859 : endif jeuvs
860 :
861 0 : if ( has_jno_i ) then
862 : do k = 1,pver
863 0 : wrk(:,k) = vmr(:,k,id_no)*rxt_rates(:,k,rid_jno_i)
864 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
865 : end do
866 : call outfld( 'PJNO_I', wrk, ncol, lchnk )
867 : endif
868 0 : if ( has_jno ) then
869 : do k = 1,pver
870 0 : wrk(:,k) = vmr(:,k,id_no)*rxt_rates(:,k,rid_jno)
871 : wrk(:,k) = wrk(:,k) * invariants(:,k,indexm)
872 : end do
873 : call outfld( 'PJNO', wrk, ncol, lchnk )
874 : endif
875 :
876 0 : call species_sums_output(vmr, mmr, ncol, lchnk)
877 :
878 0 : end subroutine chm_diags
879 :
880 0 : subroutine het_diags( het_rates, mmr, pdel, lchnk, ncol )
881 :
882 0 : use cam_history, only : outfld
883 : use phys_grid, only : get_wght_all_p
884 :
885 : integer, intent(in) :: lchnk
886 : integer, intent(in) :: ncol
887 : real(r8), intent(in) :: het_rates(ncol,pver,max(1,gas_pcnst))
888 : real(r8), intent(in) :: mmr(ncol,pver,gas_pcnst)
889 : real(r8), intent(in) :: pdel(ncol,pver)
890 :
891 0 : real(r8), dimension(ncol) :: noy_wk, sox_wk, nhx_wk, wrk_wd
892 : integer :: m, k
893 0 : real(r8) :: wght(ncol)
894 : !
895 : ! output integrated wet deposition field
896 : !
897 0 : noy_wk(:) = 0._r8
898 0 : sox_wk(:) = 0._r8
899 0 : nhx_wk(:) = 0._r8
900 :
901 0 : call get_wght_all_p(lchnk, ncol, wght)
902 :
903 : do m = 1,gas_pcnst
904 : !
905 : ! compute vertical integral
906 : !
907 : wrk_wd(:ncol) = 0._r8
908 : do k = 1,pver
909 : wrk_wd(:ncol) = wrk_wd(:ncol) + het_rates(:ncol,k,m) * mmr(:ncol,k,m) * pdel(:ncol,k)
910 : end do
911 : !
912 : wrk_wd(:ncol) = wrk_wd(:ncol) * rgrav * wght(:ncol) * rearth**2
913 : !
914 : if (gas_wetdep_method=='MOZ') then
915 : call outfld( wetdep_name(m), wrk_wd(:ncol), ncol, lchnk )
916 : call outfld( wtrate_name(m), het_rates(:ncol,:,m), ncol, lchnk )
917 :
918 : if ( any(noy_species == m ) ) then
919 : noy_wk(:ncol) = noy_wk(:ncol) + wrk_wd(:ncol)*N_molwgt/adv_mass(m)
920 : endif
921 : if ( m == id_n2o5 ) then ! 2 NOy molecules in N2O5
922 : noy_wk(:ncol) = noy_wk(:ncol) + wrk_wd(:ncol)*N_molwgt/adv_mass(m)
923 : endif
924 : if ( any(sox_species == m ) ) then
925 : sox_wk(:ncol) = sox_wk(:ncol) + wrk_wd(:ncol)*S_molwgt/adv_mass(m)
926 : endif
927 : if ( any(nhx_species == m ) ) then
928 : nhx_wk(:ncol) = nhx_wk(:ncol) + wrk_wd(:ncol)*N_molwgt/adv_mass(m)
929 : endif
930 : endif
931 : end do
932 0 : if (gas_wetdep_method=='MOZ') then
933 0 : call outfld( 'wet_deposition_NOy_as_N', noy_wk(:ncol), ncol, lchnk )
934 0 : call outfld( 'WD_SOX', sox_wk(:ncol), ncol, lchnk )
935 0 : call outfld( 'wet_deposition_NHx_as_N', nhx_wk(:ncol), ncol, lchnk )
936 : endif
937 :
938 0 : end subroutine het_diags
939 :
940 : end module mo_chm_diags
|