1 | !$Id: chemmain.F90 125 2009-04-21 10:22:57Z acosce $ |
---|
2 | !! ========================================================================= |
---|
3 | !! INCA - INteraction with Chemistry and Aerosols |
---|
4 | !! |
---|
5 | !! Copyright Laboratoire des Sciences du Climat et de l'Environnement (LSCE) |
---|
6 | !! Unite mixte CEA-CNRS-UVSQ |
---|
7 | !! |
---|
8 | !! Contributors to this INCA subroutine: |
---|
9 | !! |
---|
10 | !! Didier Hauglustaine, LSCE, hauglustaine@cea.fr |
---|
11 | !! Stacy Walters, NCAR, stacy@ucar.edu |
---|
12 | !! C. Textor, LSCE |
---|
13 | !! |
---|
14 | !! Anne Cozic, LSCE, anne.cozic@cea.fr |
---|
15 | !! Yann Meurdesoif, LSCE, yann.meurdesoif@cea.fr |
---|
16 | !! |
---|
17 | !! This software is a computer program whose purpose is to simulate the |
---|
18 | !! atmospheric gas phase and aerosol composition. The model is designed to be |
---|
19 | !! used within a transport model or a general circulation model. This version |
---|
20 | !! of INCA was designed to be coupled to the LMDz GCM. LMDz-INCA accounts |
---|
21 | !! for emissions, transport (resolved and sub-grid scale), photochemical |
---|
22 | !! transformations, and scavenging (dry deposition and washout) of chemical |
---|
23 | !! species and aerosols interactively in the GCM. Several versions of the INCA |
---|
24 | !! model are currently used depending on the envisaged applications with the |
---|
25 | !! chemistry-climate model. |
---|
26 | !! |
---|
27 | !! This software is governed by the CeCILL license under French law and |
---|
28 | !! abiding by the rules of distribution of free software. You can use, |
---|
29 | !! modify and/ or redistribute the software under the terms of the CeCILL |
---|
30 | !! license as circulated by CEA, CNRS and INRIA at the following URL |
---|
31 | !! "http://www.cecill.info". |
---|
32 | !! |
---|
33 | !! As a counterpart to the access to the source code and rights to copy, |
---|
34 | !! modify and redistribute granted by the license, users are provided only |
---|
35 | !! with a limited warranty and the software's author, the holder of the |
---|
36 | !! economic rights, and the successive licensors have only limited |
---|
37 | !! liability. |
---|
38 | !! |
---|
39 | !! In this respect, the user's attention is drawn to the risks associated |
---|
40 | !! with loading, using, modifying and/or developing or reproducing the |
---|
41 | !! software by the user in light of its specific status of free software, |
---|
42 | !! that may mean that it is complicated to manipulate, and that also |
---|
43 | !! therefore means that it is reserved for developers and experienced |
---|
44 | !! professionals having in-depth computer knowledge. Users are therefore |
---|
45 | !! encouraged to load and test the software's suitability as regards their |
---|
46 | !! requirements in conditions enabling the security of their systems and/or |
---|
47 | !! data to be ensured and, more generally, to use and operate it in the |
---|
48 | !! same conditions as regards security. |
---|
49 | !! |
---|
50 | !! The fact that you are presently reading this means that you have had |
---|
51 | !! knowledge of the CeCILL license and that you accept its terms. |
---|
52 | !! ========================================================================= |
---|
53 | |
---|
54 | #include <inca_define.h> |
---|
55 | |
---|
56 | |
---|
57 | SUBROUTINE CHEMMAIN( & |
---|
58 | mmr , & |
---|
59 | nstep , & |
---|
60 | calday , & |
---|
61 | ncdate , & !not used |
---|
62 | ncsec , & !not used |
---|
63 | lat , & |
---|
64 | delt , & |
---|
65 | ps , & |
---|
66 | pmid , & |
---|
67 | pdel , & |
---|
68 | area , & |
---|
69 | oro , & |
---|
70 | tsurf , & ! not used |
---|
71 | albs , & |
---|
72 | zma , & |
---|
73 | phis , & |
---|
74 | cldfr , & |
---|
75 | cldfr_st , & |
---|
76 | cldfr_cv , & |
---|
77 | cldtop , & |
---|
78 | cldbot , & |
---|
79 | cwat , & |
---|
80 | flxrst , & |
---|
81 | flxrcv , & |
---|
82 | flxsst , & |
---|
83 | flxscv , & |
---|
84 | flxupd , & |
---|
85 | flxmass_w , & |
---|
86 | tfld , & |
---|
87 | sh , & |
---|
88 | ql , & ! variable en prevision de la chimie strato |
---|
89 | rh , & |
---|
90 | nx , & |
---|
91 | ny , & |
---|
92 | source) |
---|
93 | |
---|
94 | !----------------------------------------------------------------------- |
---|
95 | ! |
---|
96 | ! INCA -- INteractions with Chemistry in the Atmosphere |
---|
97 | ! |
---|
98 | ! Advances the volumetric mixing ratio forward one time step via a |
---|
99 | ! combination of explicit, EBI, QSSA, fully implicit, and/or rodas |
---|
100 | ! algorithms. |
---|
101 | ! |
---|
102 | ! Didier Hauglustaine and Stacy Walters, 2000. |
---|
103 | !----------------------------------------------------------------------- |
---|
104 | |
---|
105 | USE CHEM_MODS, ONLY : nas, invariants, hrates, hrates_cv, o3_prod, o3_loss, o3_st_flx, & |
---|
106 | so2_p_h2soh, so2_p_dmsoh, so2_p_dmsno3, so2_p_dmsooh, asmsam_p_dmsooh, & |
---|
107 | dmso_p_dmsoh, asso4m_p_so2oh, wet3d_so2, wet3d_dms, wet3d_dmso, wet3d_hno3, wet3d_nh3, & |
---|
108 | wet3d_noy, wet3d_h2o2, wet3d_hono, hno3_p_g, hno3_p_a, hno3_l_g, hno3_l_a, & |
---|
109 | nh3_l_g, nh3_l_a, csno3_p_a1, csno3_p_a2, cino3_p_a, & |
---|
110 | wet3d_asso4m, wet3d_asnh4m, wet3d_asno3m, wet3d_csno3m, wet3d_cino3m, & |
---|
111 | hono_p_g, hono_p_a, hono_l_g, hono_l_a, & |
---|
112 | asno3m_p_nh3hno3, asnh4m_p_nh3hno3, hno3_p_nh3hno3, nh3_p_nh3hno3, & |
---|
113 | ASAPp1a_p, ASAPp2a_p, ASARp1a_p, ASARp2a_p, pom_p_g, & |
---|
114 | co2_basprod, co2_nmhcprod, co2_radicalprod, & |
---|
115 | hno3_prod, hno3_loss, co_prod, co_loss, ch4_loss, n2o_loss, nadv_mass, & |
---|
116 | no_daytime, day_cnt , adv_mass, h2o, h2oc, extfrc, wetloss |
---|
117 | |
---|
118 | USE SPECIES_NAMES |
---|
119 | USE CHEM_TRACNM |
---|
120 | USE CHEM_CONS, ONLY : gravit, uma |
---|
121 | ! USE MASS_DIAGS |
---|
122 | USE AIRPLANE_SRC, ONLY : itrop |
---|
123 | USE SFLX |
---|
124 | USE TRANSPORT_CONTROLS |
---|
125 | USE INCA_DIM |
---|
126 | USE MOD_INCA_PARA |
---|
127 | |
---|
128 | USE MOD_GRID_INCA, ONLY : PLON_GLO |
---|
129 | USE RATE_INDEX_MOD |
---|
130 | |
---|
131 | USE INCA_WRITE_FIELD_P |
---|
132 | USE XIOS_INCA |
---|
133 | |
---|
134 | USE OXYDANT_COM |
---|
135 | |
---|
136 | USE PARAM_CHEM |
---|
137 | USE LIGHTNING |
---|
138 | |
---|
139 | #ifdef STRAT |
---|
140 | USE HETCHEM |
---|
141 | USE LGLIVED_MOD |
---|
142 | #endif |
---|
143 | |
---|
144 | IMPLICIT NONE |
---|
145 | |
---|
146 | !----------------------------------------------------------------------- |
---|
147 | ! ... Dummy arguments |
---|
148 | !----------------------------------------------------------------------- |
---|
149 | INTEGER, INTENT(in) :: nstep ! time index |
---|
150 | INTEGER, INTENT(in) :: lat ! latitude index (always 1 and place holder in LMDz) |
---|
151 | INTEGER, INTENT(in) :: nx, ny ! dyn grid |
---|
152 | INTEGER, INTENT(in) :: ncdate ! date at end present time step |
---|
153 | INTEGER, INTENT(in) :: ncsec ! seconds relative to ncdate |
---|
154 | INTEGER, INTENT(in) :: cldtop(PLON) ! cloud top level ( 1 ... PLEV ) |
---|
155 | INTEGER, INTENT(in) :: cldbot(PLON) ! cloud bot level ( 1 ... PLEV ) |
---|
156 | REAL, INTENT(in) :: calday ! time of year in days for midpoint time |
---|
157 | REAL, INTENT(in) :: delt ! timestep in seconds |
---|
158 | REAL, INTENT(inout) :: mmr(PLON,PLEV,PCNST) ! xported species ( mmr ) |
---|
159 | REAL, INTENT(in) :: ps(PLON) ! surface press ( pascals ) |
---|
160 | REAL, INTENT(in) :: pmid(PLON,PLEV) ! midpoint press ( pascals ) |
---|
161 | REAL, INTENT(in) :: pdel(PLON,PLEV) ! delta press across midpoints |
---|
162 | REAL, INTENT(in) :: area(PLON) ! grid cell area |
---|
163 | REAL, INTENT(in) :: oro(PLON) ! surface orography flag |
---|
164 | REAL, INTENT(in) :: tsurf(PLON) ! surface temperature |
---|
165 | REAL, INTENT(in) :: zma(PLON,PLEV) ! abs geopot height at midpoints ( m ) |
---|
166 | REAL, INTENT(in) :: phis(PLON) ! surf geopot |
---|
167 | REAL, INTENT(in) :: cldfr(PLON,PLEV) ! cloud fraction (all clouds) |
---|
168 | REAL, INTENT(in) :: cldfr_cv(PLON,PLEV) ! cloud fraction (convective clouds) |
---|
169 | REAL, INTENT(in) :: cldfr_st(PLON,PLEV) ! cloud fraction (large scale clouds) |
---|
170 | REAL, INTENT(in) :: cwat(PLON,PLEV) ! total cloud water (kg/kg) |
---|
171 | REAL, INTENT(in) :: flxrst(PLON,PLEVP) !liquid water flux (stratiform) kgH2O/m2/s |
---|
172 | REAL, INTENT(in) :: flxrcv(PLON,PLEVP) !solid water flux (stratiform) kgH2O/m2/s |
---|
173 | REAL, INTENT(in) :: flxsst(PLON,PLEVP) !liquid water flux (convection) kgH2O/m2/s |
---|
174 | REAL, INTENT(in) :: flxscv(PLON,PLEVP) !solid water flux (convection) kgH2O/m2/s |
---|
175 | REAL, INTENT(in) :: flxmass_w(PLON,PLEV) ! vertical mass flux |
---|
176 | REAL, INTENT(in) :: flxupd(PLON,PLEV) ! entrainment flux kgAIR/m2/s |
---|
177 | REAL, INTENT(in) :: tfld(PLON,PLEV) ! midpoint temperature |
---|
178 | REAL, INTENT(in) :: sh(PLON,PLEV) ! specific humidity ( kg/kg ) |
---|
179 | REAL, INTENT(in) :: ql(PLON,PLEV) ! eau liquide |
---|
180 | REAL, INTENT(in) :: rh(PLON,PLEV) ! relative humidity |
---|
181 | REAL, INTENT(in) :: albs(PLON) ! surface albedo |
---|
182 | REAL, INTENT(out) :: source(PLON,nbtrac) |
---|
183 | |
---|
184 | !----------------------------------------------------------------------- |
---|
185 | ! ... Local variables |
---|
186 | !----------------------------------------------------------------------- |
---|
187 | INTEGER :: i, k, m, n, it,j |
---|
188 | INTEGER :: klev |
---|
189 | # if GRPCNT != 0 |
---|
190 | REAL :: group_ratios(PLON,PLEV,GRPCNT) |
---|
191 | # endif |
---|
192 | # if NCOL != 0 |
---|
193 | REAL :: col_dens(PLON,PLEV,NCOL) |
---|
194 | # else |
---|
195 | REAL :: col_dens(1) |
---|
196 | # endif |
---|
197 | # if RELCNT != 0 |
---|
198 | REAL :: rel_ratios(PLON,PLEV,RELCNT) |
---|
199 | # endif |
---|
200 | # if HETCNT != 0 |
---|
201 | REAL, ALLOCATABLE, SAVE :: het_rates(:,:,:) |
---|
202 | REAL, ALLOCATABLE, SAVE :: het_rates_st(:,:,:) |
---|
203 | REAL, ALLOCATABLE, SAVE :: het_rates_cv(:,:,:) |
---|
204 | !$OMP THREADPRIVATE(het_rates,het_rates_st,het_rates_cv) |
---|
205 | # else |
---|
206 | REAL :: het_rates(1) = 0.0 |
---|
207 | # endif |
---|
208 | REAL :: vmr(PLON,PLEV,PCNST) ! xported species ( vmr ) |
---|
209 | REAL :: reaction_rates(PLON,PLEV,RXNCNT) |
---|
210 | REAL, DIMENSION(PLON,PLEV) :: h2ovmr ! water vapor volume mixing ratio |
---|
211 | REAL, DIMENSION(PLON,PLEV) :: mbar ! mean wet atmospheric mass ( amu ) |
---|
212 | REAL, DIMENSION(PLON,PLEV) :: zmid ! midpoint geopotential in km |
---|
213 | REAL, DIMENSION(PCNST) :: vprodj ! volume production |
---|
214 | REAL, DIMENSION(PCNST) :: vlossj ! volume loss |
---|
215 | REAL, DIMENSION(PLON) :: zen_angle |
---|
216 | REAL, DIMENSION(PLON) :: loc_angle |
---|
217 | REAL :: sunon(PLON) |
---|
218 | REAL :: sunoff(PLON) |
---|
219 | REAL :: delt_inverse |
---|
220 | REAL :: zelev(PLON) |
---|
221 | REAL :: zalt(PLON,PLEV) |
---|
222 | LOGICAL :: polar_night(PLON) |
---|
223 | LOGICAL :: polar_day(PLON) |
---|
224 | LOGICAL :: zangtz(PLON) |
---|
225 | |
---|
226 | REAL :: tfld_lim(PLON,PLEV) |
---|
227 | REAL :: albs_lim(PLON) |
---|
228 | REAL :: tfld_glo(PLON_GLO,PLEV) |
---|
229 | REAL :: pmid_glo(PLON_GLO,PLEV) |
---|
230 | REAL, DIMENSION(PLON,PLEV,RXNCNT) :: reacflux |
---|
231 | LOGICAL,SAVE :: first = .TRUE. |
---|
232 | !$OMP THREADPRIVATE(first) |
---|
233 | |
---|
234 | |
---|
235 | !----------------------------------------------------------------------- |
---|
236 | ! ... Function interface |
---|
237 | !----------------------------------------------------------------------- |
---|
238 | REAL :: TSECND |
---|
239 | |
---|
240 | ! Exemple utilisation de writefied pour le debuggage |
---|
241 | ! DO i=1,PCNST |
---|
242 | ! CALL writefield_p('q_'//tracnam(i), mmr(:,:,i), PLEV) |
---|
243 | ! ENDDO |
---|
244 | |
---|
245 | reaction_rates(:,:,:) = 0.0 |
---|
246 | IF (first) THEN |
---|
247 | ALLOCATE(het_rates(PLON,PLEV,HETCNT)) |
---|
248 | ALLOCATE(het_rates_st(PLON,PLEV,HETCNT)) |
---|
249 | ALLOCATE(het_rates_cv(PLON,PLEV,HETCNT)) |
---|
250 | het_rates(:,:,:) = 0.0 |
---|
251 | het_rates_st(:,:,:) = 0.0 |
---|
252 | het_rates_cv(:,:,:) = hrates_cv(:,:,:) |
---|
253 | first=.FALSE. |
---|
254 | |
---|
255 | ENDIF |
---|
256 | |
---|
257 | delt_inverse = 1. / delt |
---|
258 | |
---|
259 | !----------------------------------------------------------------------- |
---|
260 | ! ... T and albedo adjustment for photorates |
---|
261 | !----------------------------------------------------------------------- |
---|
262 | tfld_lim(:,:) = MAX(tfld(:,:),226.) |
---|
263 | albs_lim(:) = MIN(albs(:),0.5) |
---|
264 | |
---|
265 | # if PHTCNT != 0 |
---|
266 | !----------------------------------------------------------------------- |
---|
267 | ! ... Calculate parameters for diurnal geometry |
---|
268 | !----------------------------------------------------------------------- |
---|
269 | CALL DIURNAL_GEOM( & |
---|
270 | lat, calday, polar_night, polar_day, & |
---|
271 | sunon, sunoff, loc_angle, zen_angle, & |
---|
272 | zangtz) |
---|
273 | |
---|
274 | # endif |
---|
275 | |
---|
276 | !----------------------------------------------------------------------- |
---|
277 | ! ... Initialize xform between mass and volume mixing ratios |
---|
278 | !----------------------------------------------------------------------- |
---|
279 | CALL INTI_MR_XFORM( sh, mbar ) |
---|
280 | |
---|
281 | !----------------------------------------------------------------------- |
---|
282 | ! ... Xform from mmr to vmr |
---|
283 | !----------------------------------------------------------------------- |
---|
284 | CALL MMR2VMR( vmr, mmr, mbar ) |
---|
285 | |
---|
286 | !----------------------------------------------------------------------- |
---|
287 | ! ... Xform water vapor from mmr to vmr and adjust in stratosphere |
---|
288 | !----------------------------------------------------------------------- |
---|
289 | #ifdef STRAT |
---|
290 | CALL ADJH2O( h2o, sh, ql, h2oc, mbar, vmr ) |
---|
291 | #else |
---|
292 | CALL ADJH2O( h2ovmr, sh,ql, h2oc, mbar, vmr ) |
---|
293 | #endif |
---|
294 | !----------------------------------------------------------------------- |
---|
295 | ! ... Xform geopotential height from m to km |
---|
296 | !----------------------------------------------------------------------- |
---|
297 | zelev(:) = 1.e-3 * phis(:PLON) / gravit |
---|
298 | zmid(:,:) = 1.e-3 * zma(:PLON,:) / gravit |
---|
299 | |
---|
300 | DO i = 1, PLON |
---|
301 | DO k = 1, PLEV |
---|
302 | zalt(i,k) = zmid(i,k) + zelev(i) |
---|
303 | END DO |
---|
304 | END DO |
---|
305 | |
---|
306 | # if NFS != 0 |
---|
307 | !----------------------------------------------------------------------- |
---|
308 | ! ... Set the "invariants" |
---|
309 | !----------------------------------------------------------------------- |
---|
310 | CALL SETINV(invariants, tfld, h2ovmr, pmid) |
---|
311 | # endif |
---|
312 | |
---|
313 | #ifndef DUSS |
---|
314 | # if NCOL != 0 |
---|
315 | !----------------------------------------------------------------------- |
---|
316 | ! ... Set the column densities at the upper boundary |
---|
317 | !----------------------------------------------------------------------- |
---|
318 | CALL SET_UB_COL(col_dens, vmr, invariants, pdel) |
---|
319 | |
---|
320 | # endif |
---|
321 | #endif |
---|
322 | |
---|
323 | !----------------------------------------------------------------------- |
---|
324 | ! ... Set rates for "tabular" and user specified reactions |
---|
325 | !----------------------------------------------------------------------- |
---|
326 | # if SETRXNCNT != 0 |
---|
327 | CALL SETRXT( reaction_rates, tfld ) |
---|
328 | # endif |
---|
329 | |
---|
330 | #ifndef DUSS |
---|
331 | # if USRRXNCNT != 0 |
---|
332 | CALL USRRXT( & |
---|
333 | reaction_rates(:,:,:) , & |
---|
334 | tfld , & |
---|
335 | invariants , & |
---|
336 | rh , & |
---|
337 | ps , & |
---|
338 | pmid , & |
---|
339 | mmr , & |
---|
340 | invariants(1,1,INDEXM) , & |
---|
341 | lat , & |
---|
342 | zen_angle , & |
---|
343 | vmr ) |
---|
344 | # endif |
---|
345 | #endif |
---|
346 | |
---|
347 | # if GASCNT != 0 |
---|
348 | # ifdef RADJFLAG |
---|
349 | ! changement d'unite des taux de reactions des reactions a 2 reactants |
---|
350 | ! ils sont calcules dans usrrxt en s-1/(mol.cm-3) et le code dans |
---|
351 | ! exp_sol et imp_sol attends des s-1 |
---|
352 | ! ce n'est pas le cas des reactions a 1 reactant dont les taux |
---|
353 | ! sont directement calcule en s-1 par usrrxt (so4aerrxt etc.) |
---|
354 | CALL ADJRXT(reaction_rates, invariants, invariants(1,1,INDEXM)) |
---|
355 | # endif |
---|
356 | # endif |
---|
357 | |
---|
358 | # if PHTCNT != 0 |
---|
359 | !----------------------------------------------------------------------- |
---|
360 | ! ... Compute the photolysis rates |
---|
361 | !----------------------------------------------------------------------- |
---|
362 | # if NCOL != 0 |
---|
363 | !----------------------------------------------------------------------- |
---|
364 | ! ... Set the column densities |
---|
365 | !----------------------------------------------------------------------- |
---|
366 | CALL SETCOL(col_dens, vmr, invariants, pdel) |
---|
367 | # endif |
---|
368 | |
---|
369 | !----------------------------------------------------------------------- |
---|
370 | ! ... Calculate the photodissociation rates |
---|
371 | !----------------------------------------------------------------------- |
---|
372 | |
---|
373 | CALL PHOTO(reaction_rates, pmid, pdel, & |
---|
374 | tfld, zalt, col_dens, zen_angle, & |
---|
375 | albs, cwat, cldfr, sunon, sunoff, & |
---|
376 | polar_night, polar_day) |
---|
377 | |
---|
378 | !----------------------------------------------------------------------- |
---|
379 | ! ... Adjust the photodissociation rates |
---|
380 | !----------------------------------------------------------------------- |
---|
381 | CALL PHTADJ(reaction_rates, invariants, invariants(1,1,INDEXM)) |
---|
382 | |
---|
383 | # endif |
---|
384 | |
---|
385 | # if HETCNT != 0 |
---|
386 | !----------------------------------------------------------------------- |
---|
387 | ! ... Compute the heterogeneous rates at time = t(n+1) |
---|
388 | !----------------------------------------------------------------------- |
---|
389 | |
---|
390 | ! ... For stratiform precipitation |
---|
391 | |
---|
392 | CALL SETHET( & |
---|
393 | het_rates_st ,& |
---|
394 | pmid ,& |
---|
395 | pdel ,& |
---|
396 | lat ,& |
---|
397 | zmid ,& |
---|
398 | tfld ,& |
---|
399 | delt ,& |
---|
400 | invariants(1,1,INDEXM) ,& |
---|
401 | flxrst ,& |
---|
402 | flxsst ,& |
---|
403 | flxupd ,& |
---|
404 | cldtop ,& |
---|
405 | cldbot ,& |
---|
406 | cldfr_st ,& |
---|
407 | 1 ,& |
---|
408 | vmr ) |
---|
409 | |
---|
410 | ! ... For convective precipitation |
---|
411 | |
---|
412 | CALL SETHET( & |
---|
413 | het_rates_cv ,& |
---|
414 | pmid ,& |
---|
415 | pdel ,& |
---|
416 | lat ,& |
---|
417 | zmid ,& |
---|
418 | tfld ,& |
---|
419 | delt ,& |
---|
420 | invariants(1,1,INDEXM) ,& |
---|
421 | flxrcv ,& |
---|
422 | flxscv ,& |
---|
423 | flxupd ,& |
---|
424 | cldtop ,& |
---|
425 | cldbot ,& |
---|
426 | cldfr_cv ,& |
---|
427 | 2 ,& |
---|
428 | vmr ) |
---|
429 | |
---|
430 | ! ... Total removal rate |
---|
431 | |
---|
432 | het_rates = het_rates_st + het_rates_cv |
---|
433 | hrates = het_rates |
---|
434 | hrates_cv = het_rates_cv |
---|
435 | |
---|
436 | # endif |
---|
437 | # if EXTCNT != 0 |
---|
438 | !----------------------------------------------------------------------- |
---|
439 | ! ... Compute the extraneous forcings |
---|
440 | !----------------------------------------------------------------------- |
---|
441 | |
---|
442 | !Initialize for this timeslice |
---|
443 | DO m = 1,EXTCNT |
---|
444 | extfrc(:,:,m) = 0. |
---|
445 | END DO |
---|
446 | |
---|
447 | #ifndef GES |
---|
448 | CALL SETEXT_LGHT( & |
---|
449 | calday, extfrc, lat, zmid, & |
---|
450 | area, cldtop, vmr, tfld, & |
---|
451 | col_dens, zen_angle, delt, & |
---|
452 | invariants) |
---|
453 | |
---|
454 | #ifdef NMHC |
---|
455 | CALL SETEXT_O3L( & |
---|
456 | calday, extfrc, lat, & |
---|
457 | zmid, area, cldtop, & |
---|
458 | vmr, tfld, col_dens, & |
---|
459 | zen_angle, delt, & |
---|
460 | invariants) |
---|
461 | |
---|
462 | CALL SETEXT_BBG( & |
---|
463 | calday, extfrc, lat, & |
---|
464 | zmid, area, cldtop, & |
---|
465 | vmr, tfld, col_dens, & |
---|
466 | zen_angle, delt, pdel, & |
---|
467 | invariants) |
---|
468 | |
---|
469 | IF (flag_plane .GE. 1) THEN |
---|
470 | CALL SETEXT_AIR( & |
---|
471 | calday, extfrc, lat, & |
---|
472 | zmid, area, cldtop, & |
---|
473 | vmr, tfld, col_dens, & |
---|
474 | zen_angle, delt, & |
---|
475 | invariants) |
---|
476 | ENDIF |
---|
477 | #endif |
---|
478 | |
---|
479 | #else |
---|
480 | DO m = 1,EXTCNT |
---|
481 | extfrc(:,:,m) = 0. |
---|
482 | END DO |
---|
483 | #endif |
---|
484 | |
---|
485 | # endif |
---|
486 | |
---|
487 | # if GRPCNT != 0 |
---|
488 | |
---|
489 | !----------------------------------------------------------------------- |
---|
490 | ! ... Set the group ratios |
---|
491 | !----------------------------------------------------------------------- |
---|
492 | CALL SET_GRP_RATIOS( & |
---|
493 | group_ratios, reaction_rates(:,:,:), & |
---|
494 | reaction_rates, vmr, mmr, nas, mbar, invariants) |
---|
495 | |
---|
496 | |
---|
497 | !----------------------------------------------------------------------- |
---|
498 | ! ... Vertical mass flux |
---|
499 | !----------------------------------------------------------------------- |
---|
500 | DO i = 1, PLON |
---|
501 | !klev = itrop(i) |
---|
502 | klev = 11 |
---|
503 | o3_st_flx(:)=flxmass_w(:,klev)*mmr(:,klev,id_OX)/uma/nadv_mass(id_O3) |
---|
504 | END DO |
---|
505 | |
---|
506 | # endif |
---|
507 | |
---|
508 | # if RELCNT != 0 |
---|
509 | !----------------------------------------------------------------------- |
---|
510 | ! ... Set the relationship ratios |
---|
511 | !----------------------------------------------------------------------- |
---|
512 | # if GRPCNT != 0 |
---|
513 | CALL SET_REL_RATIOS(& |
---|
514 | rel_ratios , PLON, PLEV, & |
---|
515 | RELCNT, group_ratios, & |
---|
516 | GRPCNT, reaction_rates(1,1,PHTCNTP1), & |
---|
517 | GASCNT, reaction_rates, PHTCNT, vmr, & |
---|
518 | PCNST, invariants, NFS) |
---|
519 | #else |
---|
520 | CALL SET_REL_RATIOS(& |
---|
521 | rel_ratios , PLON, PLEV, & |
---|
522 | RELCNT, reaction_rates(1,1,PHTCNTP1), & |
---|
523 | GASCNT, reaction_rates, PHTCNT, vmr, & |
---|
524 | PCNST, invariants, NFS) |
---|
525 | #endif |
---|
526 | #endif |
---|
527 | |
---|
528 | |
---|
529 | #if RELCNT != 0 || GRPCNT != 0 |
---|
530 | !----------------------------------------------------------------------- |
---|
531 | ! ... Modify the reaction rate of any reaction |
---|
532 | ! with group member or proportional reactant(s) |
---|
533 | !----------------------------------------------------------------------- |
---|
534 | #if RELCNT != 0 |
---|
535 | CALL RXT_MOD(reaction_rates, het_rates, rel_ratios, group_ratios) |
---|
536 | #else |
---|
537 | CALL RXT_MOD(reaction_rates, het_rates, group_ratios) |
---|
538 | #endif |
---|
539 | #endif |
---|
540 | |
---|
541 | !======================================================================= |
---|
542 | ! ... Call the class solution algorithms |
---|
543 | !======================================================================= |
---|
544 | #if CLSCNT1 != 0 |
---|
545 | !----------------------------------------------------------------------- |
---|
546 | ! ... Solve for "explicit" species |
---|
547 | !----------------------------------------------------------------------- |
---|
548 | CALL EXP_SOL(vmr, reaction_rates, het_rates, extfrc, & |
---|
549 | #ifdef NMHC |
---|
550 | co_prod, co_loss, ch4_loss, n2o_loss, invariants(:,:,INDEXM),& |
---|
551 | #endif |
---|
552 | nstep, delt) |
---|
553 | #endif |
---|
554 | |
---|
555 | #ifndef DUSS |
---|
556 | # if CLSCNT4 != 0 |
---|
557 | !----------------------------------------------------------------------- |
---|
558 | ! ... Solve for Implicit species |
---|
559 | !----------------------------------------------------------------------- |
---|
560 | |
---|
561 | CALL IMP_SOL( & |
---|
562 | vmr, & |
---|
563 | reaction_rates, & |
---|
564 | het_rates, & |
---|
565 | extfrc, & |
---|
566 | o3_prod, & |
---|
567 | o3_loss, & |
---|
568 | #ifdef NMHC |
---|
569 | co2_basprod, & |
---|
570 | co2_nmhcprod, & |
---|
571 | co2_radicalprod, & |
---|
572 | hno3_prod, & |
---|
573 | hno3_loss, & |
---|
574 | #endif |
---|
575 | #ifdef AER |
---|
576 | so2_p_h2soh, & |
---|
577 | so2_p_dmsoh, & |
---|
578 | so2_p_dmsno3, & |
---|
579 | so2_p_dmsooh, & |
---|
580 | asmsam_p_dmsooh, & |
---|
581 | dmso_p_dmsoh, & |
---|
582 | asso4m_p_so2oh, & |
---|
583 | wet3d_so2, & |
---|
584 | wet3d_dms, & |
---|
585 | wet3d_dmso, & |
---|
586 | #ifdef NMHC |
---|
587 | pom_p_g, & |
---|
588 | hno3_p_g, & |
---|
589 | hno3_p_a, & |
---|
590 | hno3_l_g, & |
---|
591 | hno3_l_a, & |
---|
592 | hono_p_g, & |
---|
593 | hono_p_a, & |
---|
594 | hono_l_g, & |
---|
595 | hono_l_a, & |
---|
596 | nh3_l_g, & |
---|
597 | nh3_l_a, & |
---|
598 | csno3_p_a1, & |
---|
599 | csno3_p_a2, & |
---|
600 | cino3_p_a, & |
---|
601 | wet3d_h2o2, & |
---|
602 | wet3d_hono, & |
---|
603 | #endif |
---|
604 | wet3d_hno3, & |
---|
605 | wet3d_nh3, & |
---|
606 | wet3d_noy, & |
---|
607 | wet3d_asso4m, & |
---|
608 | wet3d_asnh4m, & |
---|
609 | wet3d_asno3m, & |
---|
610 | wet3d_csno3m, & |
---|
611 | wet3d_cino3m, & |
---|
612 | #endif |
---|
613 | pdel, & |
---|
614 | #ifdef AERONLY |
---|
615 | invariants, & |
---|
616 | #endif |
---|
617 | invariants(:,:,INDEXM), & |
---|
618 | nstep, & |
---|
619 | delt_inverse, & |
---|
620 | lat, wetloss ) |
---|
621 | # endif |
---|
622 | #endif |
---|
623 | |
---|
624 | #ifdef STRAT |
---|
625 | !----------------------------------------------------------------------- |
---|
626 | ! ... Sedimentation in the stratosphere |
---|
627 | !----------------------------------------------------------------------- |
---|
628 | CALL STRSEDCALC( invariants(:,:,INDEXM), pdel, vmr, delt,h2o, h2oc ) |
---|
629 | #endif |
---|
630 | |
---|
631 | DO i = 1, PLON |
---|
632 | #ifdef AER |
---|
633 | #if !defined(AERONLY) |
---|
634 | wflux(i,id_ASSO4M) = SUM(wet3d_asso4m(i,:)) |
---|
635 | wflux(i,id_ASNH4M) = SUM(wet3d_asnh4m(i,:)) |
---|
636 | wflux(i,id_ASNO3M) = SUM(wet3d_asno3m(i,:)) |
---|
637 | wflux(i,id_CSNO3M) = SUM(wet3d_csno3m(i,:)) |
---|
638 | wflux(i,id_CINO3M) = SUM(wet3d_cino3m(i,:)) |
---|
639 | #endif |
---|
640 | wfluxso2(i) = SUM(wet3d_so2(i,:)) |
---|
641 | wfluxdms(i) = SUM(wet3d_dms(i,:)) |
---|
642 | wfluxdmso(i) = SUM(wet3d_dmso(i,:)) |
---|
643 | wfluxnoy(i) = SUM(wet3d_noy(i,:)) |
---|
644 | wfluxnh3(i) = SUM(wet3d_nh3(i,:)) |
---|
645 | wfluxhno3(i) = SUM(wet3d_hno3(i,:)) |
---|
646 | #endif |
---|
647 | ENDDO |
---|
648 | |
---|
649 | #if !defined(DUSS) && defined(AER) |
---|
650 | !DH Nitrates formation |
---|
651 | CALL AERTHERM( & |
---|
652 | delt , & |
---|
653 | tfld , & |
---|
654 | rh , & |
---|
655 | pmid , & |
---|
656 | invariants(:,:,INDEXM) , & |
---|
657 | asno3m_p_nh3hno3 , & |
---|
658 | asnh4m_p_nh3hno3 , & |
---|
659 | hno3_p_nh3hno3 , & |
---|
660 | nh3_p_nh3hno3 , & |
---|
661 | mmr , & |
---|
662 | vmr ) |
---|
663 | #endif |
---|
664 | !----------------------------------------------------------------------- |
---|
665 | ! ... Check for negative values and reset to zero |
---|
666 | !----------------------------------------------------------------------- |
---|
667 | CALL NEGTRC( 'After chemistry ', vmr ) |
---|
668 | |
---|
669 | |
---|
670 | IF (calcul_flux) THEN |
---|
671 | !----------------------------------------------------------------------- |
---|
672 | ! ... Compute the flux through each reaction |
---|
673 | !----------------------------------------------------------------------- |
---|
674 | CALL reac_flx(vmr, reaction_rates, invariants, reacflux) |
---|
675 | |
---|
676 | CALL xios_inca_change_context("inca") |
---|
677 | DO n=1, RXNCNT |
---|
678 | CALL xios_inca_send_field("flux_"//trim(reacname(n)), reacflux(:,:,n)) |
---|
679 | ENDDO |
---|
680 | CALL xios_inca_change_context("LMDZ") |
---|
681 | ENDIF |
---|
682 | |
---|
683 | !----------------------------------------------------------------------- |
---|
684 | ! ... Xform from vmr to mmr |
---|
685 | !----------------------------------------------------------------------- |
---|
686 | CALL VMR2MMR( vmr, mmr, & |
---|
687 | # if GRPCNT != 0 |
---|
688 | nas, & |
---|
689 | group_ratios, & |
---|
690 | # endif |
---|
691 | mbar ) |
---|
692 | |
---|
693 | #ifndef DUSS |
---|
694 | #ifdef AER |
---|
695 | CALL GASTOPARTICLE(delt,mmr,mbar,invariants(:,:,INDEXM)) |
---|
696 | |
---|
697 | #endif |
---|
698 | #endif |
---|
699 | |
---|
700 | #if defined(NMHC) && defined(AER) |
---|
701 | CALL soa_main(mmr, & |
---|
702 | mbar, & |
---|
703 | pmid, & |
---|
704 | tfld, & |
---|
705 | ASAPp1a_p, & |
---|
706 | ASAPp2a_p, & |
---|
707 | ASARp1a_p, & |
---|
708 | ASARp2a_p, & |
---|
709 | delt) |
---|
710 | #endif |
---|
711 | |
---|
712 | !DH Update the drydep flux with latest surface mixing ratio |
---|
713 | CALL DRYDEP (pmid, & |
---|
714 | tfld, & |
---|
715 | mmr) |
---|
716 | |
---|
717 | ! calcul des sources |
---|
718 | DO it=1,nbtrac |
---|
719 | DO k = 1, PLON |
---|
720 | source(k,it) = eflux(k,it)-dflux(k,it) |
---|
721 | ENDDO |
---|
722 | ENDDO |
---|
723 | |
---|
724 | |
---|
725 | flux_source(:,:) = source(:,:) |
---|
726 | |
---|
727 | #if defined(NMHC) && defined(AER) |
---|
728 | CALL xios_inca_change_context("inca") ! |
---|
729 | CALL xios_inca_send_field("O3_surf",vmr(:,1,id_O3)) ! Define surface concentration of O3 --- |
---|
730 | CALL xios_inca_send_field("SO2_surf",vmr(:,1,id_SO2)) ! Define surface concentration of SO2 --- |
---|
731 | CALL xios_inca_send_field("NO2_surf",vmr(:,1,id_NO2)) ! Define surface concentration of NO2 --- |
---|
732 | CALL xios_inca_send_field("NH3_surf",vmr(:,1,id_NH3)) ! Define surface concentration of NH3 --- |
---|
733 | CALL xios_inca_send_field("CH2O_surf",vmr(:,1,id_CH2O)) ! Define surface concentration of CH2O --- |
---|
734 | CALL xios_inca_send_field("OH_surf",vmr(:,1,id_OH)) ! Define surface concentration of OH --- |
---|
735 | CALL xios_inca_change_context("LMDZ") ! |
---|
736 | #endif |
---|
737 | |
---|
738 | |
---|
739 | END SUBROUTINE CHEMMAIN |
---|
740 | |
---|
741 | SUBROUTINE endrun |
---|
742 | !---------------------------------------------------------------------- |
---|
743 | ! ... Abort the model |
---|
744 | !---------------------------------------------------------------------- |
---|
745 | |
---|
746 | IMPLICIT NONE |
---|
747 | |
---|
748 | EXTERNAL abort |
---|
749 | |
---|
750 | CALL abort |
---|
751 | |
---|
752 | END SUBROUTINE endrun |
---|
753 | |
---|