/[lmdze]/trunk/phylmd/Mobidic/regr_pr_comb_coefoz.f
ViewVC logotype

Annotation of /trunk/phylmd/Mobidic/regr_pr_comb_coefoz.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (hide annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 3 months ago) by guez
File size: 3986 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

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 90 use netcdf, only: nf90_nowrite
50 guez 18 use netcdf95, only: nf95_open, nf95_close
51 guez 3 use phyetat0_m, only: rlat
52 guez 90 use regr_pr_av_m, only: regr_pr_av
53     use regr_pr_int_m, only: regr_pr_int
54 guez 3
55 guez 7 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
56    
57 guez 3 ! Variables local to the procedure:
58 guez 18
59 guez 3 integer ncid ! for NetCDF
60    
61 guez 7 real coefoz(klon, llm)
62     ! (temporary storage for an ozone coefficient)
63     ! (On the "physics" grid.
64     ! "coefoz(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
65     ! middle of layer "k".)
66    
67     real a6(klon, llm)
68     ! (derivative of "P_net_Mob" with respect to column-density of ozone
69 guez 3 ! above, in cm2 s-1)
70     ! (On the "physics" grid.
71 guez 7 ! "a6(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
72 guez 3 ! middle of layer "k".)
73    
74     real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg
75    
76     real, parameter:: Clx = 3.8e-9
77     ! (total chlorine content in the upper stratosphere)
78    
79 guez 7 integer k
80 guez 3
81     !------------------------------------
82    
83 guez 10 print *, "Call sequence information: regr_pr_comb_coefoz"
84 guez 3
85     call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
86    
87 guez 90 call regr_pr_av(ncid, "a2", julien, a2)
88 guez 3
89 guez 90 call regr_pr_av(ncid, "a4", julien, a4_mass)
90 guez 7 a4_mass = a4_mass * 48. / 29.
91    
92 guez 90 call regr_pr_av(ncid, "a6", julien, a6)
93 guez 7
94 guez 3 ! Compute "a6_mass" avoiding underflow, do not divide by 1e4
95     ! before dividing by molecular mass:
96     a6_mass = a6 / (1e4 * 29. * amu)
97     ! (factor 1e4: conversion from cm2 to m2)
98    
99 guez 7 ! Combine coefficients to get "c_Mob":
100     ! (We use as few local variables as possible, in order to spare
101     ! main memory.)
102 guez 3
103 guez 90 call regr_pr_av(ncid, "P_net_Mob", julien, c_Mob)
104 guez 7
105 guez 90 call regr_pr_av(ncid, "r_Mob", julien, coefoz)
106 guez 9 c_mob = c_mob - a2 * coefoz
107 guez 7
108 guez 90 call regr_pr_int(ncid, "Sigma_Mob", julien, top_value=0., v3=coefoz)
109 guez 7 c_mob = (c_mob - a6 * coefoz) * 48. / 29.
110    
111 guez 90 call regr_pr_av(ncid, "temp_Mob", julien, coefoz)
112 guez 7 c_mob = c_mob - a4_mass * coefoz
113    
114 guez 90 call regr_pr_av(ncid, "R_Het", julien, r_het_interm)
115 guez 3 ! Heterogeneous chemistry is only at high latitudes:
116 guez 7 forall (k = 1: llm)
117     where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.
118 guez 3 end forall
119 guez 11 r_het_interm = r_het_interm * (Clx / 3.8e-9)**2
120 guez 3
121     call nf95_close(ncid)
122    
123 guez 7 end subroutine regr_pr_comb_coefoz
124 guez 3
125 guez 7 end module regr_pr_comb_coefoz_m

  ViewVC Help
Powered by ViewVC 1.1.21