SCALE-RM
Functions/Subroutines
scale_atmos_bottom Module Reference

module ATMOSPHERE / Bottom boundary treatment More...

Functions/Subroutines

subroutine, public atmos_bottom_estimate (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 treatment

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 ( 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 54 of file scale_atmos_sub_bottom.F90.

References scale_const::const_grav, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ks, and lagrange_interp().

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

54  use scale_const, only: &
55  grav => const_grav
56  implicit none
57 
58  real(RP), intent(in) :: DENS (KA,IA,JA)
59  real(RP), intent(in) :: PRES (KA,IA,JA)
60  real(RP), intent(in) :: CZ (KA,IA,JA)
61  real(RP), intent(in) :: Zsfc (IA,JA)
62  real(RP), intent(in) :: Z1 (IA,JA)
63  real(RP), intent(out) :: SFC_DENS(IA,JA)
64  real(RP), intent(out) :: SFC_PRES(IA,JA)
65 
66  integer :: i, j
67  !---------------------------------------------------------------------------
68 
69  ! estimate surface density (extrapolation)
70  do j = jsb, jeb
71  do i = isb, ieb
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  do j = jsb, jeb
80  do i = isb, ieb
81  sfc_pres(i,j) = pres(ks,i,j) &
82  + 0.5_rp * ( sfc_dens(i,j) + dens(ks,i,j) ) * grav * z1(i,j)
83  enddo
84  enddo
85 
86  return
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module CONSTANT
Definition: scale_const.F90:14
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 91 of file scale_atmos_sub_bottom.F90.

Referenced by atmos_bottom_estimate().

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