source: LMDZ5/trunk/libf/phylmd/o3_chem_m.F90 @ 4011

Last change on this file since 4011 was 2346, checked in by Ehouarn Millour, 9 years ago

Physics/dynamics separation:

  • remove all references to dimensions.h from physics. nbp_lon (==iim) , nbp_lat (==jjm+1) and nbp_lev (==llm) from mod_grid_phy_lmdz should be used instead.
  • added module regular_lonlat_mod in phy_common to store information about the global (lon-lat) grid cell boundaries and centers.

EM

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 5.4 KB
Line 
1! $Id$
2module o3_chem_m
3
4  IMPLICIT none
5
6  private o3_prod
7
8contains
9
10  subroutine o3_chem(julien, gmtime, t_seri, zmasse, pdtphys, rlat, rlon, q)
11
12    ! This procedure evolves the ozone mass fraction through a time
13    ! step taking only chemistry into account.
14
15    ! All the 2-dimensional arrays are on the partial "physics" grid.
16    ! Their shape is "(/klon, nbp_lev/)".
17    ! Index "(i, :)" is for longitude "rlon(i)", latitude "rlat(i)".
18
19    use assert_m, only: assert
20    use dimphy, only: klon
21    use regr_pr_comb_coefoz_m, only: c_Mob, a4_mass, a2, r_het_interm
22    use mod_grid_phy_lmdz, only: nbp_lev
23    use nrtype, only: pi
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, nbp_lev) temperature, in K
28
29    real, intent(in):: zmasse(:, :) ! (klon, nbp_lev)
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(in):: rlat(:), rlon(:)
36    ! (longitude and latitude of each horizontal position, in degrees)
37
38    real, intent(inout):: q(:, :) ! (klon, nbp_lev) mass fraction of ozone
39    ! "q(:, k)" is at middle of layer "k".)
40
41    ! Variables local to the procedure:
42    ! (for "pi")
43    integer k
44
45    real c(klon, nbp_lev)
46    ! (constant term during a time step in the net mass production
47    ! rate of ozone by chemistry, per unit mass of air, in s-1)
48    ! "c(:, k)" is at middle of layer "k".)
49
50    real b(klon, nbp_lev)
51    ! (coefficient of "q" in the net mass production
52    ! rate of ozone by chemistry, per unit mass of air, in s-1)
53    ! "b(:, k)" is at middle of layer "k".)
54
55    real dq_o3_chem(klon, nbp_lev)
56    ! (variation of ozone mass fraction due to chemistry during a time step)
57    ! "dq_o3_chem(:, k)" is at middle of layer "k".)
58
59    real earth_long
60    ! (longitude vraie de la Terre dans son orbite solaire, par
61    ! rapport au point vernal (21 mars), en degrés)
62
63    real pmu0(klon) ! mean of cosine of solar zenith angle during "pdtphys"
64    real trash1
65    real trash2(klon)
66
67    !-------------------------------------------------------------
68
69    call assert(klon == (/size(q, 1), size(t_seri, 1), size(zmasse, 1), &
70         size(rlat), size(rlon)/), "o3_chem klon")
71    call assert(nbp_lev == (/size(q, 2), size(t_seri, 2), size(zmasse, 2)/), &
72         "o3_chem nbp_lev")
73
74    c = c_Mob + a4_mass * t_seri
75
76    ! Compute coefficient "b":
77
78    ! Heterogeneous chemistry is only at low temperature:
79    where (t_seri < 195.)
80       b = r_het_interm
81    elsewhere
82       b = 0.
83    end where
84
85    ! Heterogeneous chemistry is only during daytime:
86    call orbite(real(julien), earth_long, trash1)
87    call zenang(earth_long, gmtime, pdtphys, rlat, rlon, pmu0, trash2)
88    forall (k = 1: nbp_lev)
89       where (pmu0 <= cos(87. / 180. * pi)) b(:, k) = 0.
90    end forall
91
92    b = b + a2
93
94    ! Midpoint method:
95
96    ! Trial step to the midpoint:
97    dq_o3_chem = o3_prod(q, zmasse, c, b) * pdtphys  / 2
98    ! "Real" step across the whole interval:
99    dq_o3_chem = o3_prod(q + dq_o3_chem, zmasse, c, b) * pdtphys
100    q = q + dq_o3_chem
101
102    ! Confine the mass fraction:
103    q = min(max(q, 0.), .01)
104
105  end subroutine o3_chem
106
107  !*************************************************
108
109  function o3_prod(q, zmasse, c, b)
110
111    ! This function computes the production rate of ozone by chemistry.
112
113    ! All the 2-dimensional arrays are on the partial "physics" grid.
114    ! Their shape is "(/klon, nbp_lev/)".
115    ! Index "(i, :)" is for longitude "rlon(i)", latitude "rlat(i)".
116
117    use regr_pr_comb_coefoz_m, only: a6_mass
118    use assert_m, only: assert
119    use dimphy, only: klon
120    use mod_grid_phy_lmdz, only: nbp_lev
121
122    real, intent(in):: q(:, :) ! mass fraction of ozone
123    ! "q(:, k)" is at middle of layer "k".)
124
125    real, intent(in):: zmasse(:, :)
126    ! (column-density of mass of air in a layer, in kg m-2)
127    ! ("zmasse(:, k)" is for layer "k".)
128
129    real, intent(in):: c(:, :)
130    ! (constant term during a time step in the net mass production
131    ! rate of ozone by chemistry, per unit mass of air, in s-1)
132    ! "c(:, k)" is at middle of layer "k".)
133
134    real, intent(in):: b(:, :)
135    ! (coefficient of "q" in the net mass production rate of ozone by
136    ! chemistry, per unit mass of air, in s-1)
137    ! ("b(:, k)" is at middle of layer "k".)
138
139    real o3_prod(klon, nbp_lev)
140    ! (net mass production rate of ozone by chemistry, per unit mass
141    ! of air, in s-1)
142    ! ("o3_prod(:, k)" is at middle of layer "k".)
143
144    ! Variables local to the procedure:
145
146    real sigma_mass(klon, nbp_lev)
147    ! (mass column-density of ozone above point, in kg m-2)
148    ! ("sigma_mass(:, k)" is at middle of layer "k".)
149
150    integer k
151
152    !-------------------------------------------------------------------
153
154    call assert(klon == (/size(q, 1), size(zmasse, 1), size(c, 1), &
155         size(b, 1)/), "o3_prod 1")
156    call assert(nbp_lev == (/size(q, 2), size(zmasse, 2), size(c, 2), &
157         size(b, 2)/), "o3_prod 2")
158
159    ! Compute the column-density above the base of layer
160    ! "k", and, as a first approximation, take it as column-density
161    ! above the middle of layer "k":
162    sigma_mass(:, nbp_lev) = zmasse(:, nbp_lev) * q(:, nbp_lev) ! top layer
163    do k =  nbp_lev - 1, 1, -1
164       sigma_mass(:, k) = sigma_mass(:, k+1) + zmasse(:, k) * q(:, k)
165    end do
166
167    o3_prod = c + b * q + a6_mass * sigma_mass
168
169  end function o3_prod
170
171end module o3_chem_m
Note: See TracBrowser for help on using the repository browser.