FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
mod_mktopo_util Module Reference

module Utility for mktopo More...

Functions/Subroutines

subroutine, public mktopoutil_barocwave_global_jw2006_calc_topo (topo, u0, eta0, lat, np)
 Calculate a topography for a global test case of idealized baroclinic wave in Jablonowski and Williamson (2006)
 

Detailed Description

module Utility for mktopo

Description
subroutines useful to prepare topography data
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ mktopoutil_barocwave_global_jw2006_calc_topo()

subroutine, public mod_mktopo_util::mktopoutil_barocwave_global_jw2006_calc_topo ( real(rp), dimension(np), intent(out) topo,
real(rp), intent(in) u0,
real(rp), intent(in) eta0,
real(rp), dimension(np), intent(in) lat,
integer, intent(in) np )

Calculate a topography for a global test case of idealized baroclinic wave in Jablonowski and Williamson (2006)

Reference
  • Jablonowski, C., & Williamson, D. L., 2006: A baroclinic instability test case for atmospheric model dynamical cores. Quarterly Journal of the Royal Meteorological Society: A journal of the atmospheric sciences, applied meteorology and physical oceanography, 132(621C), 2943-2975.
Parameters
[in]u0The value of η at a reference level (position of the jet)
[in]eta0The parameter associated with zonal jet maximum amplitude [m/s]
[in]latlatitude [radian]

Definition at line 42 of file mod_mktopo_util.F90.

44
45 use scale_const, only: &
46 ohm => const_ohm, &
47 grav => const_grav, &
48 rplanet => const_radius
49
50 implicit none
51
52 integer, intent(in) :: Np
53 real(RP), intent(out) :: topo(Np)
54 real(RP), intent(in) :: U0
55 real(RP), intent(in) :: ETA0
56 real(RP), intent(in) :: lat(Np)
57
58 real(RP) :: sin_lat(Np)
59 real(RP) :: sin_lat_pow_6(Np)
60 real(RP) :: cos_lat(Np)
61
62 real(RP ) :: tmp
63 !-------------------------------------------
64
65 sin_lat(:) = sin(lat(:))
66 sin_lat_pow_6(:) = sin_lat(:)**6
67 cos_lat(:) = cos(lat(:))
68
69 tmp = cos( (1.0_rp - eta0) * 0.5_rp * pi )
70 tmp = u0 * tmp * sqrt(tmp)
71
72 ! Calc horizontal variation of geopotential height
73 topo(:) = tmp * &
74 ( tmp * ( - 2.0_rp * sin_lat_pow_6(:) * ( cos_lat(:)**2 + 1.0_rp / 3.0_rp ) + 10.0_rp / 63.0_rp ) &
75 + rplanet * ohm * ( 8.0_rp / 5.0_rp * cos_lat(:)**3 * ( sin_lat(:)**2 + 2.0_rp / 3.0_rp ) - 0.25_rp * pi ) &
76 ) / grav
77
78 return

Referenced by mod_mktopo::mktopo_write().