1 | MODULE dcmip2016_baroclinic_wave_mod |
---|
2 | |
---|
3 | !======================================================================= |
---|
4 | ! |
---|
5 | ! Date: July 29, 2015 |
---|
6 | ! |
---|
7 | ! Functions for setting up idealized initial conditions for the |
---|
8 | ! Ullrich, Melvin, Staniforth and Jablonowski baroclinic instability. |
---|
9 | ! |
---|
10 | ! SUBROUTINE baroclinic_wave_sample( |
---|
11 | ! deep,moist,pertt,X,lon,lat,p,z,zcoords,u,v,w,t,phis,ps,rho,q) |
---|
12 | ! |
---|
13 | ! Options: |
---|
14 | ! deep deep atmosphere (1 = yes or 0 = no) |
---|
15 | ! moist include moisture (1 = yes or 0 = no) |
---|
16 | ! pertt type of perturbation (0 = exponential, 1 = stream function) |
---|
17 | ! X Earth scaling factor |
---|
18 | ! |
---|
19 | ! Given a point specified by: |
---|
20 | ! lon longitude (radians) |
---|
21 | ! lat latitude (radians) |
---|
22 | ! p/z pressure (Pa) / height (m) |
---|
23 | ! zcoords 1 if z is specified, 0 if p is specified |
---|
24 | ! |
---|
25 | ! the functions will return: |
---|
26 | ! p pressure if z is specified and zcoords = 1 (Pa) |
---|
27 | ! u zonal wind (m s^-1) |
---|
28 | ! v meridional wind (m s^-1) |
---|
29 | ! t temperature (K) |
---|
30 | ! thetav virtual potential temperature (K) |
---|
31 | ! phis surface geopotential (m^2 s^-2) |
---|
32 | ! ps surface pressure (Pa) |
---|
33 | ! rho density (kj m^-3) |
---|
34 | ! q water vapor mixing ratio (kg/kg) |
---|
35 | ! |
---|
36 | ! |
---|
37 | ! Author: Paul Ullrich |
---|
38 | ! University of California, Davis |
---|
39 | ! Email: paullrich@ucdavis.edu |
---|
40 | ! |
---|
41 | !======================================================================= |
---|
42 | |
---|
43 | IMPLICIT NONE |
---|
44 | |
---|
45 | !======================================================================= |
---|
46 | ! Physical constants |
---|
47 | !======================================================================= |
---|
48 | |
---|
49 | REAL(8), PARAMETER :: & |
---|
50 | a = 6371220.0d0, & ! Reference Earth's Radius (m) |
---|
51 | Rd = 287.0d0, & ! Ideal gas const dry air (J kg^-1 K^1) |
---|
52 | g = 9.80616d0, & ! Gravity (m s^2) |
---|
53 | cp = 1004.5d0, & ! Specific heat capacity (J kg^-1 K^1) |
---|
54 | Lvap = 2.5d6, & ! Latent heat of vaporization of water |
---|
55 | Rvap = 461.5d0, & ! Ideal gas constnat for water vapor |
---|
56 | Mvap = 0.608d0, & ! Ratio of molar mass of dry air/water |
---|
57 | pi = 3.14159265358979d0, & ! pi |
---|
58 | p0 = 100000.0d0, & ! surface pressure (Pa) |
---|
59 | kappa = 2.d0/7.d0, & ! Ratio of Rd to cp |
---|
60 | omega = 7.29212d-5, & ! Reference rotation rate of the Earth (s^-1) |
---|
61 | deg2rad = pi/180.d0 ! Conversion factor of degrees to radians |
---|
62 | |
---|
63 | !======================================================================= |
---|
64 | ! Test case parameters |
---|
65 | !======================================================================= |
---|
66 | REAL(8), PARAMETER :: & |
---|
67 | T0E = 310.d0 , & ! temperature at equatorial surface (K) |
---|
68 | T0P = 240.d0 , & ! temperature at polar surface (K) |
---|
69 | B = 2.d0 , & ! jet half-width parameter |
---|
70 | K = 3.d0 , & ! jet width parameter |
---|
71 | lapse = 0.005d0 ! lapse rate parameter |
---|
72 | |
---|
73 | REAL(8), PARAMETER :: & |
---|
74 | pertu0 = 0.5d0 , & ! SF Perturbation wind velocity (m/s) |
---|
75 | pertr = 1.d0/6.d0 , & ! SF Perturbation radius (Earth radii) |
---|
76 | pertup = 1.0d0 , & ! Exp. perturbation wind velocity (m/s) |
---|
77 | pertexpr = 0.1d0 , & ! Exp. perturbation radius (Earth radii) |
---|
78 | pertlon = pi/9.d0 , & ! Perturbation longitude |
---|
79 | pertlat = 2.d0*pi/9.d0, & ! Perturbation latitude |
---|
80 | pertz = 15000.d0 , & ! Perturbation height cap |
---|
81 | dxepsilon = 1.d-5 ! Small value for numerical derivatives |
---|
82 | |
---|
83 | REAL(8), PARAMETER :: & |
---|
84 | moistqlat = 2.d0*pi/9.d0, & ! Humidity latitudinal width |
---|
85 | moistqp = 34000.d0, & ! Humidity vertical pressure width |
---|
86 | moisttr = 0.1d0, & ! Vertical cut-off pressure for humidity |
---|
87 | moistqs = 1.d-12, & ! Humidity above cut-off |
---|
88 | moistq0 = 0.018d0, & ! Maximum specific humidity |
---|
89 | moistqr = 0.9d0, & ! Maximum saturation ratio |
---|
90 | moisteps = 0.622d0, & ! Ratio of gas constants |
---|
91 | moistT0 = 273.16d0, & ! Reference temperature (K) |
---|
92 | moistE0Ast = 610.78d0 ! Saturation vapor pressure at T0 (Pa) |
---|
93 | |
---|
94 | CONTAINS |
---|
95 | |
---|
96 | !======================================================================= |
---|
97 | ! Generate the baroclinic instability initial conditions |
---|
98 | !======================================================================= |
---|
99 | SUBROUTINE baroclinic_wave_test(deep,moist,pertt,X,lon,lat,p,z,zcoords,u,v,t,thetav,phis,ps,rho,q) & |
---|
100 | BIND(c, name = "baroclinic_wave_test") |
---|
101 | |
---|
102 | IMPLICIT NONE |
---|
103 | |
---|
104 | !----------------------------------------------------------------------- |
---|
105 | ! input/output params parameters at given location |
---|
106 | !----------------------------------------------------------------------- |
---|
107 | INTEGER, INTENT(IN) :: & |
---|
108 | deep, & ! Deep (1) or Shallow (0) test case |
---|
109 | moist, & ! Moist (1) or Dry (0) test case |
---|
110 | pertt ! Perturbation type |
---|
111 | |
---|
112 | REAL(8), INTENT(IN) :: & |
---|
113 | lon, & ! Longitude (radians) |
---|
114 | lat, & ! Latitude (radians) |
---|
115 | X ! Earth scaling parameter |
---|
116 | |
---|
117 | REAL(8), INTENT(INOUT) :: & |
---|
118 | p, & ! Pressure (Pa) |
---|
119 | z ! Altitude (m) |
---|
120 | |
---|
121 | INTEGER, INTENT(IN) :: zcoords ! 1 if z coordinates are specified |
---|
122 | ! 0 if p coordinates are specified |
---|
123 | |
---|
124 | REAL(8), INTENT(OUT) :: & |
---|
125 | u, & ! Zonal wind (m s^-1) |
---|
126 | v, & ! Meridional wind (m s^-1) |
---|
127 | t, & ! Temperature (K) |
---|
128 | thetav, & ! Virtual potential temperature (K) |
---|
129 | phis, & ! Surface Geopotential (m^2 s^-2) |
---|
130 | ps, & ! Surface Pressure (Pa) |
---|
131 | rho, & ! density (kg m^-3) |
---|
132 | q ! water vapor mixing ratio (kg/kg) |
---|
133 | |
---|
134 | !------------------------------------------------ |
---|
135 | ! Local variables |
---|
136 | !------------------------------------------------ |
---|
137 | REAL(8) :: aref, omegaref |
---|
138 | REAL(8) :: T0, constH, constC, scaledZ, inttau2, rratio |
---|
139 | REAL(8) :: inttermU, bigU, rcoslat, omegarcoslat |
---|
140 | REAL(8) :: eta, qratio, qnum, qden |
---|
141 | |
---|
142 | !------------------------------------------------ |
---|
143 | ! Pressure and temperature |
---|
144 | !------------------------------------------------ |
---|
145 | if (zcoords .eq. 1) then |
---|
146 | CALL evaluate_pressure_temperature(deep, X, lon, lat, z, p, t) |
---|
147 | else |
---|
148 | CALL evaluate_z_temperature(deep, X, lon, lat, p, z, t) |
---|
149 | end if |
---|
150 | |
---|
151 | !------------------------------------------------ |
---|
152 | ! Compute test case constants |
---|
153 | !------------------------------------------------ |
---|
154 | aref = a / X |
---|
155 | omegaref = omega * X |
---|
156 | |
---|
157 | T0 = 0.5d0 * (T0E + T0P) |
---|
158 | |
---|
159 | constH = Rd * T0 / g |
---|
160 | |
---|
161 | constC = 0.5d0 * (K + 2.d0) * (T0E - T0P) / (T0E * T0P) |
---|
162 | |
---|
163 | scaledZ = z / (B * constH) |
---|
164 | |
---|
165 | inttau2 = constC * z * exp(- scaledZ**2) |
---|
166 | |
---|
167 | ! radius ratio |
---|
168 | if (deep .eq. 0) then |
---|
169 | rratio = 1.d0 |
---|
170 | else |
---|
171 | rratio = (z + aref) / aref; |
---|
172 | end if |
---|
173 | |
---|
174 | !----------------------------------------------------- |
---|
175 | ! Initialize surface pressure |
---|
176 | !----------------------------------------------------- |
---|
177 | ps = p0 |
---|
178 | |
---|
179 | !----------------------------------------------------- |
---|
180 | ! Initialize velocity field |
---|
181 | !----------------------------------------------------- |
---|
182 | inttermU = (rratio * cos(lat))**(K - 1.d0) - (rratio * cos(lat))**(K + 1.d0) |
---|
183 | bigU = g / aref * K * inttau2 * inttermU * t |
---|
184 | if (deep .eq. 0) then |
---|
185 | rcoslat = aref * cos(lat) |
---|
186 | else |
---|
187 | rcoslat = (z + aref) * cos(lat) |
---|
188 | end if |
---|
189 | |
---|
190 | omegarcoslat = omegaref * rcoslat |
---|
191 | |
---|
192 | u = - omegarcoslat + sqrt(omegarcoslat**2 + rcoslat * bigU) |
---|
193 | v = 0.d0 |
---|
194 | |
---|
195 | !----------------------------------------------------- |
---|
196 | ! Add perturbation to the velocity field |
---|
197 | !----------------------------------------------------- |
---|
198 | |
---|
199 | ! Exponential type |
---|
200 | if (pertt .eq. 0) then |
---|
201 | u = u + evaluate_exponential(lon, lat, z) |
---|
202 | |
---|
203 | ! Stream function type |
---|
204 | elseif (pertt .eq. 1) then |
---|
205 | u = u - 1.d0 / (2.d0 * dxepsilon) * & |
---|
206 | ( evaluate_streamfunction(lon, lat + dxepsilon, z) & |
---|
207 | - evaluate_streamfunction(lon, lat - dxepsilon, z)) |
---|
208 | |
---|
209 | v = v + 1.d0 / (2.d0 * dxepsilon * cos(lat)) * & |
---|
210 | ( evaluate_streamfunction(lon + dxepsilon, lat, z) & |
---|
211 | - evaluate_streamfunction(lon - dxepsilon, lat, z)) |
---|
212 | end if |
---|
213 | |
---|
214 | !----------------------------------------------------- |
---|
215 | ! Initialize surface geopotential |
---|
216 | !----------------------------------------------------- |
---|
217 | phis = 0.d0 |
---|
218 | |
---|
219 | !----------------------------------------------------- |
---|
220 | ! Initialize density |
---|
221 | !----------------------------------------------------- |
---|
222 | rho = p / (Rd * t) |
---|
223 | |
---|
224 | !----------------------------------------------------- |
---|
225 | ! Initialize specific humidity |
---|
226 | !----------------------------------------------------- |
---|
227 | if (moist .eq. 1) then |
---|
228 | eta = p/p0 |
---|
229 | |
---|
230 | if (eta .gt. moisttr) then |
---|
231 | q = moistq0 * exp(- (lat/moistqlat)**4) & |
---|
232 | * exp(- ((eta-1.d0)*p0/moistqp)**2) |
---|
233 | else |
---|
234 | q = moistqs |
---|
235 | end if |
---|
236 | |
---|
237 | ! Convert virtual temperature to temperature |
---|
238 | t = t / (1.d0 + Mvap * q) |
---|
239 | |
---|
240 | else |
---|
241 | q = 0.d0 |
---|
242 | end if |
---|
243 | |
---|
244 | !----------------------------------------------------- |
---|
245 | ! Initialize virtual potential temperature |
---|
246 | !----------------------------------------------------- |
---|
247 | thetav = t * (1.d0 + 0.61d0 * q) * (p0 / p)**(Rd / cp) |
---|
248 | |
---|
249 | END SUBROUTINE baroclinic_wave_test |
---|
250 | |
---|
251 | !----------------------------------------------------------------------- |
---|
252 | ! Calculate pointwise pressure and temperature |
---|
253 | !----------------------------------------------------------------------- |
---|
254 | SUBROUTINE evaluate_pressure_temperature(deep, X, lon, lat, z, p, t) |
---|
255 | |
---|
256 | INTEGER, INTENT(IN) :: deep ! Deep (1) or Shallow (0) test case |
---|
257 | |
---|
258 | REAL(8), INTENT(IN) :: & |
---|
259 | X, & ! Earth scaling ratio |
---|
260 | lon, & ! Longitude (radians) |
---|
261 | lat, & ! Latitude (radians) |
---|
262 | z ! Altitude (m) |
---|
263 | |
---|
264 | REAL(8), INTENT(OUT) :: & |
---|
265 | p, & ! Pressure (Pa) |
---|
266 | t ! Temperature (K) |
---|
267 | |
---|
268 | REAL(8) :: aref, omegaref |
---|
269 | REAL(8) :: T0, constA, constB, constC, constH, scaledZ |
---|
270 | REAL(8) :: tau1, tau2, inttau1, inttau2 |
---|
271 | REAL(8) :: rratio, inttermT |
---|
272 | |
---|
273 | !-------------------------------------------- |
---|
274 | ! Constants |
---|
275 | !-------------------------------------------- |
---|
276 | aref = a / X |
---|
277 | omegaref = omega * X |
---|
278 | |
---|
279 | T0 = 0.5d0 * (T0E + T0P) |
---|
280 | constA = 1.d0 / lapse |
---|
281 | constB = (T0 - T0P) / (T0 * T0P) |
---|
282 | constC = 0.5d0 * (K + 2.d0) * (T0E - T0P) / (T0E * T0P) |
---|
283 | constH = Rd * T0 / g |
---|
284 | |
---|
285 | scaledZ = z / (B * constH) |
---|
286 | |
---|
287 | !-------------------------------------------- |
---|
288 | ! tau values |
---|
289 | !-------------------------------------------- |
---|
290 | tau1 = constA * lapse / T0 * exp(lapse * z / T0) & |
---|
291 | + constB * (1.d0 - 2.d0 * scaledZ**2) * exp(- scaledZ**2) |
---|
292 | tau2 = constC * (1.d0 - 2.d0 * scaledZ**2) * exp(- scaledZ**2) |
---|
293 | |
---|
294 | inttau1 = constA * (exp(lapse * z / T0) - 1.d0) & |
---|
295 | + constB * z * exp(- scaledZ**2) |
---|
296 | inttau2 = constC * z * exp(- scaledZ**2) |
---|
297 | |
---|
298 | !-------------------------------------------- |
---|
299 | ! radius ratio |
---|
300 | !-------------------------------------------- |
---|
301 | if (deep .eq. 0) then |
---|
302 | rratio = 1.d0 |
---|
303 | else |
---|
304 | rratio = (z + aref) / aref; |
---|
305 | end if |
---|
306 | |
---|
307 | !-------------------------------------------- |
---|
308 | ! interior term on temperature expression |
---|
309 | !-------------------------------------------- |
---|
310 | inttermT = (rratio * cos(lat))**K & |
---|
311 | - K / (K + 2.d0) * (rratio * cos(lat))**(K + 2.d0) |
---|
312 | |
---|
313 | !-------------------------------------------- |
---|
314 | ! temperature |
---|
315 | !-------------------------------------------- |
---|
316 | t = 1.d0 / (rratio**2 * (tau1 - tau2 * inttermT)) |
---|
317 | |
---|
318 | !-------------------------------------------- |
---|
319 | ! hydrostatic pressure |
---|
320 | !-------------------------------------------- |
---|
321 | p = p0 * exp(- g / Rd * (inttau1 - inttau2 * inttermT)) |
---|
322 | |
---|
323 | END SUBROUTINE evaluate_pressure_temperature |
---|
324 | |
---|
325 | !----------------------------------------------------------------------- |
---|
326 | ! Calculate pointwise z and temperature given pressure |
---|
327 | !----------------------------------------------------------------------- |
---|
328 | SUBROUTINE evaluate_z_temperature(deep, X, lon, lat, p, z, t) |
---|
329 | |
---|
330 | INTEGER, INTENT(IN) :: deep ! Deep (1) or Shallow (0) test case |
---|
331 | |
---|
332 | REAL(8), INTENT(IN) :: & |
---|
333 | X, & ! Earth scaling ratio |
---|
334 | lon, & ! Longitude (radians) |
---|
335 | lat, & ! Latitude (radians) |
---|
336 | p ! Pressure (Pa) |
---|
337 | |
---|
338 | REAL(8), INTENT(OUT) :: & |
---|
339 | z, & ! Altitude (m) |
---|
340 | t ! Temperature (K) |
---|
341 | |
---|
342 | INTEGER :: ix |
---|
343 | |
---|
344 | REAL(8) :: z0, z1, z2 |
---|
345 | REAL(8) :: p0, p1, p2 |
---|
346 | |
---|
347 | z0 = 0.d0 |
---|
348 | z1 = 10000.d0 |
---|
349 | |
---|
350 | CALL evaluate_pressure_temperature(deep, X, lon, lat, z0, p0, t) |
---|
351 | CALL evaluate_pressure_temperature(deep, X, lon, lat, z1, p1, t) |
---|
352 | |
---|
353 | DO ix = 1, 100 |
---|
354 | z2 = z1 - (p1 - p) * (z1 - z0) / (p1 - p0) |
---|
355 | |
---|
356 | CALL evaluate_pressure_temperature(deep, X, lon, lat, z2, p2, t) |
---|
357 | |
---|
358 | IF (ABS((p2 - p)/p) .lt. 1.0d-13) THEN |
---|
359 | EXIT |
---|
360 | END IF |
---|
361 | |
---|
362 | z0 = z1 |
---|
363 | p0 = p1 |
---|
364 | |
---|
365 | z1 = z2 |
---|
366 | p1 = p2 |
---|
367 | END DO |
---|
368 | |
---|
369 | z = z2 |
---|
370 | |
---|
371 | CALL evaluate_pressure_temperature(deep, X, lon, lat, z, p0, t) |
---|
372 | |
---|
373 | END SUBROUTINE evaluate_z_temperature |
---|
374 | |
---|
375 | !----------------------------------------------------------------------- |
---|
376 | ! Exponential perturbation function |
---|
377 | !----------------------------------------------------------------------- |
---|
378 | REAL(8) FUNCTION evaluate_exponential(lon, lat, z) |
---|
379 | |
---|
380 | REAL(8), INTENT(IN) :: & |
---|
381 | lon, & ! Longitude (radians) |
---|
382 | lat, & ! Latitude (radians) |
---|
383 | z ! Altitude (meters) |
---|
384 | |
---|
385 | REAL(8) :: greatcircler, perttaper |
---|
386 | |
---|
387 | ! Great circle distance |
---|
388 | greatcircler = 1.d0 / pertexpr & |
---|
389 | * acos(sin(pertlat) * sin(lat) + cos(pertlat) * cos(lat) * cos(lon - pertlon)) |
---|
390 | |
---|
391 | ! Vertical tapering of stream function |
---|
392 | if (z < pertz) then |
---|
393 | perttaper = 1.d0 - 3.d0 * z**2 / pertz**2 + 2.d0 * z**3 / pertz**3 |
---|
394 | else |
---|
395 | perttaper = 0.d0 |
---|
396 | end if |
---|
397 | |
---|
398 | ! Zonal velocity perturbation |
---|
399 | if (greatcircler < 1.d0) then |
---|
400 | evaluate_exponential = pertup * perttaper * exp(- greatcircler**2) |
---|
401 | else |
---|
402 | evaluate_exponential = 0.d0 |
---|
403 | end if |
---|
404 | |
---|
405 | END FUNCTION evaluate_exponential |
---|
406 | |
---|
407 | !----------------------------------------------------------------------- |
---|
408 | ! Stream function perturbation function |
---|
409 | !----------------------------------------------------------------------- |
---|
410 | REAL(8) FUNCTION evaluate_streamfunction(lon, lat, z) |
---|
411 | |
---|
412 | REAL(8), INTENT(IN) :: & |
---|
413 | lon, & ! Longitude (radians) |
---|
414 | lat, & ! Latitude (radians) |
---|
415 | z ! Altitude (meters) |
---|
416 | |
---|
417 | REAL(8) :: greatcircler, perttaper, cospert |
---|
418 | |
---|
419 | ! Great circle distance |
---|
420 | greatcircler = 1.d0 / pertr & |
---|
421 | * acos(sin(pertlat) * sin(lat) + cos(pertlat) * cos(lat) * cos(lon - pertlon)) |
---|
422 | |
---|
423 | ! Vertical tapering of stream function |
---|
424 | if (z < pertz) then |
---|
425 | perttaper = 1.d0 - 3.d0 * z**2 / pertz**2 + 2.d0 * z**3 / pertz**3 |
---|
426 | else |
---|
427 | perttaper = 0.d0 |
---|
428 | end if |
---|
429 | |
---|
430 | ! Horizontal tapering of stream function |
---|
431 | if (greatcircler .lt. 1.d0) then |
---|
432 | cospert = cos(0.5d0 * pi * greatcircler) |
---|
433 | else |
---|
434 | cospert = 0.d0 |
---|
435 | end if |
---|
436 | |
---|
437 | evaluate_streamfunction = & |
---|
438 | (- pertu0 * pertr * perttaper * cospert**4) |
---|
439 | |
---|
440 | END FUNCTION evaluate_streamfunction |
---|
441 | |
---|
442 | END MODULE dcmip2016_baroclinic_wave_mod |
---|
443 | |
---|