FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
scale_atm_dyn_dgm_hydrostatic Module Reference

module FElib / Fluid dyn solver / Atmosphere / Common More...

Data Types

interface  hydrostaic_build_rho_xyz
 

Functions/Subroutines

subroutine, public hydrostatic_calc_basicstate_constt (dens_hyd, pres_hyd, temp0, pres_sfc, x, y, z, lcmesh3d, elem)
 Calculate density and pressure in hydrostatic balance with a constant temperature.
 
subroutine, public hydrostatic_calc_basicstate_constpt (dens_hyd, pres_hyd, pottemp0, pres_sfc, x, y, z, lcmesh3d, elem)
 Calculate density and pressure in hydrostatic balance with a constant potential temperature.
 
subroutine, public hydrostatic_calc_basicstate_constbvfreq (dens_hyd, pres_hyd, bruntvaisalafreq, pottemp0, pres_sfc, x, y, z, lcmesh3d, elem)
 Calculate density and pressure in hydrostatic balance with a constant Brunt–Väisälä frequency.
 
subroutine, public hydrostatic_calc_basicstate_consttlaps (dens_hyd, pres_hyd, tlaps, temp0, pres_sfc, x, y, z, lcmesh3d, elem)
 Calculate density and pressure in hydrostatic balance with a constant lapse rate of temperature.
 
subroutine, public hydrostatic_calc_basicstate_constptlaps (dens_hyd, pres_hyd, ptlaps, pottemp0, pres_sfc, x, y, z, lcmesh3d, elem)
 Calculate density and pressure in hydrostatic balance with a constant lapse rate of potential temperature.
 

Detailed Description

module FElib / Fluid dyn solver / Atmosphere / Common

Description
Construct hydrostatic state for atmospheric dynamical process.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ hydrostatic_calc_basicstate_constt()

subroutine, public scale_atm_dyn_dgm_hydrostatic::hydrostatic_calc_basicstate_constt ( real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) dens_hyd,
real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) pres_hyd,
real(rp), intent(in) temp0,
real(rp), intent(in) pres_sfc,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) x,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) y,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) z,
class(localmesh3d), intent(in) lcmesh3d,
class(elementbase3d), intent(in) elem )

Calculate density and pressure in hydrostatic balance with a constant temperature.

Parameters
DENS_hydhydrostatic density [kg/m3]
PRES_hydhydrostatic pressure [Pa]
Temp0Temperature of isothermal atmosphere [K]
PRES_sfcSurface pressure [Pa]
xx-coordinate
yy-coordinate
zz-coordinate [m]
lcmesh3DA object to manage a local 3D mesh
elemA object to manage a 3D finite element

Definition at line 86 of file scale_atm_dyn_dgm_hydrostatic.F90.

89
90 implicit none
91
92 class(LocalMesh3D), intent(in) :: lcmesh3D
93 class(ElementBase3D), intent(in) :: elem
94 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh3D%NeA)
95 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh3D%NeA)
96 real(RP), intent(in) :: x(elem%Np,lcmesh3D%Ne)
97 real(RP), intent(in) :: y(elem%Np,lcmesh3D%Ne)
98 real(RP), intent(in) :: z(elem%Np,lcmesh3D%Ne)
99 real(RP), intent(in) :: Temp0
100 real(RP), intent(in) :: PRES_sfc
101
102 integer :: ke
103 real(RP) :: H0
104 !-----------------------------------------------
105
106 h0 = rdry * temp0 / grav
107
108 !$omp parallel do
109 do ke=lcmesh3d%NeS, lcmesh3d%NeE
110 pres_hyd(:,ke) = pres_sfc * exp( - z(:,ke) / h0 )
111 dens_hyd(:,ke) = pres_hyd(:,ke) / ( rdry * temp0 )
112 end do
113
114 return

◆ hydrostatic_calc_basicstate_constpt()

subroutine, public scale_atm_dyn_dgm_hydrostatic::hydrostatic_calc_basicstate_constpt ( real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) dens_hyd,
real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) pres_hyd,
real(rp), intent(in) pottemp0,
real(rp), intent(in) pres_sfc,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) x,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) y,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) z,
class(localmesh3d), intent(in) lcmesh3d,
class(elementbase3d), intent(in) elem )

Calculate density and pressure in hydrostatic balance with a constant potential temperature.

Parameters
DENS_hydhydrostatic density [kg/m3]
PRES_hydhydrostatic pressure [Pa]
PotTemp0Constant potential temperature of atmosphere [K]
PRES_sfcSurface pressure [Pa]
xx-coordinate
yy-coordinate
zz-coordinate [m]
lcmesh3DA object to manage a local 3D mesh
elemA object to manage a 3D finite element

Definition at line 129 of file scale_atm_dyn_dgm_hydrostatic.F90.

