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, QV, &
48  SFC_TEMP, &
49  FZ, &
50  SFC_DENS, SFC_PRES )
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  !$acc data copyin(DENS, PRES, QV, SFC_TEMP, FZ) copyout(SFC_DENS, SFC_PRES)
78 
79  !$omp parallel do OMP_SCHEDULE_ &
80  !$omp private (Rtot,dz1,dz2,F2H1,F2H2,LogP1,LogP2,PRESH)
81  !$acc kernels
82  do j = js, je
83  do i = is, ie
84  ! interpolate half-level pressure
85  dz1 = fz(ks+1,i,j) - fz(ks,i,j)
86  dz2 = fz(ks,i,j) - fz(ks-1,i,j)
87 
88  f2h1 = dz2 / ( dz1 + dz2 )
89  f2h2 = dz1 / ( dz1 + dz2 )
90 
91  logp1 = log( pres(ks+1,i,j) )
92  logp2 = log( pres(ks,i,j) )
93 
94  presh = exp( f2h1 * logp1 + f2h2 * logp2 )
95 
96  rtot = rdry * ( 1.0_rp + epstvap * qv(ks,i,j) )
97  ! ( PRESH - SFC_PRES ) / dz2 = - GRAV * DENS(KS)
98  sfc_pres(i,j) = presh + grav * dens(ks,i,j) * dz2
99  sfc_dens(i,j) = sfc_pres(i,j) / ( rtot * sfc_temp(i,j) )
100  enddo
101  enddo
102  !$acc end kernels
103 
104  !$acc end data
105 
106  return
107  end subroutine atmos_bottom_estimate
108 
109  !-----------------------------------------------------------------------------
110  function lagrange_interp( p, x, y ) result(q)
111  !$acc routine seq
112  implicit none
113 
114  real(rp), intent(in) :: p
115  real(rp), intent(in) :: x(3), y(3)
116  real(rp) :: q
117  !---------------------------------------------------------------------------
118 
119  q = ( (p-x(2)) * (p-x(3)) ) / ( (x(1)-x(2)) * (x(1)-x(3)) ) * y(1) &
120  + ( (p-x(1)) * (p-x(3)) ) / ( (x(2)-x(1)) * (x(2)-x(3)) ) * y(2) &
121  + ( (p-x(1)) * (p-x(2)) ) / ( (x(3)-x(1)) * (x(3)-x(2)) ) * y(3)
122 
123  return
124  end function lagrange_interp
125 
126 end module scale_atmos_bottom
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
scale_const::const_epstvap
real(rp), public const_epstvap
1 / epsilon - 1
Definition: scale_const.F90:76
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_bottom::lagrange_interp
real(rp) function lagrange_interp(p, x, y)
Definition: scale_atmos_bottom.F90:111
scale_atmos_bottom
module atmosphere / bottom boundary extrapolation
Definition: scale_atmos_bottom.F90:12
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_atmos_bottom::atmos_bottom_estimate
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)
Definition: scale_atmos_bottom.F90:51