/[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 265 - (hide annotations)
Tue Mar 20 09:35:59 2018 UTC (6 years, 2 months ago) by guez
File size: 4266 byte(s)
Rename module dimens_m to dimensions.
1 guez 7 module regr_pr_comb_coefoz_m
2 guez 3
3 guez 265 use dimensions, only: llm
4 guez 3 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 162 subroutine regr_pr_comb_coefoz(julien, paprs, pplay)
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 162 real, intent(in):: paprs(:, :) ! (klon, llm + 1)
58     ! (pression pour chaque inter-couche, en Pa)
59    
60     real, intent(in):: pplay(:, :) ! (klon, llm)
61     ! (pression pour le mileu de chaque couche, en Pa)
62    
63 guez 3 ! Variables local to the procedure:
64 guez 18
65 guez 3 integer ncid ! for NetCDF
66    
67 guez 7 real coefoz(klon, llm)
68     ! (temporary storage for an ozone coefficient)
69     ! (On the "physics" grid.
70     ! "coefoz(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
71     ! middle of layer "k".)
72    
73     real a6(klon, llm)
74     ! (derivative of "P_net_Mob" with respect to column-density of ozone
75 guez 3 ! above, in cm2 s-1)
76     ! (On the "physics" grid.
77 guez 7 ! "a6(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
78 guez 3 ! middle of layer "k".)
79    
80     real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg
81    
82     real, parameter:: Clx = 3.8e-9
83     ! (total chlorine content in the upper stratosphere)
84    
85 guez 7 integer k
86 guez 3
87     !------------------------------------
88    
89 guez 10 print *, "Call sequence information: regr_pr_comb_coefoz"
90 guez 3
91     call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
92    
93 guez 162 call regr_pr_av(ncid, "a2", julien, paprs, a2)
94 guez 3
95 guez 162 call regr_pr_av(ncid, "a4", julien, paprs, a4_mass)
96 guez 7 a4_mass = a4_mass * 48. / 29.
97    
98 guez 162 call regr_pr_av(ncid, "a6", julien, paprs, a6)
99 guez 7
100 guez 3 ! Compute "a6_mass" avoiding underflow, do not divide by 1e4
101     ! before dividing by molecular mass:
102     a6_mass = a6 / (1e4 * 29. * amu)
103     ! (factor 1e4: conversion from cm2 to m2)
104    
105 guez 7 ! Combine coefficients to get "c_Mob":
106     ! (We use as few local variables as possible, in order to spare
107     ! main memory.)
108 guez 3
109 guez 162 call regr_pr_av(ncid, "P_net_Mob", julien, paprs, c_Mob)
110 guez 7
111 guez 162 call regr_pr_av(ncid, "r_Mob", julien, paprs, coefoz)
112 guez 9 c_mob = c_mob - a2 * coefoz
113 guez 7
114 guez 162 call regr_pr_int(ncid, "Sigma_Mob", julien, pplay, top_value=0., v3=coefoz)
115 guez 7 c_mob = (c_mob - a6 * coefoz) * 48. / 29.
116    
117 guez 162 call regr_pr_av(ncid, "temp_Mob", julien, paprs, coefoz)
118 guez 7 c_mob = c_mob - a4_mass * coefoz
119    
120 guez 162 call regr_pr_av(ncid, "R_Het", julien, paprs, r_het_interm)
121 guez 3 ! Heterogeneous chemistry is only at high latitudes:
122 guez 7 forall (k = 1: llm)
123     where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.
124 guez 3 end forall
125 guez 11 r_het_interm = r_het_interm * (Clx / 3.8e-9)**2
126 guez 3
127     call nf95_close(ncid)
128    
129 guez 7 end subroutine regr_pr_comb_coefoz
130 guez 3
131 guez 7 end module regr_pr_comb_coefoz_m

  ViewVC Help
Powered by ViewVC 1.1.21