/[lmdze]/trunk/libf/phylmd/Mobidic/regr_pr_comb_coefoz.f90
ViewVC logotype

Annotation of /trunk/libf/phylmd/Mobidic/regr_pr_comb_coefoz.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 9 - (hide annotations)
Mon Mar 31 13:58:05 2008 UTC (16 years, 1 month ago) by guez
File size: 5642 byte(s)
New variables "*_dir" in "g95.mk".
Corrected some bugs: "etat0_lim" works, but not "gcm".

1 guez 7 module regr_pr_comb_coefoz_m
2 guez 3
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 guez 7 real, save:: c_Mob(klon, llm)
11 guez 3 ! (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 guez 7 ! "c_Mob(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
15 guez 3 ! middle of layer "k".)
16    
17 guez 7 real, save:: a2(klon, llm)
18 guez 3 ! (derivative of mass production rate of ozone per unit mass of
19     ! air with respect to ozone mass fraction, in s-1)
20     ! (On the "physics" grid.
21 guez 7 ! "a2(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
22 guez 3 ! middle of layer "k".)
23    
24 guez 7 real, save:: a4_mass(klon, llm)
25 guez 3 ! (derivative of mass production rate of ozone per unit mass of
26     ! air with respect to temperature, in s-1 K-1)
27     ! (On the "physics" grid.
28 guez 7 ! "a4_mass(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
29 guez 3 ! middle of layer "k".)
30    
31 guez 7 real, save:: a6_mass(klon, llm)
32 guez 3 ! (derivative of mass production rate of ozone per unit mass of
33     ! air with respect to mass column-density of ozone above, in m2 s-1 kg-1)
34     ! (On the "physics" grid.
35 guez 7 ! "a6_mass(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
36 guez 3 ! middle of layer "k".)
37    
38 guez 7 real, save:: r_het_interm(klon, llm)
39 guez 3 ! (net mass production rate by heterogeneous chemistry, per unit
40     ! mass of ozone, corrected for chlorine content and latitude, but
41     ! not for temperature and sun direction, in s-1)
42     ! (On the "physics" grid.
43 guez 7 ! "r_het_interm(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
44 guez 3 ! middle of layer "k".)
45    
46     private klon, llm
47    
48     contains
49    
50 guez 7 subroutine regr_pr_comb_coefoz(julien)
51 guez 3
52 guez 7 ! "regr_pr_comb_coefoz" stands for "regrid pressure combine
53     ! coefficients ozone".
54 guez 3
55 guez 7 ! This subroutine :
56     ! -- reads from a file all eight parameters for ozone chemistry,
57     ! at the current day ;
58     ! -- regrids the parameters in pressure to the LMDZ vertical grid ;
59     ! -- packs the parameters to the "physics" horizontal grid ;
60     ! -- combines the eight parameters to define the five module variables.
61    
62     ! We assume that, in "coefoz_LMDZ.nc", the pressure levels are in hPa
63     ! and strictly increasing.
64    
65     use netcdf95, only: nf95_open, nf95_close, nf95_get_coord
66     use netcdf, only: nf90_nowrite
67     use regr_pr_coefoz, only: regr_pr_av_coefoz, regr_pr_int_coefoz
68 guez 3 use phyetat0_m, only: rlat
69    
70 guez 7 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
71    
72 guez 3 ! Variables local to the procedure:
73     integer ncid ! for NetCDF
74    
75 guez 7 real, pointer:: plev(:)
76     ! (pressure level of input data, converted to Pa, in strictly
77     ! increasing order)
78    
79     integer n_plev ! number of pressure levels in the input data
80    
81     real, allocatable:: press_in_edg(:)
82     ! (edges of pressure intervals for input data, in Pa, in strictly
83     ! increasing order)
84    
85     real coefoz(klon, llm)
86     ! (temporary storage for an ozone coefficient)
87     ! (On the "physics" grid.
88     ! "coefoz(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
89     ! middle of layer "k".)
90    
91     real a6(klon, llm)
92     ! (derivative of "P_net_Mob" with respect to column-density of ozone
93 guez 3 ! above, in cm2 s-1)
94     ! (On the "physics" grid.
95 guez 7 ! "a6(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
96 guez 3 ! middle of layer "k".)
97    
98     real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg
99    
100     real, parameter:: Clx = 3.8e-9
101     ! (total chlorine content in the upper stratosphere)
102    
103 guez 7 integer k
104 guez 3
105     !------------------------------------
106    
107     print *, "Call sequence information: read_coefoz"
108    
109     call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
110    
111 guez 7 call nf95_get_coord(ncid, "plev", plev)
112     ! Convert from hPa to Pa because "regr_pr_av" and "regr_pr_int" require so:
113     plev = plev * 100.
114     n_plev = size(plev)
115 guez 3
116 guez 7 ! Compute edges of pressure intervals:
117     allocate(press_in_edg(n_plev + 1))
118     press_in_edg(1) = 0.
119     ! We choose edges halfway in logarithm:
120     forall (k = 2:n_plev) press_in_edg(k) = sqrt(plev(k - 1) * plev(k))
121     press_in_edg(n_plev + 1) = huge(0.)
122     ! (infinity, but any value guaranteed to be greater than the
123     ! surface pressure would do)
124    
125     call regr_pr_av_coefoz(ncid, "a2", julien, press_in_edg, a2)
126    
127     call regr_pr_av_coefoz(ncid, "a4", julien, press_in_edg, a4_mass)
128     a4_mass = a4_mass * 48. / 29.
129    
130     call regr_pr_av_coefoz(ncid, "a6", julien, press_in_edg, a6)
131    
132 guez 3 ! Compute "a6_mass" avoiding underflow, do not divide by 1e4
133     ! before dividing by molecular mass:
134     a6_mass = a6 / (1e4 * 29. * amu)
135     ! (factor 1e4: conversion from cm2 to m2)
136    
137 guez 7 ! Combine coefficients to get "c_Mob":
138     ! (We use as few local variables as possible, in order to spare
139     ! main memory.)
140 guez 3
141 guez 7 call regr_pr_av_coefoz(ncid, "P_net_Mob", julien, press_in_edg, c_Mob)
142    
143     call regr_pr_av_coefoz(ncid, "r_Mob", julien, press_in_edg, coefoz)
144 guez 9 c_mob = c_mob - a2 * coefoz
145 guez 7
146     call regr_pr_int_coefoz(ncid, "Sigma_Mob", julien, plev, top_value=0., &
147 guez 9 v3=coefoz)
148 guez 7 c_mob = (c_mob - a6 * coefoz) * 48. / 29.
149    
150     call regr_pr_av_coefoz(ncid, "temp_Mob", julien, press_in_edg, coefoz)
151     c_mob = c_mob - a4_mass * coefoz
152    
153     call regr_pr_av_coefoz(ncid, "R_Het", julien, press_in_edg, r_het_interm)
154 guez 3 ! Heterogeneous chemistry is only at high latitudes:
155 guez 7 forall (k = 1: llm)
156     where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.
157 guez 3 end forall
158 guez 7 where (r_het_interm /= 0.) r_het_interm = r_het_interm * (Clx / 3.8e-9)**2
159 guez 3
160 guez 7 deallocate(plev) ! pointer
161 guez 3 call nf95_close(ncid)
162    
163 guez 7 end subroutine regr_pr_comb_coefoz
164 guez 3
165 guez 7 end module regr_pr_comb_coefoz_m

  ViewVC Help
Powered by ViewVC 1.1.21