132
133 implicit none
134
135 class(LocalMesh3D), intent(in) :: lcmesh3D
136 class(ElementBase3D), intent(in) :: elem
137 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh3D%NeA)
138 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh3D%NeA)
139 real(RP), intent(in) :: x(elem%Np,lcmesh3D%Ne)
140 real(RP), intent(in) :: y(elem%Np,lcmesh3D%Ne)
141 real(RP), intent(in) :: z(elem%Np,lcmesh3D%Ne)
142 real(RP), intent(in) :: PotTemp0
143 real(RP), intent(in) :: PRES_sfc
144
145 integer :: ke
146
147 real(RP) :: exner(elem%Np)
148 real(RP) :: exner_sfc
149 real(RP) :: RovCP
150 real(RP) :: CPovR
151 !-----------------------------------------------
152
153 rovcp = rdry / cpdry
154 cpovr = cpdry / rdry
155 exner_sfc = (pres_sfc / pres00)**rovcp
156
157 !$omp parallel do private(exner)
158 do ke=lcmesh3d%NeS, lcmesh3d%NeE
159 ! Cp * PT0 * d exner / dz = - g
160 exner(:) = exner_sfc - grav / (cpdry * pottemp0) * z(:,ke)
161 pres_hyd(:,ke) = pres00 * exner(:)**cpovr
162 dens_hyd(:,ke) = pres_hyd(:,ke) / ( rdry * exner(:) * pottemp0 )
163 end do
164
165 return

◆ hydrostatic_calc_basicstate_constbvfreq()

subroutine, public scale_atm_dyn_dgm_hydrostatic::hydrostatic_calc_basicstate_constbvfreq ( real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) dens_hyd,
real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) pres_hyd,
real(rp), intent(in) bruntvaisalafreq,
real(rp), intent(in) pottemp0,
real(rp), intent(in) pres_sfc,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) x,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) y,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) z,
class(localmesh3d), intent(in) lcmesh3d,
class(elementbase3d), intent(in) elem )

Calculate density and pressure in hydrostatic balance with a constant Brunt–Väisälä frequency.

Parameters
DENS_hydhydrostatic density [kg/m3]
PRES_hydhydrostatic pressure [Pa]
BruntVaisalaFreqConstant Brunt–Väisälä frequency [s-1]
PotTemp0Constant potential temperature of atmosphere [K]
PRES_sfcSurface pressure [Pa]
xx-coordinate
yy-coordinate
zz-coordinate [m]
lcmesh3DA object to manage a local 3D mesh
elemA object to manage a 3D finite element

Definition at line 181 of file scale_atm_dyn_dgm_hydrostatic.F90.

184
185 implicit none
186
187 class(LocalMesh3D), intent(in) :: lcmesh3D
188 class(ElementBase3D), intent(in) :: elem
189 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh3D%NeA)
190 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh3D%NeA)
191 real(RP), intent(in) :: x(elem%Np,lcmesh3D%Ne)
192 real(RP), intent(in) :: y(elem%Np,lcmesh3D%Ne)
193 real(RP), intent(in) :: z(elem%Np,lcmesh3D%Ne)
194 real(RP), intent(in) :: BruntVaisalaFreq
195 real(RP), intent(in) :: PotTemp0
196 real(RP), intent(in) :: PRES_sfc
197
198 integer :: ke
199
200 real(RP) :: PT(elem%NP)
201 real(RP) :: exner(elem%Np)
202 real(RP) :: exner_sfc
203 real(RP) :: RovCP
204 real(RP) :: CPovR
205 !-----------------------------------------------
206
207 rovcp = rdry / cpdry
208 cpovr = cpdry / rdry
209 exner_sfc = (pres_sfc / pres00)**rovcp
210
211 !$omp parallel do private(PT, exner)
212 do ke=lcmesh3d%NeS, lcmesh3d%NeE
213 ! d exner / dz = - g / ( Cp * PT0 ) * exp (- N2/g * z)
214 ! exner = exner(zs) - g^2 / (Cp * N^2) [ 1/PT (z) - 1/PT(zs) ]
215 pt(:) = pottemp0 * exp( bruntvaisalafreq**2 / grav * z(:,ke) )
216 exner(:) = exner_sfc + grav**2 / ( cpdry * bruntvaisalafreq**2 ) * ( 1.0_rp / pt(:) - 1.0_rp / pottemp0 )
217
218 pres_hyd(:,ke) = pres00 * exner(:)**cpovr
219 dens_hyd(:,ke) = pres_hyd(:,ke) / ( rdry * exner(:) * pt(:) )
220 end do
221
222 return

◆ hydrostatic_calc_basicstate_consttlaps()

subroutine, public scale_atm_dyn_dgm_hydrostatic::hydrostatic_calc_basicstate_consttlaps ( real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) dens_hyd,
real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) pres_hyd,
real(rp), intent(in) tlaps,
real(rp), intent(in) temp0,
real(rp), intent(in) pres_sfc,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) x,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) y,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) z,
class(localmesh3d), intent(in) lcmesh3d,
class(elementbase3d), intent(in) elem )

Calculate density and pressure in hydrostatic balance with a constant lapse rate of temperature.

