1 |
guez |
7 |
module regr_pr_comb_coefoz_m |
2 |
guez |
3 |
|
3 |
|
|
use dimens_m, only: llm |
4 |
|
|
use dimphy, only: klon |
5 |
|
|
|
6 |
|
|
implicit none |
7 |
|
|
|
8 |
guez |
10 |
! The five module variables declared here are on the "physics" grid. |
9 |
|
|
! The value of each variable for index "(i, k)" is at longitude |
10 |
|
|
! "rlon(i)", latitude "rlat(i)" and middle of layer "k". |
11 |
|
|
|
12 |
guez |
7 |
real, save:: c_Mob(klon, llm) |
13 |
guez |
3 |
! (sum of Mobidic terms in the net mass production rate of ozone |
14 |
|
|
! by chemistry, per unit mass of air, in s-1) |
15 |
|
|
|
16 |
guez |
7 |
real, save:: a2(klon, llm) |
17 |
guez |
3 |
! (derivative of mass production rate of ozone per unit mass of |
18 |
|
|
! air with respect to ozone mass fraction, in s-1) |
19 |
|
|
|
20 |
guez |
7 |
real, save:: a4_mass(klon, llm) |
21 |
guez |
3 |
! (derivative of mass production rate of ozone per unit mass of |
22 |
|
|
! air with respect to temperature, in s-1 K-1) |
23 |
|
|
|
24 |
guez |
7 |
real, save:: a6_mass(klon, llm) |
25 |
guez |
3 |
! (derivative of mass production rate of ozone per unit mass of |
26 |
|
|
! air with respect to mass column-density of ozone above, in m2 s-1 kg-1) |
27 |
|
|
|
28 |
guez |
7 |
real, save:: r_het_interm(klon, llm) |
29 |
guez |
3 |
! (net mass production rate by heterogeneous chemistry, per unit |
30 |
|
|
! mass of ozone, corrected for chlorine content and latitude, but |
31 |
|
|
! not for temperature and sun direction, in s-1) |
32 |
|
|
|
33 |
|
|
private klon, llm |
34 |
|
|
|
35 |
|
|
contains |
36 |
|
|
|
37 |
guez |
7 |
subroutine regr_pr_comb_coefoz(julien) |
38 |
guez |
3 |
|
39 |
guez |
7 |
! "regr_pr_comb_coefoz" stands for "regrid pressure combine |
40 |
|
|
! coefficients ozone". |
41 |
guez |
3 |
|
42 |
guez |
7 |
! This subroutine : |
43 |
guez |
10 |
! -- reads from a file all eight coefficients for ozone chemistry, |
44 |
guez |
7 |
! at the current day ; |
45 |
guez |
10 |
! -- regrids the coefficients in pressure to the LMDZ vertical grid ; |
46 |
|
|
! -- packs the coefficients to the "physics" horizontal grid ; |
47 |
|
|
! -- combines the eight coefficients to define the five module variables. |
48 |
guez |
7 |
|
49 |
guez |
18 |
use netcdf95, only: nf95_open, nf95_close |
50 |
guez |
7 |
use netcdf, only: nf90_nowrite |
51 |
|
|
use regr_pr_coefoz, only: regr_pr_av_coefoz, regr_pr_int_coefoz |
52 |
guez |
3 |
use phyetat0_m, only: rlat |
53 |
|
|
|
54 |
guez |
7 |
integer, intent(in):: julien ! jour julien, 1 <= julien <= 360 |
55 |
|
|
|
56 |
guez |
3 |
! Variables local to the procedure: |
57 |
guez |
18 |
|
58 |
guez |
3 |
integer ncid ! for NetCDF |
59 |
|
|
|
60 |
guez |
7 |
real coefoz(klon, llm) |
61 |
|
|
! (temporary storage for an ozone coefficient) |
62 |
|
|
! (On the "physics" grid. |
63 |
|
|
! "coefoz(i, k)" is at longitude "rlon(i)", latitude "rlat(i)", |
64 |
|
|
! middle of layer "k".) |
65 |
|
|
|
66 |
|
|
real a6(klon, llm) |
67 |
|
|
! (derivative of "P_net_Mob" with respect to column-density of ozone |
68 |
guez |
3 |
! above, in cm2 s-1) |
69 |
|
|
! (On the "physics" grid. |
70 |
guez |
7 |
! "a6(i, k)" is at longitude "rlon(i)", latitude "rlat(i)", |
71 |
guez |
3 |
! middle of layer "k".) |
72 |
|
|
|
73 |
|
|
real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg |
74 |
|
|
|
75 |
|
|
real, parameter:: Clx = 3.8e-9 |
76 |
|
|
! (total chlorine content in the upper stratosphere) |
77 |
|
|
|
78 |
guez |
7 |
integer k |
79 |
guez |
3 |
|
80 |
|
|
!------------------------------------ |
81 |
|
|
|
82 |
guez |
10 |
print *, "Call sequence information: regr_pr_comb_coefoz" |
83 |
guez |
3 |
|
84 |
|
|
call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) |
85 |
|
|
|
86 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "a2", julien, a2) |
87 |
guez |
3 |
|
88 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "a4", julien, a4_mass) |
89 |
guez |
7 |
a4_mass = a4_mass * 48. / 29. |
90 |
|
|
|
91 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "a6", julien, a6) |
92 |
guez |
7 |
|
93 |
guez |
3 |
! Compute "a6_mass" avoiding underflow, do not divide by 1e4 |
94 |
|
|
! before dividing by molecular mass: |
95 |
|
|
a6_mass = a6 / (1e4 * 29. * amu) |
96 |
|
|
! (factor 1e4: conversion from cm2 to m2) |
97 |
|
|
|
98 |
guez |
7 |
! Combine coefficients to get "c_Mob": |
99 |
|
|
! (We use as few local variables as possible, in order to spare |
100 |
|
|
! main memory.) |
101 |
guez |
3 |
|
102 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "P_net_Mob", julien, c_Mob) |
103 |
guez |
7 |
|
104 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "r_Mob", julien, coefoz) |
105 |
guez |
9 |
c_mob = c_mob - a2 * coefoz |
106 |
guez |
7 |
|
107 |
guez |
17 |
call regr_pr_int_coefoz(ncid, "Sigma_Mob", julien, top_value=0., v3=coefoz) |
108 |
guez |
7 |
c_mob = (c_mob - a6 * coefoz) * 48. / 29. |
109 |
|
|
|
110 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "temp_Mob", julien, coefoz) |
111 |
guez |
7 |
c_mob = c_mob - a4_mass * coefoz |
112 |
|
|
|
113 |
guez |
17 |
call regr_pr_av_coefoz(ncid, "R_Het", julien, r_het_interm) |
114 |
guez |
3 |
! Heterogeneous chemistry is only at high latitudes: |
115 |
guez |
7 |
forall (k = 1: llm) |
116 |
|
|
where (abs(rlat) <= 45.) r_het_interm(:, k) = 0. |
117 |
guez |
3 |
end forall |
118 |
guez |
11 |
r_het_interm = r_het_interm * (Clx / 3.8e-9)**2 |
119 |
guez |
3 |
|
120 |
|
|
call nf95_close(ncid) |
121 |
|
|
|
122 |
guez |
7 |
end subroutine regr_pr_comb_coefoz |
123 |
guez |
3 |
|
124 |
guez |
7 |
end module regr_pr_comb_coefoz_m |