SCALE-RM
scale_land_sfc.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: land_sfc_setup
27 
28  abstract interface
29  subroutine lndsfc( &
30  LST_t, &
31  ZMFLX, &
32  XMFLX, &
33  YMFLX, &
34  SHFLX, &
35  LHFLX, &
36  GHFLX, &
37  U10, &
38  V10, &
39  T2, &
40  Q2, &
41  TMPA, &
42  PRSA, &
43  WA, &
44  UA, &
45  VA, &
46  RHOA, &
47  QVA, &
48  Z1, &
49  PBL, &
50  PRSS, &
51  LWD, &
52  SWD, &
53  TG, &
54  LST, &
55  QVEF, &
56  ALB_LW, &
57  ALB_SW, &
58  DZG, &
59  TCS, &
60  Z0M, &
61  Z0H, &
62  Z0E, &
63  dt )
64  use scale_precision
66  implicit none
67 
68  real(RP), intent(out) :: lst_t(ia,ja) ! tendency of LST
69  real(RP), intent(out) :: zmflx(ia,ja) ! z-momentum flux at the surface [kg/m2/s]
70  real(RP), intent(out) :: xmflx(ia,ja) ! x-momentum flux at the surface [kg/m2/s]
71  real(RP), intent(out) :: ymflx(ia,ja) ! y-momentum flux at the surface [kg/m2/s]
72  real(RP), intent(out) :: shflx(ia,ja) ! sensible heat flux at the surface [J/m2/s]
73  real(RP), intent(out) :: lhflx(ia,ja) ! latent heat flux at the surface [J/m2/s]
74  real(RP), intent(out) :: ghflx(ia,ja) ! ground heat flux at the surface [J/m2/s]
75  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m [m/s]
76  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m [m/s]
77  real(RP), intent(out) :: t2 (ia,ja) ! temperature at 2m [K]
78  real(RP), intent(out) :: q2 (ia,ja) ! water vapor at 2m [kg/kg]
79 
80  real(RP), intent(in) :: tmpa(ia,ja) ! temperature at the lowest atmospheric layer [K]
81  real(RP), intent(in) :: prsa(ia,ja) ! pressure at the lowest atmospheric layer [Pa]
82  real(RP), intent(in) :: wa (ia,ja) ! velocity w at the lowest atmospheric layer [m/s]
83  real(RP), intent(in) :: ua (ia,ja) ! velocity u at the lowest atmospheric layer [m/s]
84  real(RP), intent(in) :: va (ia,ja) ! velocity v at the lowest atmospheric layer [m/s]
85  real(RP), intent(in) :: rhoa(ia,ja) ! density at the lowest atmospheric layer [kg/m3]
86  real(RP), intent(in) :: qva (ia,ja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
87  real(RP), intent(in) :: z1 (ia,ja) ! cell center height at the lowest atmospheric layer [m]
88  real(RP), intent(in) :: pbl (ia,ja) ! the top of atmospheric mixing layer [m]
89  real(RP), intent(in) :: prss(ia,ja) ! pressure at the surface [Pa]
90  real(RP), intent(in) :: lwd (ia,ja) ! downward long-wave radiation flux at the surface (upward positive) [J/m2/s]
91  real(RP), intent(in) :: swd (ia,ja) ! downward short-wave radiation flux at the surface (upward positive) [J/m2/s]
92 
93  real(RP), intent(in) :: tg (ia,ja) ! soil temperature [K]
94  real(RP), intent(in) :: lst (ia,ja) ! land surface temperature [K]
95  real(RP), intent(in) :: qvef (ia,ja) ! efficiency of evaporation [0-1]
96  real(RP), intent(in) :: alb_lw(ia,ja) ! surface albedo for LW [0-1]
97  real(RP), intent(in) :: alb_sw(ia,ja) ! surface albedo for SW [0-1]
98  real(RP), intent(in) :: dzg (ia,ja) ! soil depth [m]
99  real(RP), intent(in) :: tcs (ia,ja) ! thermal conductivity for soil [J/m/K/s]
100  real(RP), intent(in) :: z0m (ia,ja) ! roughness length for momemtum [m]
101  real(RP), intent(in) :: z0h (ia,ja) ! roughness length for heat [m]
102  real(RP), intent(in) :: z0e (ia,ja) ! roughness length for vapor [m]
103  real(DP), intent(in) :: dt ! delta time
104  end subroutine lndsfc
105  end interface
106  procedure(lndsfc), pointer :: land_sfc => null()
107  public :: land_sfc
108 
109  !-----------------------------------------------------------------------------
110  !
111  !++ Public parameters & variables
112  !
113  !-----------------------------------------------------------------------------
114  !
115  !++ Private procedure
116  !
117  !-----------------------------------------------------------------------------
118  !
119  !++ Private parameters & variables
120  !
121  !-----------------------------------------------------------------------------
122 contains
123 
124  subroutine land_sfc_setup( LAND_TYPE )
125  use scale_process, only: &
127  use scale_land_sfc_slab, only: &
130  implicit none
131 
132  character(len=*), intent(in) :: LAND_TYPE
133  !---------------------------------------------------------------------------
134 
135  select case( land_type )
136  case ( 'CONST' )
137  call land_sfc_slab_setup( land_type )
139  case ( 'SLAB' )
140  call land_sfc_slab_setup( land_type )
142  end select
143 
144  return
145  end subroutine land_sfc_setup
146 
147 end module scale_land_sfc
subroutine, public prc_mpistop
Abort MPI.
subroutine, public land_sfc_slab(LST_t, ZMFLX, XMFLX, YMFLX, SHFLX, LHFLX, GHFLX, U10, V10, T2, Q2, TMPA, PRSA, WA, UA, VA, RHOA, QVA, Z1, PBL, PRSS, LWD, SWD, TG, LST, QVEF, ALB_LW, ALB_SW, DZG, TCS, Z0M, Z0H, Z0E, dt)
module STDIO
Definition: scale_stdio.F90:12
module grid index
integer, public ia
of x whole cells (local, with HALO)
module PROCESS
module LAND / Surface fluxes
subroutine, public land_sfc_setup(LAND_TYPE)
procedure(lndsfc), pointer, public land_sfc
subroutine, public land_sfc_slab_setup(LAND_TYPE)
Setup.
module LAND / Surface fluxes with slab land model
module PRECISION
integer, public ja
of y whole cells (local, with HALO)