SCALE-RM
scale_atmos_bottom.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
27  public :: atmos_bottom_estimate
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  !-----------------------------------------------------------------------------
34  !
35  !++ Private procedure
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private parameters & variables
40  !
41  !-----------------------------------------------------------------------------
42 contains
43  !-----------------------------------------------------------------------------
45  subroutine atmos_bottom_estimate( &
46  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
47  DENS, PRES, &
48  CZ, Zsfc, Z1, &
49  SFC_DENS, SFC_PRES )
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
88  end subroutine atmos_bottom_estimate
89 
90  !-----------------------------------------------------------------------------
91  function lagrange_interp( p, x, y ) result(q)
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
104  end function lagrange_interp
105 
106 end module scale_atmos_bottom
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)
module atmosphere / bottom boundary extrapolation
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
module profiler
Definition: scale_prof.F90:11
real(rp) function lagrange_interp(p, x, y)
module PRECISION
module STDIO
Definition: scale_io.F90:10