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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 17 - (show annotations)
Tue Aug 5 13:31:32 2008 UTC (15 years, 9 months ago) by guez
File size: 4224 byte(s)
Created rule for "compare_sampl_*" files in
"Documentation/Manuel_LMDZE.texfol/Graphiques/GNUmakefile".

Extracted "qcheck", "radiornpb", "minmaxqfi" into separate files.

Read pressure coordinate of ozone coefficients once per run instead of
every day.

Added some "intent" attributes.

Added argument "nq" to "ini_histday". Replaced calls to "gr_fi_ecrit"
by calls to "gr_phy_write_2d". "Sigma_O3_Royer" is written to
"histday.nc" only if "nq >= 4". Moved "ini_histrac" to module
"ini_hist".

Compute "zmasse" in "physiq", pass it to "phytrac".

Removed computations of "pftsol*" and "ppsrf*" in "phytrac".

Do not use variable "rg" from module "YOMCST" in "TLIFT".

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

  ViewVC Help
Powered by ViewVC 1.1.21