1 |
guez |
3 |
module read_coefoz_m |
2 |
|
|
|
3 |
|
|
! This module is clean: no C preprocessor directive, no include line. |
4 |
|
|
|
5 |
|
|
use dimens_m, only: llm |
6 |
|
|
use dimphy, only: klon |
7 |
|
|
|
8 |
|
|
implicit none |
9 |
|
|
|
10 |
|
|
real, save:: c_Mob(klon, llm, 12) |
11 |
|
|
! (sum of Mobidic terms in the net mass production rate of ozone |
12 |
|
|
! by chemistry, per unit mass of air, in s-1) |
13 |
|
|
! (On the "physics" grid. |
14 |
|
|
! Third dimension is the number of the month in the year. |
15 |
|
|
! "c_Mob(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
16 |
|
|
! middle of layer "k".) |
17 |
|
|
|
18 |
|
|
real, save:: a2(klon, llm, 12) |
19 |
|
|
! (derivative of mass production rate of ozone per unit mass of |
20 |
|
|
! air with respect to ozone mass fraction, in s-1) |
21 |
|
|
! (On the "physics" grid. |
22 |
|
|
! Third dimension is the number of the month in the year. |
23 |
|
|
! "a2(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
24 |
|
|
! middle of layer "k".) |
25 |
|
|
|
26 |
|
|
real, save:: a4_mass(klon, llm, 12) |
27 |
|
|
! (derivative of mass production rate of ozone per unit mass of |
28 |
|
|
! air with respect to temperature, in s-1 K-1) |
29 |
|
|
! (On the "physics" grid. |
30 |
|
|
! Third dimension is the number of the month in the year. |
31 |
|
|
! "a4_mass(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
32 |
|
|
! middle of layer "k".) |
33 |
|
|
|
34 |
|
|
real, save:: a6_mass(klon, llm, 12) |
35 |
|
|
! (derivative of mass production rate of ozone per unit mass of |
36 |
|
|
! air with respect to mass column-density of ozone above, in m2 s-1 kg-1) |
37 |
|
|
! (On the "physics" grid. |
38 |
|
|
! Third dimension is the number of the month in the year. |
39 |
|
|
! "a6_mass(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
40 |
|
|
! middle of layer "k".) |
41 |
|
|
|
42 |
|
|
real, save:: r_het_interm(klon, llm, 12) |
43 |
|
|
! (net mass production rate by heterogeneous chemistry, per unit |
44 |
|
|
! mass of ozone, corrected for chlorine content and latitude, but |
45 |
|
|
! not for temperature and sun direction, in s-1) |
46 |
|
|
! (On the "physics" grid. |
47 |
|
|
! Third dimension is the number of the month in the year. |
48 |
|
|
! "r_het_interm(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
49 |
|
|
! middle of layer "k".) |
50 |
|
|
|
51 |
|
|
private klon, llm |
52 |
|
|
|
53 |
|
|
contains |
54 |
|
|
|
55 |
|
|
subroutine read_coefoz |
56 |
|
|
|
57 |
|
|
! This subroutine reads from a file all eight parameters for ozone |
58 |
|
|
! chemistry, at all months. |
59 |
|
|
! The parameters are packed to the "physics" grid. |
60 |
|
|
! Finally, the eight parameters are combined to define the five |
61 |
|
|
! module variables. |
62 |
|
|
|
63 |
|
|
use netcdf95, only: nf95_open, nf95_close, nf90_nowrite |
64 |
|
|
use o3_Mob_ph_m, only: o3_Mob_ph |
65 |
|
|
use phyetat0_m, only: rlat |
66 |
|
|
|
67 |
|
|
! Variables local to the procedure: |
68 |
|
|
integer ncid ! for NetCDF |
69 |
|
|
|
70 |
|
|
real a6(klon, llm, 12) |
71 |
|
|
! (derivative of P_net_Mob with respect to column-density of ozone |
72 |
|
|
! above, in cm2 s-1) |
73 |
|
|
! (On the "physics" grid. |
74 |
|
|
! Third dimension is the number of the month in the year. |
75 |
|
|
! "a6(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)", |
76 |
|
|
! middle of layer "k".) |
77 |
|
|
|
78 |
|
|
real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg |
79 |
|
|
|
80 |
|
|
real, parameter:: Clx = 3.8e-9 |
81 |
|
|
! (total chlorine content in the upper stratosphere) |
82 |
|
|
|
83 |
|
|
integer k, month |
84 |
|
|
|
85 |
|
|
!------------------------------------ |
86 |
|
|
|
87 |
|
|
print *, "Call sequence information: read_coefoz" |
88 |
|
|
|
89 |
|
|
call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) |
90 |
|
|
|
91 |
|
|
a2 = o3_Mob_ph(ncid, "a2") |
92 |
|
|
a4_mass = o3_Mob_ph(ncid, "a4") * 48. / 29. |
93 |
|
|
a6 = o3_Mob_ph(ncid, "a6") |
94 |
|
|
|
95 |
|
|
! Compute "a6_mass" avoiding underflow, do not divide by 1e4 |
96 |
|
|
! before dividing by molecular mass: |
97 |
|
|
a6_mass = a6 / (1e4 * 29. * amu) |
98 |
|
|
! (factor 1e4: conversion from cm2 to m2) |
99 |
|
|
|
100 |
|
|
c_Mob = 48. / 29. * (o3_Mob_ph(ncid, "P_net_Mob") & |
101 |
|
|
- a2 * o3_Mob_ph(ncid, "r_Mob") - a6 * o3_Mob_ph(ncid, "Sigma_Mob")) & |
102 |
|
|
- a4_mass * o3_Mob_ph(ncid, "temp_Mob") |
103 |
|
|
|
104 |
|
|
r_het_interm = o3_Mob_ph(ncid, "R_Het") * (Clx / 3.8e-9)**2 |
105 |
|
|
! Heterogeneous chemistry is only at high latitudes: |
106 |
|
|
forall (k = 1: llm, month = 1: 12) |
107 |
|
|
where (abs(rlat) <= 45.) r_het_interm(:, k, month) = 0. |
108 |
|
|
end forall |
109 |
|
|
|
110 |
|
|
call nf95_close(ncid) |
111 |
|
|
|
112 |
|
|
end subroutine read_coefoz |
113 |
|
|
|
114 |
|
|
end module read_coefoz_m |