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, &
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(inout) :: SFC_Z0M (IA,JA) ! surface roughness length (momentum) [m]
82  real(RP), intent(inout) :: SFC_Z0H (IA,JA) ! surface roughness length (heat) [m]
83  real(RP), intent(inout) :: SFC_Z0E (IA,JA) ! surface roughness length (vapor) [m]
84  real(RP), intent(out) :: SFLX_MW (IA,JA) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
85  real(RP), intent(out) :: SFLX_MU (IA,JA) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
86  real(RP), intent(out) :: SFLX_MV (IA,JA) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
87  real(RP), intent(out) :: SFLX_SH (IA,JA) ! surface flux for sensible heat (area center) [J/m2/s]
88  real(RP), intent(out) :: SFLX_LH (IA,JA) ! surface flux for latent heat (area center) [J/m2/s]
89  real(RP), intent(out) :: SFLX_QTRC (IA,JA,QA) ! surface flux for tracer mass (area center) [kg/m2/s]
90  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m height
91  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m height
92  real(RP), intent(out) :: T2 (IA,JA) ! temperature t at 2m height
93  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor q at 2m height
94  end subroutine sf
95  end interface
96  procedure(sf), pointer :: atmos_phy_sf => null()
97  public :: atmos_phy_sf
98 
99 contains
100 
101  !-----------------------------------------------------------------------------
102  !
103  !-----------------------------------------------------------------------------
104  subroutine atmos_phy_sf_setup( SF_TYPE )
105  use scale_process, only: &
107  use scale_atmos_phy_sf_const, only: &
110  use scale_atmos_phy_sf_bulk, only: &
113  implicit none
114 
115  character(len=*), intent(in) :: sf_type
116  !---------------------------------------------------------------------------
117 
118  if( io_l ) write(io_fid_log,*) '*** => ', trim(sf_type), ' is selected.'
119 
120  select case( sf_type )
121  case('CONST')
122  call atmos_phy_sf_const_setup( sf_type )
124  case('BULK')
125  call atmos_phy_sf_bulk_setup( sf_type )
127  case('OFF')
128  ! do nothing
129  case('COUPLE')
130  if( io_l ) write(io_fid_log,*) '*** do nothing here'
131  case default
132  write(*,*) 'xxx invalid Surface flux type(', trim(sf_type), '). CHECK!'
133  call prc_mpistop
134  end select
135 
136  return
137  end subroutine atmos_phy_sf_setup
138 
139 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.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module STDIO
Definition: scale_stdio.F90:12
module grid index
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_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Constant flux.
module TRACER
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_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_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
Calculate surface flux.
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public atmos_phy_sf_const_setup(ATMOS_PHY_SF_TYPE)
Setup.