SCALE-RM
Functions/Subroutines
scale_land_phy_slab Module Reference

module LAND / Physics Slab model More...

Functions/Subroutines

subroutine, public land_phy_slab_setup (LAND_TYPE)
 Setup. More...
 
subroutine, public land_phy_slab (TEMP_t, WATER_t, TEMP, WATER, WaterLimit, ThermalCond, HeatCapacity, WaterDiff, SFLX_GH, SFLX_prec, SFLX_evap, CDZ, dt)
 Physical processes for land submodel. More...
 

Detailed Description

module LAND / Physics Slab model

Description
slab-type land physics module
Author
Team SCALE
NAMELIST
  • PARAM_LAND_PHY_SLAB
    nametypedefault valuecomment
    LAND_PHY_UPDATE_BOTTOM_TEMP logical .false. Is LAND_TEMP updated in the lowest level?
    LAND_PHY_UPDATE_BOTTOM_WATER logical .false. Is LAND_WATER updated in the lowest level?

History Output
No history output

Function/Subroutine Documentation

◆ land_phy_slab_setup()

subroutine, public scale_land_phy_slab::land_phy_slab_setup ( character(len=*), intent(in)  LAND_TYPE)

Setup.

Definition at line 53 of file scale_land_phy_slab.F90.

References scale_const::const_cl, scale_const::const_dwatr, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, and scale_process::prc_mpistop().

Referenced by scale_land_phy::land_phy_setup().

53  use scale_process, only: &
55  use scale_const, only: &
56  dwatr => const_dwatr, &
57  cl => const_cl
58  implicit none
59 
60  character(len=*), intent(in) :: LAND_TYPE
61 
62  namelist / param_land_phy_slab / &
63  land_phy_update_bottom_temp, &
64  land_phy_update_bottom_water
65 
66  integer :: ierr
67  !---------------------------------------------------------------------------
68 
69  if( io_l ) write(io_fid_log,*)
70  if( io_l ) write(io_fid_log,*) '++++++ Module[SLAB] / Categ[LAND PHY] / Origin[SCALElib]'
71 
72  !--- read namelist
73  rewind(io_fid_conf)
74  read(io_fid_conf,nml=param_land_phy_slab,iostat=ierr)
75  if( ierr < 0 ) then !--- missing
76  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
77  elseif( ierr > 0 ) then !--- fatal error
78  write(*,*) 'xxx Not appropriate names in namelist PARAM_LAND_PHY_SLAB. Check!'
79  call prc_mpistop
80  endif
81  if( io_nml ) write(io_fid_nml,nml=param_land_phy_slab)
82 
83  water_denscs = dwatr * cl
84 
85  if( io_l ) write(io_fid_log,*)
86  if( io_l ) write(io_fid_log,*) '*** Update soil temperature of bottom layer? : ', land_phy_update_bottom_temp
87  if( io_l ) write(io_fid_log,*) '*** Update soil moisture of bottom layer? : ', land_phy_update_bottom_water
88 
89  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_phy_slab()

subroutine, public scale_land_phy_slab::land_phy_slab ( real(rp), dimension (lkmax,ia,ja), intent(out)  TEMP_t,
real(rp), dimension (lkmax,ia,ja), intent(out)  WATER_t,
real(rp), dimension (lkmax,ia,ja), intent(in)  TEMP,
real(rp), dimension (lkmax,ia,ja), intent(in)  WATER,
real(rp), dimension (ia,ja), intent(in)  WaterLimit,
real(rp), dimension (ia,ja), intent(in)  ThermalCond,
real(rp), dimension(ia,ja), intent(in)  HeatCapacity,
real(rp), dimension (ia,ja), intent(in)  WaterDiff,
real(rp), dimension (ia,ja), intent(in)  SFLX_GH,
real(rp), dimension (ia,ja), intent(in)  SFLX_prec,
real(rp), dimension (ia,ja), intent(in)  SFLX_evap,
real(rp), dimension (lkmax), intent(in)  CDZ,
real(dp), intent(in)  dt 
)

Physical processes for land submodel.

Definition at line 108 of file scale_land_phy_slab.F90.

