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

Annotation of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 266 - (hide annotations)
Thu Apr 19 17:54:55 2018 UTC (6 years, 1 month ago) by guez
File size: 6921 byte(s)
Define macros of the preprocessor CPP_IIM, CPP_JJM, CPP_LLM so we can
control the resolution from the compilation command, and automate
compilation for several resolutions.

In module yoethf_m, transform variables into named constants. So we do
not need procedure yoethf any longer.

Bug fix in program test_inter_barxy, missing calls to fyhyp and fxhyp,
and definition of rlatu.

Remove variable iecri of module conf_gcm_m. The files dyn_hist*.nc are
written every time step. We are simplifying the output system, pending
replacement by a whole new system.

Modify possible value of vert_sampling from "param" to
"strato_custom", following LMDZ. Default values of corresponding
namelist variables are now the values used for LMDZ CMIP6.

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

  ViewVC Help
Powered by ViewVC 1.1.21