Parameters
DENS_hydhydrostatic density [kg/m3]
PRES_hydhydrostatic pressure [Pa]
TLAPSConstant lapse rate of atmosphere [K/m]
Temp0Temperature at z=0 [K]
PRES_sfcSurface pressure [Pa]
xx-coordinate
yy-coordinate
zz-coordinate [m]
lcmesh3DA object to manage a local 3D mesh
elemA object to manage a 3D finite element

Definition at line 238 of file scale_atm_dyn_dgm_hydrostatic.F90.

241
242 implicit none
243
244 class(LocalMesh3D), intent(in) :: lcmesh3D
245 class(ElementBase3D), intent(in) :: elem
246 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh3D%NeA)
247 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh3D%NeA)
248 real(RP), intent(in) :: x(elem%Np,lcmesh3D%Ne)
249 real(RP), intent(in) :: y(elem%Np,lcmesh3D%Ne)
250 real(RP), intent(in) :: z(elem%Np,lcmesh3D%Ne)
251 real(RP), intent(in) :: TLAPS
252 real(RP), intent(in) :: Temp0
253 real(RP), intent(in) :: PRES_sfc
254
255 integer :: ke
256
257 real(RP) :: fac
258 real(RP) :: TEMP(elem%Np)
259 !-----------------------------------------------
260
261 fac = grav / ( rdry * tlaps )
262
263 !$omp parallel do private(TEMP)
264 do ke=lcmesh3d%NeS, lcmesh3d%NeE
265 temp(:) = temp0 * ( 1.0_rp - tlaps / temp0 * z(:,ke) )
266 pres_hyd(:,ke) = pres00 * ( temp(:) / temp0 )**fac
267 dens_hyd(:,ke) = pres_hyd(:,ke) / ( rdry * temp(:) )
268 end do
269
270 return

◆ hydrostatic_calc_basicstate_constptlaps()

subroutine, public scale_atm_dyn_dgm_hydrostatic::hydrostatic_calc_basicstate_constptlaps ( real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) dens_hyd,
real(rp), dimension(elem%np,lcmesh3d%nea), intent(out) pres_hyd,
real(rp), intent(in) ptlaps,
real(rp), intent(in) pottemp0,
real(rp), intent(in) pres_sfc,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) x,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) y,
real(rp), dimension(elem%np,lcmesh3d%ne), intent(in) z,
class(localmesh3d), intent(in) lcmesh3d,
class(elementbase3d), intent(in) elem )

Calculate density and pressure in hydrostatic balance with a constant lapse rate of potential temperature.

Parameters
DENS_hydhydrostatic density [kg/m3]
PRES_hydhydrostatic pressure [Pa]
PTLAPSConstant lapse rate of potential temperature in atmosphere [K/m]
PotTemp0Potentital temperature at z=0 [K]
PRES_sfcSurface pressure [Pa]
xx-coordinate
yy-coordinate
zz-coordinate [m]
lcmesh3DA object to manage a local 3D mesh
elemA object to manage a 3D finite element

Definition at line 286 of file scale_atm_dyn_dgm_hydrostatic.F90.

289
290 implicit none
291
292 class(LocalMesh3D), intent(in) :: lcmesh3D
293 class(ElementBase3D), intent(in) :: elem
294 real(RP), intent(out) :: DENS_hyd(elem%Np,lcmesh3D%NeA)
295 real(RP), intent(out) :: PRES_hyd(elem%Np,lcmesh3D%NeA)
296 real(RP), intent(in) :: x(elem%Np,lcmesh3D%Ne)
297 real(RP), intent(in) :: y(elem%Np,lcmesh3D%Ne)
298 real(RP), intent(in) :: z(elem%Np,lcmesh3D%Ne)
299 real(RP), intent(in) :: PTLAPS
300 real(RP), intent(in) :: PotTemp0
301 real(RP), intent(in) :: PRES_sfc
302
303 integer :: ke
304
305 real(RP) :: PT(elem%NP)
306 real(RP) :: exner(elem%Np)
307 real(RP) :: exner_sfc
308 real(RP) :: RovCP
309 real(RP) :: CPovR
310 !-----------------------------------------------
311
312 rovcp = rdry / cpdry
313 cpovr = cpdry / rdry
314 exner_sfc = (pres_sfc / pres00)**rovcp
315
316 !$omp parallel do private(PT, exner)
317 do ke=lcmesh3d%NeS, lcmesh3d%NeE
318 ! d exner / dz = - g / ( Cp * PT0 ) / (1 + PTLAPS/PT0 * z)
319 ! exner = exner(zs) - g / (Cp * PTLAPS ) * log[ 1 + PTLAPS/PT0 * z ]
320 pt(:) = pottemp0 + ptlaps * z(:,ke)
321 exner(:) = exner_sfc - grav / ( cpdry * ptlaps ) * log( 1.0_rp + ptlaps / pottemp0 * z(:,ke) )
322
323 pres_hyd(:,ke) = pres00 * exner(:)**cpovr
324 dens_hyd(:,ke) = pres_hyd(:,ke) / ( rdry * exner(:) * pt(:) )
325 end do
326
327 return