SCALE-RM
Functions/Subroutines
scale_atmos_bottom Module Reference

module atmosphere / bottom boundary extrapolation More...

Functions/Subroutines

subroutine, public atmos_bottom_estimate (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, CZ, Zsfc, Z1, SFC_DENS, SFC_PRES)
 Calc bottom boundary of atmosphere (just above surface) More...
 
real(rp) function lagrange_interp (p, x, y)
 

Detailed Description

module atmosphere / bottom boundary extrapolation

Description
Bottom boundary treatment of model domain Extrapolation of DENS, PRES, TEMP
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_bottom_estimate()

subroutine, public scale_atmos_bottom::atmos_bottom_estimate ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(rp), dimension (ka,ia,ja), intent(in)  CZ,
real(rp), dimension (ia,ja), intent(in)  Zsfc,
real(rp), dimension (ia,ja), intent(in)  Z1,
real(rp), dimension(ia,ja), intent(out)  SFC_DENS,
real(rp), dimension(ia,ja), intent(out)  SFC_PRES 
)

Calc bottom boundary of atmosphere (just above surface)

Definition at line 50 of file scale_atmos_bottom.F90.

References scale_const::const_grav, and lagrange_interp().

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency(), mod_atmos_driver::atmos_surface_set(), and mod_atmos_vars::atmos_vars_history_setpres().

50  use scale_const, only: &
51  grav => const_grav
52  implicit none
53  integer, intent(in) :: ka, ks, ke
54  integer, intent(in) :: ia, is, ie
55  integer, intent(in) :: ja, js, je
56 
57  real(RP), intent(in) :: dens (ka,ia,ja)
58  real(RP), intent(in) :: pres (ka,ia,ja)
59  real(RP), intent(in) :: cz (ka,ia,ja)
60  real(RP), intent(in) :: zsfc (ia,ja)
61  real(RP), intent(in) :: z1 (ia,ja)
62  real(RP), intent(out) :: sfc_dens(ia,ja)
63  real(RP), intent(out) :: sfc_pres(ia,ja)
64 
65  integer :: i, j
66  !---------------------------------------------------------------------------
67 
68  ! estimate surface density (extrapolation)
69  !$omp parallel do OMP_SCHEDULE_
70  do j = js, je
71  do i = is, ie
72  sfc_dens(i,j) = lagrange_interp( zsfc(i,j), & ! [IN]
73  cz(ks:ks+2,i,j), & ! [IN]
74  dens(ks:ks+2,i,j) ) ! [IN]
75  enddo
76  enddo
77 
78  ! estimate surface pressure (hydrostatic balance)
79  !$omp parallel do OMP_SCHEDULE_
80  do j = js, je
81  do i = is, ie
82  sfc_pres(i,j) = pres(ks,i,j) &
83  + 0.5_rp * ( sfc_dens(i,j) + dens(ks,i,j) ) * grav * z1(i,j)
84  enddo
85  enddo
86 
87  return
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
Here is the call graph for this function:
Here is the caller graph for this function:

◆ lagrange_interp()

real(rp) function scale_atmos_bottom::lagrange_interp ( real(rp), intent(in)  p,
real(rp), dimension(3), intent(in)  x,
real(rp), dimension(3), intent(in)  y 
)

Definition at line 92 of file scale_atmos_bottom.F90.

Referenced by atmos_bottom_estimate().

92  implicit none
93 
94  real(RP), intent(in) :: p
95  real(RP), intent(in) :: x(3), y(3)
96  real(RP) :: q
97  !---------------------------------------------------------------------------
98 
99  q = ( (p-x(2)) * (p-x(3)) ) / ( (x(1)-x(2)) * (x(1)-x(3)) ) * y(1) &
100  + ( (p-x(1)) * (p-x(3)) ) / ( (x(2)-x(1)) * (x(2)-x(3)) ) * y(2) &
101  + ( (p-x(1)) * (p-x(2)) ) / ( (x(3)-x(1)) * (x(3)-x(2)) ) * y(3)
102 
103  return
Here is the caller graph for this function: