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

Annotation of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 253 - (hide annotations)
Tue Jan 23 15:49:10 2018 UTC (6 years, 3 months ago) by guez
Original Path: trunk/Sources/dyn3d/disvert.f
File size: 6845 byte(s)
No need for intermediary variables rlong and rlat in inithist.
1 guez 66 module disvert_m
2 guez 3
3     use dimens_m, only: llm
4    
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     ! pressions approximatives des milieux de couches, en Pa
18    
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     ! Allowed values: "tropo", "param", "strato", "read_hybrid", "read_pressure"
49 guez 3
50 guez 99 ! These variables are used only in the case vert_sampling ==
51     ! "param", and all are in km:
52     real:: vert_scale = 7. ! scale height
53     real:: vert_dzmin = 0.02 ! width of first layer
54     real:: vert_dzlow = 1. ! dz in the low atmosphere
55     real:: vert_z0low = 8. ! height at which resolution reaches dzlow
56     real:: vert_dzmid = 3. ! dz in the mid atmosphere
57     real:: vert_z0mid = 70. ! height at which resolution reaches dzmid
58     real:: vert_h_mid = 20. ! width of the transition
59     real:: vert_dzhig = 11. ! dz in the high atmosphere
60     real:: vert_z0hig = 80. ! height at which resolution reaches dz
61     real:: vert_h_hig = 20. ! width of the transition
62 guez 3
63 guez 253 real, allocatable:: p(:) ! (llm + 1) pressure (in hPa)
64 guez 3
65 guez 99 namelist /disvert_nml/vert_sampling, vert_scale, vert_dzmin, vert_dzlow, &
66     vert_z0low, vert_dzmid, vert_z0mid, vert_h_mid, vert_dzhig, &
67     vert_z0hig, vert_h_hig, dsigmin
68 guez 3
69     !-----------------------------------------------------------------------
70    
71     print *, "Call sequence information: disvert"
72    
73 guez 99 write(unit=*, nml=disvert_nml)
74 guez 3 print *, "Enter namelist 'disvert_nml'."
75     read(unit=*, nml=disvert_nml)
76 guez 57 write(unit_nml, nml=disvert_nml)
77 guez 3
78 guez 99 s(1) = 1.
79     s(llm+1) = 0.
80    
81 guez 78 select case (vert_sampling)
82 guez 99
83 guez 65 case ("tropo")
84 guez 99 ! with llm = 19 for CMIP 3
85    
86 guez 3 forall (l = 1: llm) ds(l) &
87 guez 99 = dsigmin + 7. * SIN(pi * (REAL(l) - 0.5) / real(llm + 1))**2
88 guez 3 ds = ds / sum(ds)
89    
90     DO l = llm, 2, -1
91     s(l) = s(l+1) + ds(l)
92     ENDDO
93    
94 guez 78 call compute_ab
95 guez 99
96 guez 65 case ("strato")
97 guez 99 ! with llm = 39 and dsigmin = 0.3 for CMIP5
98    
99 guez 3 forall (l = 1: llm) x(l) = pi * (l - 0.5) / (llm + 1)
100    
101 guez 99 ds = (dsigmin + 7. * SIN(x)**2) * (1. - tanh(2 * x / pi - 1.))**2 / 4.
102 guez 3 ds = ds / sum(ds)
103    
104     DO l = llm, 2, -1
105     s(l) = s(l+1) + ds(l)
106     ENDDO
107 guez 78
108     call compute_ab
109 guez 99
110     case ("param")
111     ! with llm = 79 for CMIP 6
112    
113     zz(1) = 0.
114     DO l = 1, llm - 1
115     zz(l + 1) = zz(l) + vert_dzmin + vert_dzlow &
116     * TANH(zz(l) / vert_z0low) + (vert_dzmid - vert_dzlow) &
117     * (TANH((zz(l) - vert_z0mid) / vert_h_mid) &
118     - TANH(- vert_z0mid / vert_h_mid)) &
119     + (vert_dzhig - vert_dzmid - vert_dzlow) &
120     * (TANH((zz(l) - vert_z0hig) / vert_h_hig) &
121     - TANH(- vert_z0hig / vert_h_hig))
122     ENDDO
123    
124     allocate(p(2: llm))
125     p = preff * EXP(- zz(2:) / vert_scale)
126     ya = pa / preff
127     s(2: llm) = hybrid(p)
128 guez 253
129 guez 99 call compute_ab
130    
131     case("read_hybrid")
132 guez 78 ! Read "ap" and "bp". First line is skipped (title line). "ap"
133     ! should be in Pa. First couple of values should correspond to
134     ! the surface, that is : "bp" should be in descending order.
135 guez 42 call new_unit(unit)
136     open(unit, file="hybrid.csv", status="old", action="read", &
137     position="rewind")
138     read(unit, fmt=*) ! skip title line
139     do l = 1, llm + 1
140 guez 78 read(unit, fmt=*) ap(l), bp(l)
141 guez 42 end do
142     close(unit)
143     ! Quick check:
144 guez 78 call assert(ap(1) == 0., ap(llm + 1) == 0., bp(1) == 1., &
145     bp(llm + 1) == 0., "disvert: bad ap or bp values")
146 guez 99 s(2: llm) = ap(2: llm) / pa + bp(2: llm)
147    
148     case("read_pressure")
149     ! Read pressure values, in Pa, in descending order, from preff
150     ! to 0. First line is skipped (title line).
151     call read_column("pressure.txt", p, first=2)
152     call assert(size(p) == llm + 1, "disvert: bad number of pressure values")
153     ! Quick check:
154     call assert(p(1) == preff, p(llm + 1) == 0., &
155     "disvert: bad pressure values")
156     ya = pa / preff
157     s(2: llm) = hybrid(p(2: llm))
158     call compute_ab
159    
160 guez 3 case default
161 guez 78 print *, 'Wrong value for "vert_sampling"'
162 guez 3 stop 1
163 guez 99
164 guez 3 END select
165    
166 guez 32 forall (l = 1: llm) presnivs(l) = 0.5 &
167     * (ap(l) + bp(l) * preff + ap(l+1) + bp(l+1) * preff)
168 guez 3
169 guez 99 END SUBROUTINE disvert
170 guez 78
171 guez 99 !**********************************************************
172 guez 78
173 guez 99 subroutine compute_ab
174 guez 78
175 guez 99 ! Calcul de "ap" et "bp".
176 guez 78
177 guez 99 where (s >= 1. / sqrt(1. - log(tiny(0.))))
178     bp = exp(1. - 1. / s**2)
179     elsewhere
180     bp = 0.
181     end where
182 guez 3
183 guez 99 ap = pa * (s - bp)
184    
185     end subroutine compute_ab
186    
187     !**********************************************************
188    
189     function hybrid(p)
190    
191     ! This procedure computes the hybrid sigma-pressure coordinate
192     ! from pressure values, assuming some reference surface
193     ! pressure. The procedure assumes, and does not check, that
194     ! pressure values are given in descending order.
195    
196     use numer_rec_95, only: rtsafe
197    
198     real, intent(in):: p(:) ! pressure (in hPa)
199     real hybrid(size(p)) ! hybrid sigma-pressure coordinate
200    
201     ! Local:
202     integer l
203    
204     !-------------------------------------------------------
205    
206     y = p(1) / preff
207     hybrid(1) = rtsafe(funcd, x1 = 0., x2 = 1., xacc = 1e-4)
208    
209     do l = 2, size(p)
210     y = p(l) / preff
211     ! Assuming descending order in pressure:
212     hybrid(l) = rtsafe(funcd, x1 = 0., x2 = hybrid(l - 1), &
213     xacc = hybrid(l - 1) * 1e-4)
214     end do
215    
216     end function hybrid
217    
218     !**********************************************************
219    
220     SUBROUTINE funcd(s, fval, fderiv)
221    
222     REAL, INTENT(IN):: s
223     REAL, INTENT(OUT):: fval, fderiv
224    
225     ! Local:
226     real b
227    
228     !------------------------------------
229    
230     if (s**3 > 1. / huge(0.)) then
231     b = exp(1. - 1. / s**2)
232     fval = ya * s + (1. - ya) * b - y
233     fderiv = ya + 2 * (1. - ya) * b / s**3
234     else
235     fval = ya * s - y
236     fderiv = ya
237     end if
238    
239     END SUBROUTINE funcd
240    
241 guez 66 end module disvert_m

  ViewVC Help
Powered by ViewVC 1.1.21