SCALE-RM
scale_atmos_phy_sf.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_tracer
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: atmos_phy_sf_setup
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  abstract interface
49  subroutine sf( &
50  ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, &
51  ATM_DENS, &
52  ATM_QTRC, &
53  ATM_Z1, dt, &
54  SFC_DENS, SFC_PRES, &
55  SFLX_LW_dn, SFLX_SW_dn, &
56  SFC_TEMP, SFC_albedo, SFC_beta, &
57  SFC_Z0M, SFC_Z0H, SFC_Z0E, &
58  SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, &
59  SFLX_QTRC, &
60  U10, V10, T2, Q2 )
61  use scale_precision
63  use scale_tracer
64  implicit none
65 
66  real(RP), intent(in) :: atm_temp (ia,ja) ! temperature at the lowermost layer (cell center) [K]
67  real(RP), intent(in) :: atm_pres (ia,ja) ! pressure at the lowermost layer (cell center) [Pa]
68  real(RP), intent(in) :: atm_w (ia,ja) ! velocity w at the lowermost layer (cell center) [m/s]
69  real(RP), intent(in) :: atm_u (ia,ja) ! velocity u at the lowermost layer (cell center) [m/s]
70  real(RP), intent(in) :: atm_v (ia,ja) ! velocity v at the lowermost layer (cell center) [m/s]
71  real(RP), intent(in) :: atm_dens (ia,ja) ! density at the lowermost layer (cell center) [kg/m3]
72  real(RP), intent(in) :: atm_qtrc (ia,ja,qa) ! tracer at the lowermost layer (cell center) [kg/kg]
73  real(RP), intent(in) :: atm_z1 (ia,ja) ! height of the lowermost grid from surface (cell center) [m]
74  real(DP), intent(in) :: dt ! delta time
75  real(RP), intent(in) :: sfc_dens (ia,ja) ! density at the surface atmosphere [kg/m3]
76  real(RP), intent(in) :: sfc_pres (ia,ja) ! pressure at the surface atmosphere [Pa]
77  real(RP), intent(in) :: sflx_lw_dn(ia,ja) ! downward longwave radiation flux at the surface [J/m2/s]
78  real(RP), intent(in) :: sflx_sw_dn(ia,ja) ! downward shortwave radiation flux at the surface [J/m2/s]
79  real(RP), intent(in) :: sfc_temp (ia,ja) ! temperature at the surface skin [K]
80  real(RP), intent(in) :: sfc_albedo(ia,ja,2) ! surface albedo (LW/SW) [0-1]
81  real(RP), intent(in) :: sfc_beta (ia,ja) ! evaporation efficiency [0-1]
82  real(RP), intent(inout) :: sfc_z0m (ia,ja) ! surface roughness length (momentum) [m]
83  real(RP), intent(inout) :: sfc_z0h (ia,ja) ! surface roughness length (heat) [m]
84  real(RP), intent(inout) :: sfc_z0e (ia,ja) ! surface roughness length (vapor) [m]
85  real(RP), intent(out) :: sflx_mw (ia,ja) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
86  real(RP), intent(out) :: sflx_mu (ia,ja) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
87  real(RP), intent(out) :: sflx_mv (ia,ja) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
88  real(RP), intent(out) :: sflx_sh (ia,ja) ! surface flux for sensible heat (area center) [J/m2/s]
89  real(RP), intent(out) :: sflx_lh (ia,ja) ! surface flux for latent heat (area center) [J/m2/s]
90  real(RP), intent(out) :: sflx_qtrc (ia,ja,qa) ! surface flux for tracer mass (area center) [kg/m2/s]
91  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m height
92  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m height
93  real(RP), intent(out) :: t2 (ia,ja) ! temperature t at 2m height
94  real(RP), intent(out) :: q2 (ia,ja) ! water vapor q at 2m height
95  end subroutine sf
96  end interface
97  procedure(sf), pointer :: atmos_phy_sf => null()
98  public :: atmos_phy_sf
99 
100 contains
101 
102  !-----------------------------------------------------------------------------
103  !
104  !-----------------------------------------------------------------------------
105  subroutine atmos_phy_sf_setup( SF_TYPE )
106  use scale_process, only: &
108 #define EXTM(pre, name, post) pre ## name ## post
109 #define NAME(pre, name, post) EXTM(pre, name, post)
110 #ifdef SF
111  use name(scale_atmos_phy_, sf,), only: &
112  name(atmos_phy_, sf, _setup), &
113  name(atmos_phy_, sf,)
114 #else
115  use scale_atmos_phy_sf_const, only: &
118  use scale_atmos_phy_sf_bulk, only: &
121 #endif
122  implicit none
123 
124  character(len=*), intent(in) :: SF_TYPE
125  !---------------------------------------------------------------------------
126 
127  select case( sf_type )
128  case('CONST')
129  call atmos_phy_sf_const_setup( sf_type )
131  case('BULK')
132  call atmos_phy_sf_bulk_setup( sf_type )
134  case('OFF', 'COUPLE')
135  ! do nothing
136  case default
137  write(*,*) 'xxx ATMPS_PHY_SF_TYPE is invalid'
138  call prc_mpistop
139  end select
140 
141  return
142  end subroutine atmos_phy_sf_setup
143 
144 end module scale_atmos_phy_sf
subroutine, public atmos_phy_sf_setup(SF_TYPE)
subroutine, public atmos_phy_sf_bulk_setup(ATMOS_PHY_SF_TYPE)
Setup.
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_sf_const(ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, ATM_DENS, ATM_QTRC, ATM_Z1, dt, SFC_DENS, SFC_PRES, SFLX_LW_dn, SFLX_SW_dn, SFC_TEMP, SFC_albedo, SFC_beta, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Constant flux.
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
module grid index
subroutine, public atmos_phy_sf_bulk(ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, ATM_DENS, ATM_QTRC, ATM_Z1, dt, SFC_DENS, SFC_PRES, SFLX_LW_dn, SFLX_SW_dn, SFC_TEMP, SFC_albedo, SFC_beta, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Calculate surface flux.
module TRACER
integer, public ia
of x whole cells (local, with HALO)
procedure(sf), pointer, public atmos_phy_sf
module ATMOSPHERE / Physics Surface fluxes
module PROCESS
module ATMOSPHERE / Physics Surface fluxes
module profiler
Definition: scale_prof.F90:10
module ATMOSPHERE / Physics Surface fluxes
module PRECISION
subroutine, public atmos_phy_sf_const_setup(ATMOS_PHY_SF_TYPE)
Setup.
integer, public ja
of y whole cells (local, with HALO)