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