FV3 Bundle
dcmip_initial_conditions_test_1_2_3_v5.f90
Go to the documentation of this file.
2 
3  !=======================================================================
4  !
5  ! Functions for setting up initial conditions for the dynamical core tests:
6  !
7  ! 11 - Deformational Advection Test
8  ! 12 - Hadley Cell Advection Test
9  ! 13 - Orography Advection Test
10  ! 20 - Impact of orography on a steady-state at rest
11  ! 21 and 22 - Non-Hydrostatic Mountain Waves Over A Schaer-Type Mountain without and with vertical wind shear
12  ! 31 - Non-Hydrostatic Gravity Waves
13  !
14  ! Given a point specified by:
15  ! lon longitude (radians)
16  ! lat latitude (radians)
17  ! p/z pressure/height
18  ! the functions will return:
19  ! u zonal wind (m s^-1)
20  ! v meridional wind (m s^-1)
21  ! w vertical velocity (m s^-1)
22  ! t temperature (K)
23  ! phis surface geopotential (m^2 s^-2)
24  ! ps surface pressure (Pa)
25  ! rho density (kj m^-3)
26  ! q specific humidity (kg/kg)
27  ! qi tracers (kg/kg)
28  ! p pressure if height based (Pa)
29  !
30  !
31  ! Authors: James Kent, Paul Ullrich, Christiane Jablonowski
32  ! (University of Michigan, dcmip@ucar.edu)
33  ! version 4
34  ! July/8/2012
35  !
36  ! Change log: (v3, June/8/2012, v4 July/8/2012, v5 July/20/2012)
37  !
38  ! v2: bug fixes in the tracer initialization for height-based models
39  ! v3: test 3-1: the density is now initialized with the unperturbed background temperature (not the perturbed temperature)
40  ! v3: constants converted to double precision
41  ! v4: modified tracers in test 1-1, now with cutoff altitudes. Outside of the vertical domain all tracers are set to 0
42  ! v4: modified cos-term in vertical velocity (now cos(2 pi t/tau)) in test 1-1, now completing one full up and down cycle
43  ! v4: added subroutine test1_advection_orography for test 1-3
44  ! v4: added subroutine test2_steady_state_mountain for test 2-0
45  ! v4: modified parameter list for tests 2-1 and 2-2 (routine test2_schaer_mountain): addition of parameters hybrid_eta, hyam, hybm
46  ! if the logical flag hybrid_eta is true then the pressure in pressure-based model with hybrid sigma-p (eta) coordinates is
47  ! computed internally. In that case the hybrid coefficients hyam and hybm need to be supplied via the parameter list,
48  ! otherwise they are not used.
49  ! v5: Change in test 11 - change to u and w, cutoff altitudes (introduced in v4) are removed again
50  ! v5: Change in test 12 - velocities multiplies by rho0/rho, different w0 and vertical location of the initial tracer
51  !
52  !
53  !=======================================================================
54 
55  IMPLICIT NONE
56 
57 !-----------------------------------------------------------------------
58 ! Physical Parameters
59 !-----------------------------------------------------------------------
60 
61  real(8), parameter :: a = 6371220.0d0, & ! Earth's Radius (m)
62  rd = 287.0d0, & ! Ideal gas const dry air (J kg^-1 K^1)
63  g = 9.80616d0, & ! Gravity (m s^2)
64  cp = 1004.5d0, & ! Specific heat capacity (J kg^-1 K^1)
65  pi = 4.d0*atan(1.d0) ! pi
66 
67 !-----------------------------------------------------------------------
68 ! Additional constants
69 !-----------------------------------------------------------------------
70 
71  real(8), parameter :: p0 = 100000.d0 ! reference pressure (Pa)
72 
73 
74 CONTAINS
75 
76 !==========================================================================================
77 ! TEST CASE 11 - PURE ADVECTION - 3D DEFORMATIONAL FLOW
78 !==========================================================================================
79 
80 ! The 3D deformational flow test is based on the deformational flow test of Nair and Lauritzen (JCP 2010),
81 ! with a prescribed vertical wind velocity which makes the test truly 3D. An unscaled planet (with scale parameter
82 ! X = 1) is selected.
83 
84 SUBROUTINE test1_advection_deformation (lon,lat,p,z,zcoords,u,v,w,t,phis,ps,rho,q,q1,q2,q3,q4)
85 
86 IMPLICIT NONE
87 !-----------------------------------------------------------------------
88 ! input/output params parameters at given location
89 !-----------------------------------------------------------------------
90 
91  real(8), intent(in) :: lon, & ! Longitude (radians)
92  lat, & ! Latitude (radians)
93  z ! Height (m)
94 
95  real(8), intent(inout) :: p ! Pressure (Pa)
96 
97  integer, intent(in) :: zcoords ! 0 or 1 see below
98 
99  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
100  v, & ! Meridional wind (m s^-1)
101  w, & ! Vertical Velocity (m s^-1)
102  t, & ! Temperature (K)
103  phis, & ! Surface Geopotential (m^2 s^-2)
104  ps, & ! Surface Pressure (Pa)
105  rho, & ! density (kg m^-3)
106  q, & ! Specific Humidity (kg/kg)
107  q1, & ! Tracer q1 (kg/kg)
108  q2, & ! Tracer q2 (kg/kg)
109  q3, & ! Tracer q3 (kg/kg)
110  q4 ! Tracer q4 (kg/kg)
111 
112  ! if zcoords = 1, then we use z and output p
113  ! if zcoords = 0, then we use p
114 
115 !-----------------------------------------------------------------------
116 ! test case parameters
117 !-----------------------------------------------------------------------
118  real(8), parameter :: tau = 12.d0 * 86400.d0, & ! period of motion 12 days
119  u0 = (2.d0*pi*a)/tau, & ! 2 pi a / 12 days
120  k0 = (10.d0*a)/tau, & ! Velocity Magnitude
121  omega0 = (23000.d0*pi)/tau, & ! Velocity Magnitude
122  t0 = 300.d0, & ! temperature
123  h = rd * t0 / g, & ! scale height
124  rr = 1.d0/2.d0, & ! horizontal half width divided by 'a'
125  zz = 1000.d0, & ! vertical half width
126  z0 = 5000.d0, & ! center point in z
127  lambda0 = 5.d0*pi/6.d0, & ! center point in longitudes
128  lambda1 = 7.d0*pi/6.d0, & ! center point in longitudes
129  phi0 = 0.d0, & ! center point in latitudes
130  phi1 = 0.d0
131 
132  real(8) :: height ! The height of the model levels
133  real(8) :: ptop ! Model top in p
134  real(8) :: sin_tmp, cos_tmp, sin_tmp2, cos_tmp2 ! Calculate great circle distances
135  real(8) :: d1, d2, r, r2, d3, d4 ! For tracer calculations
136  real(8) :: s, bs ! Shape function, and parameter
137  real(8) :: lonp ! Translational longitude, depends on time
138  real(8) :: ud ! Divergent part of u
139  real(8) :: time ! Initially set to zero seconds, needs
140  ! to be modified when used in dycore
141 
142 !-----------------------------------------------------------------------
143 ! HEIGHT AND PRESSURE
144 !-----------------------------------------------------------------------
145 
146  ! Height and pressure are aligned (p = p0 exp(-z/H))
147 
148  if (zcoords .eq. 1) then
149 
150  height = z
151  p = p0 * exp(-z/h)
152 
153  else
154 
155  height = h * log(p0/p)
156 
157  endif
158 
159  ! Model top in p
160 
161  ptop = p0*exp(-12000.d0/h)
162 
163 !-----------------------------------------------------------------------
164 ! THE VELOCITIES ARE TIME DEPENDENT AND THEREFORE MUST BE UPDATED
165 ! IN THE DYNAMICAL CORE
166 !-----------------------------------------------------------------------
167 
168  ! These are initial conditions hence time = 0
169 
170  time = 0.d0
171 
172  ! Translational longitude = longitude when time = 0
173 
174  lonp = lon - 2.d0*pi*time/tau
175 
176  ! Shape function
177 !********
178 ! change in version 5: shape function
179 !********
180  bs = 0.2
181  s = 1.0 + exp( (ptop-p0)/(bs*ptop) ) - exp( (p-p0)/(bs*ptop)) - exp( (ptop-p)/(bs*ptop))
182 
183  ! Zonal Velocity
184 !********
185 ! change in version 5: ud
186 !********
187 
188  ud = (omega0*a)/(bs*ptop) * cos(lonp) * (cos(lat)**2.0) * cos(2.0*pi*time/tau) * &
189  ( - exp( (p-p0)/(bs*ptop)) + exp( (ptop-p)/(bs*ptop)) )
190 
191  u = k0*sin(lonp)*sin(lonp)*sin(2.d0*lat)*cos(pi*time/tau) + u0*cos(lat) + ud
192 
193  ! Meridional Velocity
194 
195  v = k0*sin(2.d0*lonp)*cos(lat)*cos(pi*time/tau)
196 
197  ! Vertical Velocity - can be changed to vertical pressure velocity by
198  ! omega = -(g*p)/(Rd*T0)*w
199  !
200 !********
201 ! change in version 4: cos(2.0*pi*time/tau) is now used instead of cos(pi*time/tau)
202 !********
203 !********
204 ! change in version 5: shape function in w
205 !********
206  w = -((rd*t0)/(g*p))*omega0*sin(lonp)*cos(lat)*cos(2.0*pi*time/tau)*s
207 
208 !-----------------------------------------------------------------------
209 ! TEMPERATURE IS CONSTANT 300 K
210 !-----------------------------------------------------------------------
211 
212  t = t0
213 
214 !-----------------------------------------------------------------------
215 ! PHIS (surface geopotential)
216 !-----------------------------------------------------------------------
217 
218  phis = 0.d0
219 
220 !-----------------------------------------------------------------------
221 ! PS (surface pressure)
222 !-----------------------------------------------------------------------
223 
224  ps = p0
225 
226 !-----------------------------------------------------------------------
227 ! RHO (density)
228 !-----------------------------------------------------------------------
229 
230  rho = p/(rd*t)
231 
232 !-----------------------------------------------------------------------
233 ! initialize Q, set to zero
234 !-----------------------------------------------------------------------
235 
236  q = 0.d0
237 
238 !-----------------------------------------------------------------------
239 ! initialize tracers
240 !-----------------------------------------------------------------------
241 
242  ! Tracer 1 - Cosine Bells
243 
244  ! To calculate great circle distance
245 
246  sin_tmp = sin(lat) * sin(phi0)
247  cos_tmp = cos(lat) * cos(phi0)
248  sin_tmp2 = sin(lat) * sin(phi1)
249  cos_tmp2 = cos(lat) * cos(phi1)
250 
251  ! great circle distance without 'a'
252 
253  r = acos(sin_tmp + cos_tmp*cos(lon-lambda0))
254  r2 = acos(sin_tmp2 + cos_tmp2*cos(lon-lambda1))
255  d1 = min( 1.d0, (r/rr)**2 + ((height-z0)/zz)**2 )
256  d2 = min( 1.d0, (r2/rr)**2 + ((height-z0)/zz)**2 )
257 
258  q1 = 0.5d0 * (1.d0 + cos(pi*d1)) + 0.5d0 * (1.d0 + cos(pi*d2))
259 
260  ! Tracer 2 - Correlated Cosine Bells
261 
262  q2 = 0.9d0 - 0.8d0*q1**2
263 
264  ! Tracer 3 - Slotted Ellipse
265 
266  ! Make the ellipse
267 
268  if (d1 .le. rr) then
269  q3 = 1.d0
270  elseif (d2 .le. rr) then
271  q3 = 1.d0
272  else
273  q3 = 0.1d0
274  endif
275 
276  ! Put in the slot
277 
278  if (height .gt. z0 .and. abs(lat) .lt. 0.125d0) then
279 
280  q3 = 0.1d0
281 
282  endif
283 
284  ! Tracer 4: q4 is chosen so that, in combination with the other three tracer
285  ! fields with weight (3/10), the sum is equal to one
286 
287  q4 = 1.d0 - 0.3d0*(q1+q2+q3)
288 
289 !************
290 ! change in version 4: added cutoff altitudes, tracers are set to zero outside this region
291 ! prevents tracers from being trapped near the bottom and top of the domain
292 !************
293  ! use a cutoff altitude for
294  ! tracer 2 3 and 4
295  ! Set them to zero outside `buffer zone'
296 !************
297 ! change in version 5: change from v4 reversed, no cutoff altitudes due to continuity equation being satisfied
298 ! commented out below
299 !************
300 
301  !if (height .gt. (z0+1.25*ZZ) .or. height .lt. (z0-1.25*ZZ)) then
302 
303  ! q2 = 0.0
304  ! q3 = 0.0
305  ! q4 = 0.0
306 
307  !endif
308 
309 
310 
311 
312 END SUBROUTINE test1_advection_deformation
313 
314 
315 
316 
317 
318 !==========================================================================================
319 ! TEST CASE 12 - PURE ADVECTION - 3D HADLEY-LIKE FLOW
320 !==========================================================================================
321 
322 SUBROUTINE test1_advection_hadley (lon,lat,p,z,zcoords,u,v,w,t,phis,ps,rho,q,q1)
324 IMPLICIT NONE
325 !-----------------------------------------------------------------------
326 ! input/output params parameters at given location
327 !-----------------------------------------------------------------------
328 
329  real(8), intent(in) :: lon, & ! Longitude (radians)
330  lat, & ! Latitude (radians)
331  z ! Height (m)
332 
333  real(8), intent(inout) :: p ! Pressure (Pa)
334 
335  integer, intent(in) :: zcoords ! 0 or 1 see below
336 
337  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
338  v, & ! Meridional wind (m s^-1)
339  w, & ! Vertical Velocity (m s^-1)
340  t, & ! Temperature (K)
341  phis, & ! Surface Geopotential (m^2 s^-2)
342  ps, & ! Surface Pressure (Pa)
343  rho, & ! density (kg m^-3)
344  q, & ! Specific Humidity (kg/kg)
345  q1 ! Tracer q1 (kg/kg)
346 
347  ! if zcoords = 1, then we use z and output p
348  ! if zcoords = 0, then we use p
349 
350 !-----------------------------------------------------------------------
351 ! test case parameters
352 !-----------------------------------------------------------------------
353  real(8), parameter :: tau = 1.d0 * 86400.d0, & ! period of motion 1 day (in s)
354  u0 = 40.d0, & ! Zonal velocity magnitude (m/s)
355  w0 = 0.15d0, & ! Vertical velocity magnitude (m/s), changed in v5
356  t0 = 300.d0, & ! temperature (K)
357  h = rd * t0 / g, & ! scale height
358  k = 5.d0, & ! number of Hadley-like cells
359  z1 = 2000.d0, & ! position of lower tracer bound (m), changed in v5
360  z2 = 5000.d0, & ! position of upper tracer bound (m), changed in v5
361  z0 = 0.5d0*(z1+z2), & ! midpoint (m)
362  ztop = 12000.d0 ! model top (m)
363 
364  real(8) :: rho0 ! reference density at z=0 m
365  real(8) :: height ! Model level heights
366  real(8) :: time ! Initially set to zero seconds, needs
367  ! to be modified when used in dycore
368 
369 !-----------------------------------------------------------------------
370 ! HEIGHT AND PRESSURE
371 !-----------------------------------------------------------------------
372 
373  ! Height and pressure are aligned (p = p0 exp(-z/H))
374 
375  if (zcoords .eq. 1) then
376 
377  height = z
378  p = p0 * exp(-z/h)
379 
380  else
381 
382  height = h * log(p0/p)
383 
384  endif
385 
386 !-----------------------------------------------------------------------
387 ! TEMPERATURE IS CONSTANT 300 K
388 !-----------------------------------------------------------------------
389 
390  t = t0
391 
392 !-----------------------------------------------------------------------
393 ! PHIS (surface geopotential)
394 !-----------------------------------------------------------------------
395 
396  phis = 0.d0
397 
398 !-----------------------------------------------------------------------
399 ! PS (surface pressure)
400 !-----------------------------------------------------------------------
401 
402  ps = p0
403 
404 !-----------------------------------------------------------------------
405 ! RHO (density)
406 !-----------------------------------------------------------------------
407 
408  rho = p/(rd*t)
409  rho0 = p0/(rd*t)
410 
411 !-----------------------------------------------------------------------
412 ! THE VELOCITIES ARE TIME DEPENDENT AND THEREFORE MUST BE UPDATED
413 ! IN THE DYNAMICAL CORE
414 !-----------------------------------------------------------------------
415 
416  time = 0.d0
417 
418  ! Zonal Velocity
419 
420  u = u0*cos(lat)
421 
422  ! Meridional Velocity
423 
424 !************
425 ! change in version 5: multiply v and w by rho0/rho
426 !************
427 
428  v = -(rho0/rho) * (a*w0*pi)/(k*ztop) *cos(lat)*sin(k*lat)*cos(pi*height/ztop)*cos(pi*time/tau)
429 
430  ! Vertical Velocity - can be changed to vertical pressure velocity by
431  ! omega = -g*rho*w
432 
433  w = (rho0/rho) *(w0/k)*(-2.d0*sin(k*lat)*sin(lat) + k*cos(lat)*cos(k*lat)) &
434  *sin(pi*height/ztop)*cos(pi*time/tau)
435 
436 
437 !-----------------------------------------------------------------------
438 ! initialize Q, set to zero
439 !-----------------------------------------------------------------------
440 
441  q = 0.d0
442 
443 !-----------------------------------------------------------------------
444 ! initialize tracers
445 !-----------------------------------------------------------------------
446 
447  ! Tracer 1 - Layer
448 
449  if (height .lt. z2 .and. height .gt. z1) then
450 
451  q1 = 0.5d0 * (1.d0 + cos( 2.d0*pi*(height-z0)/(z2-z1) ) )
452 
453  else
454 
455  q1 = 0.d0
456 
457  endif
458 
459 END SUBROUTINE test1_advection_hadley
460 
461 
462 
463 !==========================================================================================
464 ! TEST CASE 13 - HORIZONTAL ADVECTION OF THIN CLOUD-LIKE TRACERS IN THE PRESENCE OF OROGRAPHY
465 !==========================================================================================
466 
467 SUBROUTINE test1_advection_orography (lon,lat,p,z,zcoords,cfv,hybrid_eta,hyam,hybm,gc,u,v,w,t,phis,ps,rho,q,q1,q2,q3,q4)
469 IMPLICIT NONE
470 !-----------------------------------------------------------------------
471 ! input/output params parameters at given location
472 !-----------------------------------------------------------------------
473 
474  real(8), intent(in) :: lon, & ! Longitude (radians)
475  lat, & ! Latitude (radians)
476  z, & ! Height (m)
477  hyam, & ! A coefficient for hybrid-eta coordinate, at model level midpoint
478  hybm, & ! B coefficient for hybrid-eta coordinate, at model level midpoint
479  gc ! bar{z} for Gal-Chen coordinate
480 
481  logical, intent(in) :: hybrid_eta ! flag to indicate whether the hybrid sigma-p (eta) coordinate is used
482  ! if set to .true., then the pressure will be computed via the
483  ! hybrid coefficients hyam and hybm, they need to be initialized
484  ! if set to .false. (for pressure-based models): the pressure is already pre-computed
485  ! and is an input value for this routine
486  ! for height-based models: pressure will always be computed based on the height and
487  ! hybrid_eta is not used
488 
489  ! Note that we only use hyam and hybm for the hybrid-eta coordiantes, and we only use
490  ! gc for the Gal-Chen coordinates. If not required then they become dummy variables
491 
492  real(8), intent(inout) :: p ! Pressure (Pa)
493 
494  integer, intent(in) :: zcoords ! 0 or 1 see below
495  integer, intent(in) :: cfv ! 0, 1 or 2 see below
496 
497  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
498  v, & ! Meridional wind (m s^-1)
499  w, & ! Vertical Velocity (m s^-1)
500  t, & ! Temperature (K)
501  phis, & ! Surface Geopotential (m^2 s^-2)
502  ps, & ! Surface Pressure (Pa)
503  rho, & ! density (kg m^-3)
504  q, & ! Specific Humidity (kg/kg)
505  q1, & ! Tracer q1 (kg/kg)
506  q2, & ! Tracer q2 (kg/kg)
507  q3, & ! Tracer q3 (kg/kg)
508  q4 ! Tracer q4 (kg/kg)
509 
510  ! if zcoords = 1, then we use z and output p
511  ! if zcoords = 0, then we use p
512 
513  ! if cfv = 0 we assume that our horizontal velocities are not coordinate following
514  ! if cfv = 1 then our velocities follow hybrid eta coordinates and we need to specify w
515  ! if cfv = 2 then our velocities follow Gal-Chen coordinates and we need to specify w
516 
517  ! In hybrid-eta coords: p = hyam p0 + hybm ps
518  ! In Gal-Chen coords: z = zs + (gc/ztop)*(ztop - zs)
519 
520  ! if other orography-following coordinates are used, the w wind needs to be newly derived for them
521 
522 !-----------------------------------------------------------------------
523 ! test case parameters
524 !-----------------------------------------------------------------------
525  real(8), parameter :: tau = 12.d0 * 86400.d0, & ! period of motion 12 days (s)
526  u0 = 2.d0*pi*a/tau, & ! Velocity Magnitude (m/s)
527  t0 = 300.d0, & ! temperature (K)
528  h = rd * t0 / g, & ! scale height (m)
529  alpha = pi/6.d0, & ! rotation angle (radians), 30 degrees
530  lambdam = 3.d0*pi/2.d0, & ! mountain longitude center point (radians)
531  phim = 0.d0, & ! mountain latitude center point (radians)
532  h0 = 2000.d0, & ! peak height of the mountain range (m)
533  rm = 3.d0*pi/4.d0, & ! mountain radius (radians)
534  zetam = pi/16.d0, & ! mountain oscillation half-width (radians)
535  lambdap = pi/2.d0, & ! cloud-like tracer longitude center point (radians)
536  phip = 0.d0, & ! cloud-like tracer latitude center point (radians)
537  rp = pi/4.d0, & ! cloud-like tracer radius (radians)
538  zp1 = 3050.d0, & ! midpoint of first (lowermost) tracer (m)
539  zp2 = 5050.d0, & ! midpoint of second tracer (m)
540  zp3 = 8200.d0, & ! midpoint of third (topmost) tracer (m)
541  dzp1 = 1000.d0, & ! thickness of first (lowermost) tracer (m)
542  dzp2 = 1000.d0, & ! thickness of second tracer (m)
543  dzp3 = 400.d0, & ! thickness of third (topmost) tracer (m)
544  ztop = 12000.d0 ! model top (m)
545 
546  real(8) :: height ! Model level heights (m)
547  real(8) :: r ! Great circle distance (radians)
548  real(8) :: rz ! height differences
549  real(8) :: zs ! Surface elevation (m)
550 
551 
552 
553 !-----------------------------------------------------------------------
554 ! PHIS (surface geopotential)
555 !-----------------------------------------------------------------------
556 
557  r = acos( sin(phim)*sin(lat) + cos(phim)*cos(lat)*cos(lon - lambdam) )
558 
559  if (r .lt. rm) then
560 
561  zs = (h0/2.d0)*(1.d0+cos(pi*r/rm))*cos(pi*r/zetam)**2.d0
562 
563  else
564 
565  zs = 0.d0
566 
567  endif
568 
569  phis = g*zs
570 
571 !-----------------------------------------------------------------------
572 ! PS (surface pressure)
573 !-----------------------------------------------------------------------
574 
575  ps = p0 * exp(-zs/h)
576 
577 
578 !-----------------------------------------------------------------------
579 ! HEIGHT AND PRESSURE
580 !-----------------------------------------------------------------------
581 
582  ! Height and pressure are aligned (p = p0 exp(-z/H))
583 
584  if (zcoords .eq. 1) then
585 
586  height = z
587  p = p0 * exp(-z/h)
588 
589  else
590 
591  if (hybrid_eta) p = hyam*p0 + hybm*ps ! compute the pressure based on the surface pressure and hybrid coefficients
592  height = h * log(p0/p)
593 
594  endif
595 
596 !-----------------------------------------------------------------------
597 ! THE VELOCITIES ARE TIME INDEPENDENT
598 !-----------------------------------------------------------------------
599 
600  ! Zonal Velocity
601 
602  u = u0*(cos(lat)*cos(alpha)+sin(lat)*cos(lon)*sin(alpha))
603 
604  ! Meridional Velocity
605 
606  v = -u0*(sin(lon)*sin(alpha))
607 
608 !-----------------------------------------------------------------------
609 ! TEMPERATURE IS CONSTANT 300 K
610 !-----------------------------------------------------------------------
611 
612  t = t0
613 
614 !-----------------------------------------------------------------------
615 ! RHO (density)
616 !-----------------------------------------------------------------------
617 
618  rho = p/(rd*t)
619 
620 !-----------------------------------------------------------------------
621 ! initialize Q, set to zero
622 !-----------------------------------------------------------------------
623 
624  q = 0.d0
625 
626 !-----------------------------------------------------------------------
627 ! VERTICAL VELOCITY IS TIME INDEPENDENT
628 !-----------------------------------------------------------------------
629 
630  ! Vertical Velocity - can be changed to vertical pressure velocity by
631  ! omega = -(g*p)/(Rd*T0)*w
632 
633  ! NOTE that if orography-following coordinates are used then the vertical
634  ! velocity needs to be translated into the new coordinate system due to
635  ! the variation of the height along coordinate surfaces
636  ! See section 1.3 and the appendix of the test case document
637 
638  if (cfv .eq. 0) then
639 
640  ! if the horizontal velocities do not follow the vertical coordinate
641 
642  w = 0.d0
643 
644  elseif (cfv .eq. 1) then
645 
646  ! if the horizontal velocities follow hybrid eta coordinates then
647  ! the perceived vertical velocity is
648 
650 
651  elseif (cfv .eq. 2) then
652 
653  ! if the horizontal velocities follow Gal Chen coordinates then
654  ! the perceived vertical velocity is
655 
657 
658 ! else
659 ! compute your own vertical velocity if other orography-following
660 ! vertical coordinate is used
661 !
662  endif
663 
664 
665 
666 !-----------------------------------------------------------------------
667 ! initialize tracers
668 !-----------------------------------------------------------------------
669 
670  ! Tracer 1 - Cloud Layer
671 
672  r = acos( sin(phip)*sin(lat) + cos(phip)*cos(lat)*cos(lon - lambdap) )
673 
674  rz = abs(height - zp1)
675 
676  if (rz .lt. 0.5d0*dzp1 .and. r .lt. rp) then
677 
678  q1 = 0.25d0*(1.d0+cos(2.d0*pi*rz/dzp1))*(1.d0+cos(pi*r/rp))
679 
680  else
681 
682  q1 = 0.d0
683 
684  endif
685 
686  rz = abs(height - zp2)
687 
688  if (rz .lt. 0.5d0*dzp2 .and. r .lt. rp) then
689 
690  q2 = 0.25d0*(1.d0+cos(2.d0*pi*rz/dzp2))*(1.d0+cos(pi*r/rp))
691 
692  else
693 
694  q2 = 0.d0
695 
696  endif
697 
698  rz = abs(height - zp3)
699 
700  if (rz .lt. 0.5d0*dzp3 .and. r .lt. rp) then
701 
702  q3 = 1.d0
703 
704  else
705 
706  q3 = 0.d0
707 
708  endif
709 
710  q4 = q1 + q2 + q3
711 
712 
713  CONTAINS
714 
715  !-----------------------------------------------------------------------
716  ! SUBROUTINE TO CALCULATE THE PERCEIVED VERTICAL VELOCITY
717  ! UNDER HYBRID-ETA COORDINATES
718  !-----------------------------------------------------------------------
719 
721  IMPLICIT NONE
722  real(8), intent(out) :: w
723 
724  real(8) :: press, & ! hyam *p0 + hybm *ps
725  r, & ! Great Circle Distance
726  dzsdx, & ! Part of surface height derivative
727  dzsdlambda, & ! Derivative of zs w.r.t lambda
728  dzsdphi, & ! Derivative of zs w.r.t phi
729  dzdlambda, & ! Derivative of z w.r.t lambda
730  dzdphi, & ! Derivative of z w.r.t phi
731  dpsdlambda, & ! Derivative of ps w.r.t lambda
732  dpsdphi ! Derivative of ps w.r.t phi
733 
734  ! Calculate pressure and great circle distance to mountain center
735 
736  press = hyam*p0 + hybm*ps
737 
738  r = acos( sin(phim)*sin(lat) + cos(phim)*cos(lat)*cos(lon - lambdam) )
739 
740  ! Derivatives of surface height
741 
742  if (r .lt. rm) then
743  dzsdx = -h0*pi/(2.d0*rm)*sin(pi*r/rm)*cos(pi*r/zetam)**2 - &
744  (h0*pi/zetam)*(1.d0+cos(pi*r/rm))*cos(pi*r/zetam)*sin(pi*r/zetam)
745  else
746  dzsdx = 0.d0
747  endif
748 
749  ! Prevent division by zero
750 
751  if (1.d0-cos(r)**2 .gt. 0.d0) then
752  dzsdlambda = dzsdx * (cos(phim)*cos(lat)*sin(lon-lambdam)) &
753  /sqrt(1.d0-cos(r)**2)
754  dzsdphi = dzsdx * (-sin(phim)*cos(lat) + cos(phim)*sin(lat)*cos(lon-lambdam)) &
755  /sqrt(1.d0-cos(r)**2)
756  else
757  dzsdlambda = 0.d0
758  dzsdphi = 0.d0
759  endif
760 
761  ! Derivatives of surface pressure
762 
763  dpsdlambda = -(g*p0/(rd*t0))*exp(-g*zs/(rd*t0))*dzsdlambda
764  dpsdphi = -(g*p0/(rd*t0))*exp(-g*zs/(rd*t0))*dzsdphi
765 
766  ! Derivatives of coordinate height
767 
768  dzdlambda = -(rd*t0/(g*press))*hybm*dpsdlambda
769  dzdphi = -(rd*t0/(g*press))*hybm*dpsdphi
770 
771  ! Prevent division by zero
772 
773  if (abs(lat) .lt. pi/2.d0) then
774  w = - (u/(a*cos(lat)))*dzdlambda - (v/a)*dzdphi
775  else
776  w = 0.d0
777  endif
778 
780 
781 
782 
783  !-----------------------------------------------------------------------
784  ! SUBROUTINE TO CALCULATE THE PERCEIVED VERTICAL VELOCITY
785  ! UNDER GAL-CHEN COORDINATES
786  !-----------------------------------------------------------------------
787 
789  IMPLICIT NONE
790  real(8), intent(out) :: w
791 
792  real(8) :: r, & ! Great Circle Distance
793  dzsdx, & ! Part of surface height derivative
794  dzsdlambda, & ! Derivative of zs w.r.t lambda
795  dzsdphi, & ! Derivative of zs w.r.t phi
796  dzdlambda, & ! Derivative of z w.r.t lambda
797  dzdphi ! Derivative of z w.r.t phi
798 
799  ! Calculate great circle distance to mountain center
800 
801  r = acos( sin(phim)*sin(lat) + cos(phim)*cos(lat)*cos(lon - lambdam) )
802 
803  ! Derivatives of surface height
804 
805  if (r .lt. rm) then
806  dzsdx = -h0*pi/(2.d0*rm)*sin(pi*r/rm)*cos(pi*r/zetam)**2 - &
807  (h0*pi/zetam)*(1.d0+cos(pi*r/rm))*cos(pi*r/zetam)*sin(pi*r/zetam)
808  else
809  dzsdx = 0.d0
810  endif
811 
812  ! Prevent division by zero
813 
814  if (1.d0-cos(r)**2 .gt. 0.d0) then
815  dzsdlambda = dzsdx * (cos(phim)*cos(lat)*sin(lon-lambdam)) &
816  /sqrt(1.d0-cos(r)**2)
817  dzsdphi = dzsdx * (-sin(phim)*cos(lat) + cos(phim)*sin(lat)*cos(lon-lambdam)) &
818  /sqrt(1.d0-cos(r)**2)
819  else
820  dzsdlambda = 0.d0
821  dzsdphi = 0.d0
822  endif
823 
824  ! Derivatives of coordinate height
825 
826  dzdlambda = (1.d0-gc/ztop)*dzsdlambda
827  dzdphi = (1.d0-gc/ztop)*dzsdphi
828 
829  ! Prevent division by zero
830 
831  if (abs(lat) .lt. pi/2.d0) then
832  w = - (u/(a*cos(lat)))*dzdlambda - (v/a)*dzdphi
833  else
834  w = 0.d0
835  endif
836 
838 
839 END SUBROUTINE test1_advection_orography
840 
841 
842 
843 
844 !==========================================================================================
845 ! TEST CASE 2X - IMPACT OF OROGRAPHY ON A NON-ROTATING PLANET
846 !==========================================================================================
847 ! The tests in section 2-x examine the impact of 3D Schaer-like circular mountain profiles on an
848 ! atmosphere at rest (2-0), and on flow fields with wind shear (2-1) and without vertical wind shear (2-2).
849 ! A non-rotating planet is used for all configurations. Test 2-0 is conducted on an unscaled regular-size
850 ! planet and primarily examines the accuracy of the pressure gradient calculation in a steady-state
851 ! hydrostatically-balanced atmosphere at rest. This test is especially appealing for models with
852 ! orography-following vertical coordinates. It increases the complexity of test 1-3, that investigated
853 ! the impact of the same Schaer-type orographic profile on the accuracy of purely-horizontal passive
854 ! tracer advection.
855 !
856 ! Tests 2-1 and 2-2 increase the complexity even further since non-zero flow fields are now prescribed
857 ! with and without vertical wind shear. In order to trigger non-hydrostatic responses the two tests are
858 ! conducted on a reduced-size planet with reduction factor $X=500$ which makes the horizontal and
859 ! vertical grid spacing comparable. This test clearly discriminates between non-hydrostatic and hydrostatic
860 ! models since the expected response is in the non-hydrostatic regime. Therefore, the flow response is
861 ! captured differently by hydrostatic models.
862 
863 
864 
865 
866 !=========================================================================
867 ! Test 2-0: Steady-State Atmosphere at Rest in the Presence of Orography
868 !=========================================================================
869 SUBROUTINE test2_steady_state_mountain (lon,lat,p,z,zcoords,hybrid_eta,hyam,hybm,u,v,w,t,phis,ps,rho,q)
871 IMPLICIT NONE
872 !-----------------------------------------------------------------------
873 ! input/output params parameters at given location
874 !-----------------------------------------------------------------------
875 
876  real(8), intent(in) :: lon, & ! Longitude (radians)
877  lat, & ! Latitude (radians)
878  z, & ! Height (m)
879  hyam, & ! A coefficient for hybrid-eta coordinate, at model level midpoint
880  hybm ! B coefficient for hybrid-eta coordinate, at model level midpoint
881 
882  logical, intent(in) :: hybrid_eta ! flag to indicate whether the hybrid sigma-p (eta) coordinate is used
883  ! if set to .true., then the pressure will be computed via the
884  ! hybrid coefficients hyam and hybm, they need to be initialized
885  ! if set to .false. (for pressure-based models): the pressure is already pre-computed
886  ! and is an input value for this routine
887  ! for height-based models: pressure will always be computed based on the height and
888  ! hybrid_eta is not used
889 
890  real(8), intent(inout) :: p ! Pressure (Pa)
891 
892  integer, intent(in) :: zcoords ! 0 or 1 see below
893 
894  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
895  v, & ! Meridional wind (m s^-1)
896  w, & ! Vertical Velocity (m s^-1)
897  t, & ! Temperature (K)
898  phis, & ! Surface Geopotential (m^2 s^-2)
899  ps, & ! Surface Pressure (Pa)
900  rho, & ! density (kg m^-3)
901  q ! Specific Humidity (kg/kg)
902 
903  ! if zcoords = 1, then we use z and output p
904  ! if zcoords = 0, then we compute or use p
905  !
906  ! In hybrid-eta coords: p = hyam p0 + hybm ps
907  !
908  ! The grid-point based initial data are computed in this routine.
909 
910 !-----------------------------------------------------------------------
911 ! test case parameters
912 !-----------------------------------------------------------------------
913  real(8), parameter :: T0 = 300.d0, & ! temperature (K)
914  gamma = 0.0065d0, & ! temperature lapse rate (K/m)
915  lambdam = 3.d0*pi/2.d0, & ! mountain longitude center point (radians)
916  phim = 0.d0, & ! mountain latitude center point (radians)
917  h0 = 2000.d0, & ! peak height of the mountain range (m)
918  rm = 3.d0*pi/4.d0, & ! mountain radius (radians)
919  zetam = pi/16.d0, & ! mountain oscillation half-width (radians)
920  ztop = 12000.d0 ! model top (m)
921 
922  real(8) :: height ! Model level heights (m)
923  real(8) :: r ! Great circle distance (radians)
924  real(8) :: zs ! Surface elevation (m)
925  real(8) :: exponent ! exponent: g/(Rd * gamma)
926  real(8) :: exponent_rev ! reversed exponent
927 
928 
929 !-----------------------------------------------------------------------
930 ! compute exponents
931 !-----------------------------------------------------------------------
932  exponent = g/(rd*gamma)
933  exponent_rev = 1.d0/exponent
934 
935 !-----------------------------------------------------------------------
936 ! PHIS (surface geopotential)
937 !-----------------------------------------------------------------------
938 
939  r = acos( sin(phim)*sin(lat) + cos(phim)*cos(lat)*cos(lon - lambdam) )
940 
941  if (r .lt. rm) then
942 
943  zs = (h0/2.d0)*(1.d0+cos(pi*r/rm))*cos(pi*r/zetam)**2.d0 ! mountain height
944 
945  else
946 
947  zs = 0.d0
948 
949  endif
950 
951  phis = g*zs
952 
953 !-----------------------------------------------------------------------
954 ! PS (surface pressure)
955 !-----------------------------------------------------------------------
956 
957  ps = p0 * (1.d0 - gamma/t0*zs)**exponent
958 
959 
960 !-----------------------------------------------------------------------
961 ! HEIGHT AND PRESSURE
962 !-----------------------------------------------------------------------
963 
964  ! Height and pressure are aligned (p = p0 * (1.d0 - gamma/T0*z)**exponent)
965 
966  if (zcoords .eq. 1) then
967 
968  height = z
969  p = p0 * (1.d0 - gamma/t0*z)**exponent
970 
971  else
972 
973  if (hybrid_eta) p = hyam*p0 + hybm*ps ! compute the pressure based on the surface pressure and hybrid coefficients
974  height = t0/gamma * (1.d0 - (p/p0)**exponent_rev) ! compute the height at this pressure
975 
976  endif
977 
978 !-----------------------------------------------------------------------
979 ! THE VELOCITIES ARE ZERO (STATE AT REST)
980 !-----------------------------------------------------------------------
981 
982  ! Zonal Velocity
983 
984  u = 0.d0
985 
986  ! Meridional Velocity
987 
988  v = 0.d0
989 
990  ! Vertical Velocity
991 
992  w = 0.d0
993 
994 !-----------------------------------------------------------------------
995 ! TEMPERATURE WITH CONSTANT LAPSE RATE
996 !-----------------------------------------------------------------------
997 
998  t = t0 - gamma*height
999 
1000 !-----------------------------------------------------------------------
1001 ! RHO (density)
1002 !-----------------------------------------------------------------------
1003 
1004  rho = p/(rd*t)
1005 
1006 !-----------------------------------------------------------------------
1007 ! initialize Q, set to zero
1008 !-----------------------------------------------------------------------
1009 
1010  q = 0.d0
1011 
1012 END SUBROUTINE test2_steady_state_mountain
1013 
1014 
1015 
1016 !=====================================================================================
1017 ! Tests 2-1 and 2-2: Non-hydrostatic Mountain Waves over a Schaer-type Mountain
1018 !=====================================================================================
1019 
1020 SUBROUTINE test2_schaer_mountain (lon,lat,p,z,zcoords,hybrid_eta,hyam,hybm,shear,u,v,w,t,phis,ps,rho,q)
1022 IMPLICIT NONE
1023 !-----------------------------------------------------------------------
1024 ! input/output params parameters at given location
1025 !-----------------------------------------------------------------------
1026 
1027  real(8), intent(in) :: lon, & ! Longitude (radians)
1028  lat, & ! Latitude (radians)
1029  z, & ! Height (m)
1030  hyam, & ! A coefficient for hybrid-eta coordinate, at model level midpoint
1031  hybm ! B coefficient for hybrid-eta coordinate, at model level midpoint
1032 
1033  logical, intent(in) :: hybrid_eta ! flag to indicate whether the hybrid sigma-p (eta) coordinate is used
1034  ! if set to .true., then the pressure will be computed via the
1035  ! hybrid coefficients hyam and hybm, they need to be initialized
1036  ! if set to .false. (for pressure-based models): the pressure is already pre-computed
1037  ! and is an input value for this routine
1038  ! for height-based models: pressure will always be computed based on the height and
1039  ! hybrid_eta is not used
1040 
1041  real(8), intent(inout) :: p ! Pressure (Pa)
1042 
1043 
1044  integer, intent(in) :: zcoords, & ! 0 or 1 see below
1045  shear ! 0 or 1 see below
1046 
1047  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
1048  v, & ! Meridional wind (m s^-1)
1049  w, & ! Vertical Velocity (m s^-1)
1050  t, & ! Temperature (K)
1051  phis, & ! Surface Geopotential (m^2 s^-2)
1052  ps, & ! Surface Pressure (Pa)
1053  rho, & ! density (kg m^-3)
1054  q ! Specific Humidity (kg/kg)
1055 
1056  ! if zcoords = 1, then we use z and output p
1057  ! if zcoords = 0, then we either compute or use p
1058 
1059  ! if shear = 1, then we use shear flow
1060  ! if shear = 0, then we use constant u
1061 
1062 !-----------------------------------------------------------------------
1063 ! test case parameters
1064 !-----------------------------------------------------------------------
1065  real(8), parameter :: X = 500.d0, & ! Reduced Earth reduction factor
1066  om = 0.d0, & ! Rotation Rate of Earth
1067  as = a/x, & ! New Radius of small Earth
1068  ueq = 20.d0, & ! Reference Velocity
1069  teq = 300.d0, & ! Temperature at Equator
1070  peq = 100000.d0, & ! Reference PS at Equator
1071  ztop = 30000.d0, & ! Model Top
1072  lambdac = pi/4.d0, & ! Lon of Schar Mountain Center
1073  phic = 0.d0, & ! Lat of Schar Mountain Center
1074  h0 = 250.d0, & ! Height of Mountain
1075  d = 5000.d0, & ! Mountain Half-Width
1076  xi = 4000.d0, & ! Mountain Wavelength
1077  cs = 0.00025d0 ! Wind Shear (shear=1)
1078 
1079  real(8) :: height ! Model level heights
1080  real(8) :: sin_tmp, cos_tmp ! Calculation of great circle distance
1081  real(8) :: r ! Great circle distance
1082  real(8) :: zs ! Surface height
1083  real(8) :: c ! Shear
1084 
1085 !-----------------------------------------------------------------------
1086 ! PHIS (surface geopotential)
1087 !-----------------------------------------------------------------------
1088 
1089  sin_tmp = sin(lat) * sin(phic)
1090  cos_tmp = cos(lat) * cos(phic)
1091 
1092  ! great circle distance with 'a/X'
1093 
1094  r = as * acos(sin_tmp + cos_tmp*cos(lon-lambdac))
1095  zs = h0 * exp(-(r**2)/(d**2))*(cos(pi*r/xi)**2)
1096  phis = g*zs
1097 
1098 !-----------------------------------------------------------------------
1099 ! SHEAR FLOW OR CONSTANT FLOW
1100 !-----------------------------------------------------------------------
1101 
1102  if (shear .eq. 1) then
1103 
1104  c = cs
1105 
1106  else
1107 
1108  c = 0.d0
1109 
1110  endif
1111 
1112 !-----------------------------------------------------------------------
1113 ! TEMPERATURE
1114 !-----------------------------------------------------------------------
1115 
1116  t = teq *(1.d0 - (c*ueq*ueq/(g))*(sin(lat)**2) )
1117 
1118 !-----------------------------------------------------------------------
1119 ! PS (surface pressure)
1120 !-----------------------------------------------------------------------
1121 
1122  ps = peq*exp( -(ueq*ueq/(2.d0*rd*teq))*(sin(lat)**2) - phis/(rd*t) )
1123 
1124 !-----------------------------------------------------------------------
1125 ! HEIGHT AND PRESSURE
1126 !-----------------------------------------------------------------------
1127 
1128  if (zcoords .eq. 1) then
1129 
1130  height = z
1131  p = peq*exp( -(ueq*ueq/(2.d0*rd*teq))*(sin(lat)**2) - g*height/(rd*t) )
1132 
1133  else
1134  if (hybrid_eta) p = hyam*p0 + hybm*ps ! compute the pressure based on the surface pressure and hybrid coefficients
1135  height = (rd*t/(g))*log(peq/p) - (t*ueq*ueq/(2.d0*teq*g))*(sin(lat)**2)
1136 
1137  endif
1138 
1139 !-----------------------------------------------------------------------
1140 ! THE VELOCITIES
1141 !-----------------------------------------------------------------------
1142 
1143  ! Zonal Velocity
1144 
1145  u = ueq * cos(lat) * sqrt( (2.d0*teq/(t))*c*height + t/(teq) )
1146 
1147  ! Meridional Velocity
1148 
1149  v = 0.d0
1150 
1151  ! Vertical Velocity = Vertical Pressure Velocity = 0
1152 
1153  w = 0.d0
1154 
1155 !-----------------------------------------------------------------------
1156 ! RHO (density)
1157 !-----------------------------------------------------------------------
1158 
1159  rho = p/(rd*t)
1160 
1161 !-----------------------------------------------------------------------
1162 ! initialize Q, set to zero
1163 !-----------------------------------------------------------------------
1164 
1165  q = 0.d0
1166 
1167 END SUBROUTINE test2_schaer_mountain
1168 
1169 
1170 
1171 
1172 
1173 
1174 
1175 
1176 
1177 
1178 !==========================================================================================
1179 ! TEST CASE 3 - GRAVITY WAVES
1180 !==========================================================================================
1181 
1182 ! The non-hydrostatic gravity wave test examines the response of models to short time-scale wavemotion
1183 ! triggered by a localized perturbation. The formulation presented in this document is new,
1184 ! but is based on previous approaches by Skamarock et al. (JAS 1994), Tomita and Satoh (FDR 2004), and
1185 ! Jablonowski et al. (NCAR Tech Report 2008)
1186 
1187 
1188 !==========
1189 ! Test 3-1
1190 !==========
1191 SUBROUTINE test3_gravity_wave (lon,lat,p,z,zcoords,u,v,w,t,phis,ps,rho,q)
1193 IMPLICIT NONE
1194 !-----------------------------------------------------------------------
1195 ! input/output params parameters at given location
1196 !-----------------------------------------------------------------------
1197 
1198  real(8), intent(in) :: lon, & ! Longitude (radians)
1199  lat, & ! Latitude (radians)
1200  z ! Height (m)
1201 
1202  real(8), intent(inout) :: p ! Pressure (Pa)
1203 
1204 
1205  integer, intent(in) :: zcoords ! 0 or 1 see below
1206 
1207  real(8), intent(out) :: u, & ! Zonal wind (m s^-1)
1208  v, & ! Meridional wind (m s^-1)
1209  w, & ! Vertical Velocity (m s^-1)
1210  t, & ! Temperature (K)
1211  phis, & ! Surface Geopotential (m^2 s^-2)
1212  ps, & ! Surface Pressure (Pa)
1213  rho, & ! density (kg m^-3)
1214  q ! Specific Humidity (kg/kg)
1215 
1216  ! if zcoords = 1, then we use z and output z
1217  ! if zcoords = 0, then we use p
1218 
1219 !-----------------------------------------------------------------------
1220 ! test case parameters
1221 !-----------------------------------------------------------------------
1222  real(8), parameter :: X = 125.d0, & ! Reduced Earth reduction factor
1223  om = 0.d0, & ! Rotation Rate of Earth
1224  as = a/x, & ! New Radius of small Earth
1225  u0 = 20.d0, & ! Reference Velocity
1226  teq = 300.d0, & ! Temperature at Equator
1227  peq = 100000.d0, & ! Reference PS at Equator
1228  ztop = 10000.d0, & ! Model Top
1229  lambdac = 2.d0*pi/3.d0, & ! Lon of Pert Center
1230  d = 5000.d0, & ! Width for Pert
1231  phic = 0.d0, & ! Lat of Pert Center
1232  delta_theta = 1.d0, & ! Max Amplitude of Pert
1233  lz = 20000.d0, & ! Vertical Wavelength of Pert
1234  n = 0.01d0, & ! Brunt-Vaisala frequency
1235  n2 = n*n, & ! Brunt-Vaisala frequency Squared
1236  bigg = (g*g)/(n2*cp) ! Constant
1237 
1238  real(8) :: height ! Model level height
1239  real(8) :: sin_tmp, cos_tmp ! Calculation of great circle distance
1240  real(8) :: r, s ! Shape of perturbation
1241  real(8) :: TS ! Surface temperature
1242  real(8) :: t_mean, t_pert ! Mean and pert parts of temperature
1243  real(8) :: theta_pert ! Pot-temp perturbation
1244 
1245 !-----------------------------------------------------------------------
1246 ! THE VELOCITIES
1247 !-----------------------------------------------------------------------
1248 
1249  ! Zonal Velocity
1250 
1251  u = u0 * cos(lat)
1252 
1253  ! Meridional Velocity
1254 
1255  v = 0.d0
1256 
1257  ! Vertical Velocity = Vertical Pressure Velocity = 0
1258 
1259  w = 0.d0
1260 
1261 !-----------------------------------------------------------------------
1262 ! PHIS (surface geopotential)
1263 !-----------------------------------------------------------------------
1264 
1265  phis = 0.d0
1266 
1267 !-----------------------------------------------------------------------
1268 ! SURFACE TEMPERATURE
1269 !-----------------------------------------------------------------------
1270 
1271  ts = bigg + (teq-bigg)*exp( -(u0*n2/(4.d0*g*g))*(u0+2.d0*om*as)*(cos(2.d0*lat)-1.d0) )
1272 
1273 !-----------------------------------------------------------------------
1274 ! PS (surface pressure)
1275 !-----------------------------------------------------------------------
1276 
1277  ps = peq*exp( (u0/(4.0*bigg*rd))*(u0+2.0*om*as)*(cos(2.0*lat)-1.0) ) &
1278  * (ts/teq)**(cp/rd)
1279 
1280 !-----------------------------------------------------------------------
1281 ! HEIGHT AND PRESSURE AND MEAN TEMPERATURE
1282 !-----------------------------------------------------------------------
1283 
1284  if (zcoords .eq. 1) then
1285 
1286  height = z
1287  p = ps*( (bigg/ts)*exp(-n2*height/g)+1.d0 - (bigg/ts) )**(cp/rd)
1288 
1289  else
1290 
1291  height = (-g/n2)*log( (ts/bigg)*( (p/ps)**(rd/cp) - 1.d0 ) + 1.d0 )
1292 
1293  endif
1294 
1295  t_mean = bigg*(1.d0 - exp(n2*height/g))+ ts*exp(n2*height/g)
1296 
1297 !-----------------------------------------------------------------------
1298 ! rho (density), unperturbed using the background temperature t_mean
1299 !-----------------------------------------------------------------------
1300 
1301 !***********
1302 ! change in version 3: density is now initialized with unperturbed background temperature,
1303 ! temperature perturbation is added afterwards
1304 !***********
1305  rho = p/(rd*t_mean)
1306 
1307 !-----------------------------------------------------------------------
1308 ! POTENTIAL TEMPERATURE PERTURBATION,
1309 ! here: converted to temperature and added to the temperature field
1310 ! models with a prognostic potential temperature field can utilize
1311 ! the potential temperature perturbation theta_pert directly and add it
1312 ! to the background theta field (not included here)
1313 !-----------------------------------------------------------------------
1314 
1315  sin_tmp = sin(lat) * sin(phic)
1316  cos_tmp = cos(lat) * cos(phic)
1317 
1318  ! great circle distance with 'a/X'
1319 
1320  r = as * acos(sin_tmp + cos_tmp*cos(lon-lambdac))
1321 
1322  s = (d**2)/(d**2 + r**2)
1323 
1324  theta_pert = delta_theta*s*sin(2.d0*pi*height/lz)
1325 
1326  t_pert = theta_pert*(p/p0)**(rd/cp)
1327 
1328  t = t_mean + t_pert
1329 
1330 !-----------------------------------------------------------------------
1331 ! initialize Q, set to zero
1332 !-----------------------------------------------------------------------
1333 
1334  q = 0.d0
1335 
1336 END SUBROUTINE test3_gravity_wave
1337 
1338 
1340 
subroutine test2_steady_state_mountain(lon, lat, p, z, zcoords, hybrid_eta, hyam, hybm, u, v, w, t, phis, ps, rho, q)
subroutine test1_advection_orography(lon, lat, p, z, zcoords, cfv, hybrid_eta, hyam, hybm, gc, u, v, w, t, phis, ps, rho, q, q1, q2, q3, q4)
subroutine test1_advection_deformation(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q, q1, q2, q3, q4)
subroutine test1_advection_orograph_gal_chen_velocity(w)
real, parameter, public pi
Ratio of circle circumference to diameter [N/A].
Definition: constants.F90:74
subroutine test1_advection_orograph_hybrid_eta_velocity(w)
subroutine test2_schaer_mountain(lon, lat, p, z, zcoords, hybrid_eta, hyam, hybm, shear, u, v, w, t, phis, ps, rho, q)
subroutine test1_advection_hadley(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q, q1)
subroutine test3_gravity_wave(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q)
#define min(a, b)
Definition: mosaic_util.h:32