/[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 18 - (show annotations)
Thu Aug 7 12:29:13 2008 UTC (15 years, 9 months ago) by guez
File size: 4103 byte(s)
In module "regr_pr", rewrote scanning of horizontal positions as a
single set of loops, using a mask.

Added some "intent" attributes.

In "dynredem0", replaced calls to Fortran 77 interface of NetCDF by
calls to NetCDF95. Removed calls to "nf_redef", regrouped all writing
operations. In "dynredem1", replaced some calls to Fortran 77
interface of NetCDF by calls to Fortran 90 interface.

Renamed variable "nqmax" to "nq_phys".

In "physiq", if "nq >= 5" then "wo" is computed from the
parameterization of "Cariolle".

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

  ViewVC Help
Powered by ViewVC 1.1.21