/[lmdze]/trunk/dyn3d/disvert.f
ViewVC logotype

Contents of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 313 - (show annotations)
Mon Dec 10 15:54:30 2018 UTC (5 years, 5 months ago) by guez
File size: 7012 byte(s)
Remove module temps. Move variable itau_dyn from module temps to
module dynetat0_m, where it is defined.

Split module dynetat0_m into dynetat0_m and dynetat0_chosen_m. The
motivation is to create smaller modules. Procedures principal_cshift
and invert_zoomx had to stay in dynetat0_m because of circular
dependency. Now we will be able to move them away. Module variables
which are chosen by the user, not computed, in program ce0l go to
dynetat0_chosen_m: day_ref, annee_ref, clon, clat, grossismx,
grossismy, dzoomx, dzoomy, taux, tauy.

Move variable "pa" from module disvert_m to module
dynetat0_chosen_m. Define "pa" in dynetat0_chosen rather than etat0.

Define day_ref and annee_ref in procedure read_serre rather than
etat0.

1 module disvert_m
2
3 use dimensions, only: llm
4
5 implicit none
6
7 private llm, hybrid, funcd, y, ya, compute_ab
8
9 real, save:: ap(llm+1) ! in Pa
10 real, save:: bp(llm+1)
11
12 REAL s(llm+1)
13 ! "s(l)" is the atmospheric hybrid sigma-pressure coordinate at
14 ! half-level, between layers "l" and "l-1"
15
16 real, save:: presnivs(llm)
17 ! approximate full level pressure for a reference surface pressure, in Pa
18
19 real, parameter:: preff = 101325. ! in Pa
20 real y, ya ! for the hybrid function
21
22 contains
23
24 SUBROUTINE disvert
25
26 ! From dyn3d/disvert.F, version 1.1.1.1, 2004/05/19 12:53:05
27 ! Author: P. Le Van
28
29 ! This procedure sets the vertical grid. It defines the host
30 ! variables "ap", "bp", "presnivs". "pa" should be defined before
31 ! this procedure is called.
32
33 ! Libraries:
34 use jumble, only: read_column, new_unit
35 use nr_util, only: pi, assert
36
37 use dynetat0_chosen_m, only: pa
38 use unit_nml_m, only: unit_nml
39
40 ! Local:
41
42 real ds(llm)
43 ! "ds(l)" : \'epaisseur de la couche "l" dans la coordonn\'ee "s"
44
45 INTEGER l, unit
46 REAL x(llm)
47 real:: dsigmin = 1.
48 real zz(llm) ! in km
49
50 character(len=20):: vert_sampling = "tropo"
51 ! Allowed values: "tropo", "strato_custom", "strato",
52 ! "read_hybrid", "read_pressure".
53
54 ! These variables are used only in the case vert_sampling ==
55 ! "strato_custom", and all are in km:
56 real:: vert_scale = 7. ! scale height
57 real:: vert_dzmin = 0.017 ! width of first layer
58 real:: vert_dzlow = 1. ! dz in the low atmosphere
59 real:: vert_z0low = 8.7 ! height at which resolution reaches dzlow
60 real:: vert_dzmid = 2. ! dz in the mid atmosphere
61 real:: vert_z0mid = 70. ! height at which resolution reaches dzmid
62 real:: vert_h_mid = 20. ! width of the transition
63 real:: vert_dzhig = 11. ! dz in the high atmosphere
64 real:: vert_z0hig = 75. ! height at which resolution reaches dz
65 real:: vert_h_hig = 20. ! width of the transition
66
67 real, allocatable:: p(:) ! (2:llm or llm + 1) pressure (in hPa)
68
69 namelist /disvert_nml/vert_sampling, vert_scale, vert_dzmin, vert_dzlow, &
70 vert_z0low, vert_dzmid, vert_z0mid, vert_h_mid, vert_dzhig, &
71 vert_z0hig, vert_h_hig, dsigmin
72
73 !-----------------------------------------------------------------------
74
75 print *, "Call sequence information: disvert"
76
77 write(unit=*, nml=disvert_nml)
78 print *, "Enter namelist 'disvert_nml'."
79 read(unit=*, nml=disvert_nml)
80 write(unit_nml, nml=disvert_nml)
81
82 s(1) = 1.
83 s(llm+1) = 0.
84
85 select case (vert_sampling)
86
87 case ("tropo")
88 ! with llm = 19 and dsigmin = 1 for CMIP 3
89
90 forall (l = 1: llm) ds(l) &
91 = dsigmin + 7. * SIN(pi * (REAL(l) - 0.5) / real(llm + 1))**2
92 ds = ds / sum(ds)
93
94 DO l = llm, 2, -1
95 s(l) = s(l+1) + ds(l)
96 ENDDO
97
98 call compute_ab
99
100 case ("strato")
101 ! with llm = 39 and dsigmin = 0.3 for CMIP5
102
103 forall (l = 1: llm) x(l) = pi * (l - 0.5) / (llm + 1)
104
105 ds = (dsigmin + 7. * SIN(x)**2) * (1. - tanh(2 * x / pi - 1.))**2 / 4.
106 ds = ds / sum(ds)
107
108 DO l = llm, 2, -1
109 s(l) = s(l+1) + ds(l)
110 ENDDO
111
112 call compute_ab
113
114 case ("strato_custom")
115 ! with llm = 79 for CMIP 6
116
117 zz(1) = 0.
118 DO l = 1, llm - 1
119 zz(l + 1) = zz(l) + vert_dzmin + vert_dzlow &
120 * TANH(zz(l) / vert_z0low) + (vert_dzmid - vert_dzlow) &
121 * (TANH((zz(l) - vert_z0mid) / vert_h_mid) &
122 - TANH(- vert_z0mid / vert_h_mid)) &
123 + (vert_dzhig - vert_dzmid - vert_dzlow) &
124 * (TANH((zz(l) - vert_z0hig) / vert_h_hig) &
125 - TANH(- vert_z0hig / vert_h_hig))
126 ENDDO
127
128 allocate(p(2: llm))
129 p = preff * EXP(- zz(2:) / vert_scale)
130 ya = pa / preff
131 s(2: llm) = hybrid(p)
132
133 call compute_ab
134
135 case("read_hybrid")
136 ! Read "ap" and "bp". First line is skipped (title line). "ap"
137 ! should be in Pa. First couple of values should correspond to
138 ! the surface, that is : "bp" should be in descending order.
139 call new_unit(unit)
140 open(unit, file="hybrid.csv", status="old", action="read", &
141 position="rewind")
142 read(unit, fmt=*) ! skip title line
143 do l = 1, llm + 1
144 read(unit, fmt=*) ap(l), bp(l)
145 end do
146 close(unit)
147 ! Quick check:
148 call assert(ap(1) == 0., ap(llm + 1) == 0., bp(1) == 1., &
149 bp(llm + 1) == 0., "disvert: bad ap or bp values")
150 s(2: llm) = ap(2: llm) / pa + bp(2: llm)
151
152 case("read_pressure")
153 ! Read pressure values, in Pa, in descending order, from preff
154 ! to 0. First line is skipped (title line).
155 call read_column("pressure.txt", p, first=2)
156 call assert(size(p) == llm + 1, "disvert: bad number of pressure values")
157 ! Quick check:
158 call assert(p(1) == preff, p(llm + 1) == 0., &
159 "disvert: bad pressure values")
160 ya = pa / preff
161 s(2: llm) = hybrid(p(2: llm))
162 call compute_ab
163
164 case default
165 print *, 'Wrong value for "vert_sampling"'
166 stop 1
167
168 END select
169
170 forall (l = 1: llm) presnivs(l) = 0.5 &
171 * (ap(l) + bp(l) * preff + ap(l+1) + bp(l+1) * preff)
172
173 END SUBROUTINE disvert
174
175 !**********************************************************
176
177 subroutine compute_ab
178
179 use dynetat0_chosen_m, only: pa
180
181 ! Calcul de "ap" et "bp".
182
183 where (s >= 1. / sqrt(1. - log(tiny(0.))))
184 bp = exp(1. - 1. / s**2)
185 elsewhere
186 bp = 0.
187 end where
188
189 ap = pa * (s - bp)
190
191 end subroutine compute_ab
192
193 !**********************************************************
194
195 function hybrid(p)
196
197 ! This procedure computes the hybrid sigma-pressure coordinate
198 ! from pressure values, assuming some reference surface
199 ! pressure. The procedure assumes, and does not check, that
200 ! pressure values are given in descending order.
201
202 use numer_rec_95, only: rtsafe
203
204 real, intent(in):: p(:) ! pressure (in hPa)
205 real hybrid(size(p)) ! hybrid sigma-pressure coordinate
206
207 ! Local:
208 integer l
209
210 !-------------------------------------------------------
211
212 y = p(1) / preff
213 hybrid(1) = rtsafe(funcd, x1 = 0., x2 = 1., xacc = 1e-4)
214
215 do l = 2, size(p)
216 y = p(l) / preff
217 ! Assuming descending order in pressure:
218 hybrid(l) = rtsafe(funcd, x1 = 0., x2 = hybrid(l - 1), &
219 xacc = hybrid(l - 1) * 1e-4)
220 end do
221
222 end function hybrid
223
224 !**********************************************************
225
226 SUBROUTINE funcd(s, fval, fderiv)
227
228 REAL, INTENT(IN):: s
229 REAL, INTENT(OUT):: fval, fderiv
230
231 ! Local:
232 real b
233
234 !------------------------------------
235
236 if (s**3 > 1. / huge(0.)) then
237 b = exp(1. - 1. / s**2)
238 fval = ya * s + (1. - ya) * b - y
239 fderiv = ya + 2 * (1. - ya) * b / s**3
240 else
241 fval = ya * s - y
242 fderiv = ya
243 end if
244
245 END SUBROUTINE funcd
246
247 end module disvert_m

  ViewVC Help
Powered by ViewVC 1.1.21