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, QV, SFC_TEMP, FZ, 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)  QV,
real(rp), dimension(ia,ja), intent(in)  SFC_TEMP,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
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 51 of file scale_atmos_bottom.F90.

51  use scale_const, only: &
52  grav => const_grav, &
53  rdry => const_rdry, &
54  epstvap => const_epstvap
55  implicit none
56  integer, intent(in) :: KA, KS, KE
57  integer, intent(in) :: IA, IS, IE
58  integer, intent(in) :: JA, JS, JE
59 
60  real(RP), intent(in) :: DENS (KA,IA,JA)
61  real(RP), intent(in) :: PRES (KA,IA,JA)
62  real(RP), intent(in) :: QV (KA,IA,JA)
63  real(RP), intent(in) :: SFC_TEMP(IA,JA)
64  real(RP), intent(in) :: FZ (0:KA,IA,JA)
65  real(RP), intent(out) :: SFC_DENS(IA,JA)
66  real(RP), intent(out) :: SFC_PRES(IA,JA)
67 
68  real(RP) :: Rtot
69  real(RP) :: dz1, dz2
70  real(RP) :: F2H1, F2H2
71  real(RP) :: LogP1, LogP2
72  real(RP) :: PRESH
73 
74  integer :: i, j
75  !---------------------------------------------------------------------------
76 
77  !$omp parallel do OMP_SCHEDULE_ &
78  !$omp private (Rtot,dz1,dz2,F2H1,F2H2,LogP1,LogP2,PRESH)
79  do j = js, je
80  do i = is, ie
81  ! interpolate half-level pressure
82  dz1 = fz(ks+1,i,j) - fz(ks,i,j)
83  dz2 = fz(ks,i,j) - fz(ks-1,i,j)
84 
85  f2h1 = dz2 / ( dz1 + dz2 )
86  f2h2 = dz1 / ( dz1 + dz2 )
87 
88  logp1 = log( pres(ks+1,i,j) )
89  logp2 = log( pres(ks,i,j) )
90 
91  presh = exp( f2h1 * logp1 + f2h2 * logp2 )
92 
93  rtot = rdry * ( 1.0_rp + epstvap * qv(ks,i,j) )
94  ! ( PRESH - SFC_PRES ) / dz2 = - GRAV * DENS(KS)
95  sfc_pres(i,j) = presh + grav * dens(ks,i,j) * dz2
96  sfc_dens(i,j) = sfc_pres(i,j) / ( rtot * sfc_temp(i,j) )
97  enddo
98  enddo
99 
100  return

References scale_const::const_epstvap, scale_const::const_grav, and scale_const::const_rdry.

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().

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 105 of file scale_atmos_bottom.F90.

105  implicit none
106 
107  real(RP), intent(in) :: p
108  real(RP), intent(in) :: x(3), y(3)
109  real(RP) :: q
110  !---------------------------------------------------------------------------
111 
112  q = ( (p-x(2)) * (p-x(3)) ) / ( (x(1)-x(2)) * (x(1)-x(3)) ) * y(1) &
113  + ( (p-x(1)) * (p-x(3)) ) / ( (x(2)-x(1)) * (x(2)-x(3)) ) * y(2) &
114  + ( (p-x(1)) * (p-x(2)) ) / ( (x(3)-x(1)) * (x(3)-x(2)) ) * y(3)
115 
116  return
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
scale_const::const_epstvap
real(rp), public const_epstvap
1 / epsilon - 1
Definition: scale_const.F90:70
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55