1 |
module o3_chem_m |
2 |
|
3 |
! This module is clean: no C preprocessor directive, no include line. |
4 |
|
5 |
IMPLICIT none |
6 |
|
7 |
private o3_prod |
8 |
|
9 |
contains |
10 |
|
11 |
subroutine o3_chem(julien, gmtime, t_seri, zmasse, pdtphys, q) |
12 |
|
13 |
! This procedure evolves the ozone mass fraction through a time |
14 |
! step taking only chemistry into account. |
15 |
|
16 |
! All the 2-dimensional arrays are on the "physics" grid. |
17 |
! Their shape is "(/klon, llm/)". |
18 |
! Index "(i, :)" is for longitude "rlon(i)", latitude "rlat(i)". |
19 |
|
20 |
use numer_rec, only: assert, pi |
21 |
use dimphy, only: klon |
22 |
use dimens_m, only: llm |
23 |
use regr_pr_comb_coefoz_m, only: c_Mob, a4_mass, a2, r_het_interm |
24 |
use orbite_m, only: orbite, zenang |
25 |
|
26 |
integer, intent(in):: julien ! jour julien, 1 <= julien <= 360 |
27 |
real, intent(in):: gmtime ! heure de la journée en fraction de jour |
28 |
real, intent(in):: t_seri(:, :) ! (klon, llm) temperature, in K |
29 |
|
30 |
real, intent(in):: zmasse(:, :) ! (klon, llm) |
31 |
! (column-density of mass of air in a cell, in kg m-2) |
32 |
! "zmasse(:, k)" is for layer "k".) |
33 |
|
34 |
real, intent(in):: pdtphys ! time step for physics, in s |
35 |
|
36 |
real, intent(inout):: q(:, :) ! (klon, llm) mass fraction of ozone |
37 |
! "q(:, k)" is at middle of layer "k".) |
38 |
|
39 |
! Variables local to the procedure: |
40 |
integer k |
41 |
|
42 |
real c(klon, llm) |
43 |
! (constant term during a time step in the net mass production |
44 |
! rate of ozone by chemistry, per unit mass of air, in s-1) |
45 |
! "c(:, k)" is at middle of layer "k".) |
46 |
|
47 |
real b(klon, llm) |
48 |
! (coefficient of "q" in the net mass production |
49 |
! rate of ozone by chemistry, per unit mass of air, in s-1) |
50 |
! "b(:, k)" is at middle of layer "k".) |
51 |
|
52 |
real dq_o3_chem(klon, llm) |
53 |
! (variation of ozone mass fraction due to chemistry during a time step) |
54 |
! "dq_o3_chem(:, k)" is at middle of layer "k".) |
55 |
|
56 |
real earth_long |
57 |
! (longitude vraie de la Terre dans son orbite solaire, par |
58 |
! rapport au point vernal (21 mars), en degrés) |
59 |
|
60 |
real pmu0(klon) ! mean of cosine of solar zenith angle during "pdtphys" |
61 |
|
62 |
!------------------------------------------------------------- |
63 |
|
64 |
call assert(klon == (/size(q, 1), size(t_seri, 1), size(zmasse, 1)/), & |
65 |
"o3_chem klon") |
66 |
call assert(llm == (/size(q, 2), size(t_seri, 2), size(zmasse, 2)/), & |
67 |
"o3_chem llm") |
68 |
|
69 |
c = c_Mob + a4_mass * t_seri |
70 |
|
71 |
! Compute coefficient "b": |
72 |
|
73 |
! Heterogeneous chemistry is only at low temperature: |
74 |
where (t_seri < 195.) |
75 |
b = r_het_interm |
76 |
elsewhere |
77 |
b = 0. |
78 |
end where |
79 |
|
80 |
! Heterogeneous chemistry is only during daytime: |
81 |
call orbite(real(julien), earth_long) |
82 |
call zenang(earth_long, gmtime, pdtphys, pmu0) |
83 |
forall (k = 1: llm) |
84 |
where (pmu0 <= cos(87. / 180. * pi)) b(:, k) = 0. |
85 |
end forall |
86 |
|
87 |
b = b + a2 |
88 |
|
89 |
! Midpoint method: |
90 |
|
91 |
! Trial step to the midpoint: |
92 |
dq_o3_chem = o3_prod(q, zmasse, c, b) * pdtphys / 2 |
93 |
! "Real" step across the whole interval: |
94 |
dq_o3_chem = o3_prod(q + dq_o3_chem, zmasse, c, b) * pdtphys |
95 |
q = q + dq_o3_chem |
96 |
|
97 |
! Confine the mass fraction: |
98 |
q = min(max(q, 0.), .01) |
99 |
|
100 |
end subroutine o3_chem |
101 |
|
102 |
!************************************************* |
103 |
|
104 |
function o3_prod(q, zmasse, c, b) |
105 |
|
106 |
! This function computes the production rate of ozone by chemistry. |
107 |
|
108 |
! All the 2-dimensional arrays are on the "physics" grid. |
109 |
! Their shape is "(/klon, llm/)". |
110 |
! Index "(i, :)" is for longitude "rlon(i)", latitude "rlat(i)". |
111 |
|
112 |
use regr_pr_comb_coefoz_m, only: a6_mass |
113 |
use numer_rec, only: assert |
114 |
use dimens_m, only: llm |
115 |
use dimphy, only: klon |
116 |
|
117 |
real, intent(in):: q(:, :) ! mass fraction of ozone |
118 |
! "q(:, k)" is at middle of layer "k".) |
119 |
|
120 |
real, intent(in):: zmasse(:, :) |
121 |
! (column-density of mass of air in a layer, in kg m-2) |
122 |
! ("zmasse(:, k)" is for layer "k".) |
123 |
|
124 |
real, intent(in):: c(:, :) |
125 |
! (constant term during a time step in the net mass production |
126 |
! rate of ozone by chemistry, per unit mass of air, in s-1) |
127 |
! "c(:, k)" is at middle of layer "k".) |
128 |
|
129 |
real, intent(in):: b(:, :) |
130 |
! (coefficient of "q" in the net mass production rate of ozone by |
131 |
! chemistry, per unit mass of air, in s-1) |
132 |
! ("b(:, k)" is at middle of layer "k".) |
133 |
|
134 |
real o3_prod(klon, llm) |
135 |
! (net mass production rate of ozone by chemistry, per unit mass |
136 |
! of air, in s-1) |
137 |
! ("o3_prod(:, k)" is at middle of layer "k".) |
138 |
|
139 |
! Variables local to the procedure: |
140 |
|
141 |
real sigma_mass(klon, llm) |
142 |
! (mass column-density of ozone above point, in kg m-2) |
143 |
! ("sigma_mass(:, k)" is at middle of layer "k".) |
144 |
|
145 |
integer k |
146 |
|
147 |
!------------------------------------------------------------------- |
148 |
|
149 |
call assert(klon == (/size(q, 1), size(zmasse, 1), size(c, 1), & |
150 |
size(b, 1)/), "o3_prod 1") |
151 |
call assert(llm == (/size(q, 2), size(zmasse, 2), size(c, 2), & |
152 |
size(b, 2)/), "o3_prod 2") |
153 |
|
154 |
! Compute the column-density above the base of layer |
155 |
! "k", and, as a first approximation, take it as column-density |
156 |
! above the middle of layer "k": |
157 |
sigma_mass(:, llm) = zmasse(:, llm) * q(:, llm) ! top layer |
158 |
do k = llm - 1, 1, -1 |
159 |
sigma_mass(:, k) = sigma_mass(:, k+1) + zmasse(:, k) * q(:, k) |
160 |
end do |
161 |
|
162 |
o3_prod = c + b * q + a6_mass * sigma_mass |
163 |
|
164 |
end function o3_prod |
165 |
|
166 |
end module o3_chem_m |