SCALE-RM
scale_atmos_sub_bottom.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
18  use scale_prof
20  use scale_tracer
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: atmos_bottom_estimate
29 
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public parameters & variables
33  !
34  !-----------------------------------------------------------------------------
35  !
36  !++ Private procedure
37  !
38  !-----------------------------------------------------------------------------
39  !
40  !++ Private parameters & variables
41  !
42  !-----------------------------------------------------------------------------
43 contains
44  !-----------------------------------------------------------------------------
46  subroutine atmos_bottom_estimate( &
47  DENS, &
48  PRES, &
49  CZ, &
50  Zsfc, &
51  Z1, &
52  SFC_DENS, &
53  SFC_PRES )
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
87  end subroutine atmos_bottom_estimate
88 
89  !-----------------------------------------------------------------------------
90  function lagrange_interp( p, x, y ) result(q)
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
103  end function lagrange_interp
104 
105 end module scale_atmos_bottom
integer, public jeb
module ATMOSPHERE / Bottom boundary treatment
module STDIO
Definition: scale_stdio.F90:12
integer, public ieb
module grid index
module TRACER
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module CONSTANT
Definition: scale_const.F90:14
subroutine, public atmos_bottom_estimate(DENS, PRES, CZ, Zsfc, Z1, SFC_DENS, SFC_PRES)
Calc bottom boundary of atmosphere (just above surface)
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
real(rp) function lagrange_interp(p, x, y)
module PRECISION
integer, public isb
integer, public jsb
integer, public ja
of whole cells: y, local, with HALO