References scale_const::const_dwatr, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_landuse::landuse_fact_land, scale_land_grid_index::lke, scale_land_grid_index::lkmax, and scale_land_grid_index::lks.

Referenced by scale_land_phy::land_phy_setup().

108  use scale_const, only: &
109  dwatr => const_dwatr
110  use scale_landuse, only: &
112  use scale_matrix, only: &
113  matrix_solver_tridiagonal
114  implicit none
115 
116  ! arguments
117  real(RP), intent(out) :: TEMP_t (LKMAX,IA,JA)
118  real(RP), intent(out) :: WATER_t (LKMAX,IA,JA)
119 
120  real(RP), intent(in) :: TEMP (LKMAX,IA,JA)
121  real(RP), intent(in) :: WATER (LKMAX,IA,JA)
122  real(RP), intent(in) :: WaterLimit (IA,JA)
123  real(RP), intent(in) :: ThermalCond (IA,JA)
124  real(RP), intent(in) :: HeatCapacity(IA,JA)
125  real(RP), intent(in) :: WaterDiff (IA,JA)
126  real(RP), intent(in) :: SFLX_GH (IA,JA)
127  real(RP), intent(in) :: SFLX_prec (IA,JA)
128  real(RP), intent(in) :: SFLX_evap (IA,JA)
129  real(RP), intent(in) :: CDZ (LKMAX)
130  real(DP), intent(in) :: dt
131 
132  ! work
133  real(RP) :: TEMP1 (LKMAX,IA,JA)
134  real(RP) :: WATER1(LKMAX,IA,JA)
135 
136  real(RP) :: LAND_DENSCS(LKMAX,IA,JA)
137  real(RP) :: ThermalDiff(LKMAX,IA,JA)
138 
139 ! real(RP) :: RUNOFF(IA,JA)
140 
141  real(RP) :: U(LKMAX,IA,JA)
142  real(RP) :: M(LKMAX,IA,JA)
143  real(RP) :: L(LKMAX,IA,JA)
144  real(RP) :: V(LKMAX,IA,JA)
145 
146  integer :: k, i, j
147  !---------------------------------------------------------------------------
148 
149  if( io_l ) write(io_fid_log,*) '*** Land physics step: Slab'
150 
151  ! Solve diffusion of soil moisture (tridiagonal matrix)
152  do j = js, je
153  do i = is, ie
154  l(lks,i,j) = 0.0_rp
155  u(lks,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(lks) * ( cdz(lks) + cdz(lks+1) ) ) * dt
156  l(lke,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(lke) * ( cdz(lke) + cdz(lke-1) ) ) * dt
157  u(lke,i,j) = 0.0_rp
158 
159  m(lks,i,j) = 1.0_rp - l(lks,i,j) - u(lks,i,j)
160  m(lke,i,j) = 1.0_rp - l(lke,i,j) - u(lke,i,j)
161  end do
162  end do
163 
164  do j = js, je
165  do i = is, ie
166  do k = lks+1, lke-1
167  l(k,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(k) * ( cdz(k) + cdz(k-1) ) ) * dt
168  u(k,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(k) * ( cdz(k) + cdz(k+1) ) ) * dt
169  m(k,i,j) = 1.0_rp - l(k,i,j) - u(k,i,j)
170  end do
171  end do
172  end do
173 
174  ! input from atmosphere
175  do j = js, je
176  do i = is, ie
177  v(lks,i,j) = water(lks,i,j) + ( sflx_prec(i,j) - sflx_evap(i,j) ) / ( cdz(lks) * dwatr ) * dt
178  end do
179  end do
180 
181  do j = js, je
182  do i = is, ie
183  do k = lks+1, lke
184  v(k,i,j) = water(k,i,j)
185  end do
186  end do
187  end do
188 
189  call matrix_solver_tridiagonal( lkmax, & ! [IN]
190  ia, is, ie, & ! [IN]
191  ja, js, je, & ! [IN]
192  u(:,:,:), & ! [IN]
193  m(:,:,:), & ! [IN]
194  l(:,:,:), & ! [IN]
195  v(:,:,:), & ! [IN]
196  water1(:,:,:) ) ! [OUT]
197 
198  if ( .not. land_phy_update_bottom_water ) then
199  do j = js, je
200  do i = is, ie
201  water1(lke,i,j) = water(lke,i,j)
202  end do
203  end do
204  endif
205 
206  ! runoff of soil moisture (vertical sum)
207  do j = js, je
208  do i = is, ie
209 ! RUNOFF(i,j) = 0.0_RP
210  do k = lks, lke
211 ! RUNOFF(i,j) = RUNOFF(i,j) + max( WATER1(k,i,j) - WaterLimit(i,j), 0.0_RP ) * CDZ(k) * DWATR
212  water1(k,i,j) = min( water1(k,i,j), waterlimit(i,j) )
213  end do
214  end do
215  end do
216 
217  ! estimate thermal diffusivity
218  do j = js, je
219  do i = is, ie
220  do k = lks, lke
221  land_denscs(k,i,j) = ( 1.0_rp - waterlimit(i,j) ) * heatcapacity(i,j) + water_denscs * water1(k,i,j)
222  thermaldiff(k,i,j) = thermalcond(i,j) / land_denscs(k,i,j)
223  end do
224  end do
225  end do
226 
227  ! Solve diffusion of soil temperature (tridiagonal matrix)
228  do j = js, je
229  do i = is, ie
230  l(lks,i,j) = 0.0_rp
231  u(lks,i,j) = -2.0_rp * thermaldiff(lks,i,j) / ( cdz(lks) * ( cdz(lks) + cdz(lks+1) ) ) * dt
232  l(lke,i,j) = -2.0_rp * thermaldiff(lke,i,j) / ( cdz(lke) * ( cdz(lke) + cdz(lke-1) ) ) * dt
233  u(lke,i,j) = 0.0_rp
234 
235  m(lks,i,j) = 1.0_rp - l(lks,i,j) - u(lks,i,j)
236  m(lke,i,j) = 1.0_rp - l(lke,i,j) - u(lke,i,j)
237  end do
238  end do
239 
240  do j = js, je
241  do i = is, ie
242  do k = lks+1, lke-1
243  l(k,i,j) = -2.0_rp * thermaldiff(k,i,j) / ( cdz(k) * ( cdz(k) + cdz(k-1) ) ) * dt
244  u(k,i,j) = -2.0_rp * thermaldiff(k,i,j) / ( cdz(k) * ( cdz(k) + cdz(k+1) ) ) * dt
245  m(k,i,j) = 1.0_rp - l(k,i,j) - u(k,i,j)
246  end do
247  end do
248  end do
249 
250  ! input from atmosphere
251  do j = js, je
252  do i = is, ie
253  v(lks,i,j) = temp(lks,i,j) - sflx_gh(i,j) / ( land_denscs(lks,i,j) * cdz(lks) ) * dt
254  end do
255  end do
256 
257  do j = js, je
258  do i = is, ie
259  do k = lks+1, lke
260  v(k,i,j) = temp(k,i,j)
261  end do
262  end do
263  end do
264 
265  call matrix_solver_tridiagonal( lkmax, & ! [IN]
266  ia, is, ie, & ! [IN]
267  ja, js, je, & ! [IN]
268  u(:,:,:), & ! [IN]
269  m(:,:,:), & ! [IN]
270  l(:,:,:), & ! [IN]
271  v(:,:,:), & ! [IN]
272  temp1(:,:,:) ) ! [OUT]
273 
274  if ( .not. land_phy_update_bottom_temp ) then
275  do j = js, je
276  do i = is, ie
277  temp1(lke,i,j) = temp(lke,i,j)
278  end do
279  end do
280  endif
281 
282  ! calculate tendency
283  do j = js, je
284  do i = is, ie
285  do k = lks, lke
286  if( landuse_fact_land(i,j) > 0.0_rp ) then
287  temp_t(k,i,j) = ( temp1(k,i,j) - temp(k,i,j) ) / dt
288  water_t(k,i,j) = ( water1(k,i,j) - water(k,i,j) ) / dt
289  else
290  temp_t(k,i,j) = 0.0_rp
291  water_t(k,i,j) = 0.0_rp
292  end if
293  end do
294  end do
295  end do
296 
297  return
module LANDUSE
module MATRIX
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
Here is the caller graph for this function: