SCALE-RM
scale_land_phy_slab.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
11  !-----------------------------------------------------------------------------
12  !
13  !++ used modules
14  !
15  use scale_precision
16  use scale_stdio
17  use scale_prof
18  use scale_debug
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: land_phy_slab_setup
29  public :: land_phy_slab
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  !-----------------------------------------------------------------------------
36  !
37  !++ Private procedure
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private parameters & variables
42  !
43  logical, private :: land_phy_slab_const = .false. ! constant condition?
44 
45  logical, private :: land_phy_update_bottom_temp = .false. ! Is LAND_TEMP updated in the lowest level?
46  logical, private :: land_phy_update_bottom_water = .false. ! Is LAND_WATER updated in the lowest level?
47 
48  real(RP), private :: water_denscs
49 
50  logical, allocatable, private :: is_lnd(:,:)
51 
52  !-----------------------------------------------------------------------------
53 contains
54  !-----------------------------------------------------------------------------
56  subroutine land_phy_slab_setup( LAND_TYPE )
57  use scale_process, only: &
59  use scale_const, only: &
60  dwatr => const_dwatr, &
61  cl => const_cl
62  use scale_landuse, only: &
64  implicit none
65 
66  character(len=*), intent(in) :: LAND_TYPE
67 
68  namelist / param_land_phy_slab / &
69  land_phy_update_bottom_temp, &
70  land_phy_update_bottom_water
71 
72  integer :: i, j
73  integer :: ierr
74  !---------------------------------------------------------------------------
75 
76  if( io_l ) write(io_fid_log,*)
77  if( io_l ) write(io_fid_log,*) '++++++ Module[SLAB] / Categ[LAND PHY] / Origin[SCALElib]'
78 
79  !--- read namelist
80  rewind(io_fid_conf)
81  read(io_fid_conf,nml=param_land_phy_slab,iostat=ierr)
82  if( ierr < 0 ) then !--- missing
83  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
84  elseif( ierr > 0 ) then !--- fatal error
85  write(*,*) 'xxx Not appropriate names in namelist PARAM_LAND_PHY_SLAB. Check!'
86  call prc_mpistop
87  endif
88  if( io_lnml ) write(io_fid_log,nml=param_land_phy_slab)
89 
90  if( land_type == 'CONST' ) then
91  land_phy_slab_const = .true.
92  else if( land_type == 'SLAB' ) then
93  land_phy_slab_const = .false.
94  else
95  write(*,*) 'xxx wrong LAND_TYPE. Check!'
96  call prc_mpistop
97  end if
98 
99  water_denscs = dwatr * cl
100 
101  if( io_l ) write(io_fid_log,*)
102  if( io_l ) write(io_fid_log,*) '*** Update soil temperature of bottom layer? : ', land_phy_update_bottom_temp
103  if( io_l ) write(io_fid_log,*) '*** Update soil moisture of bottom layer? : ', land_phy_update_bottom_water
104 
105  ! judge to run slab land model
106  allocate( is_lnd(ia,ja) )
107 
108  do j = js, je
109  do i = is, ie
110  if( landuse_fact_land(i,j) > 0.0_rp ) then
111  is_lnd(i,j) = .true.
112  else
113  is_lnd(i,j) = .false.
114  end if
115  end do
116  end do
117 
118  return
119  end subroutine land_phy_slab_setup
120 
121  !-----------------------------------------------------------------------------
123  subroutine land_phy_slab( &
124  TEMP_t, &
125  WATER_t, &
126  TEMP, &
127  WATER, &
128  WaterLimit, &
129  ThermalCond, &
130  HeatCapacity, &
131  WaterDiff, &
132  SFLX_GH, &
133  SFLX_prec, &
134  SFLX_evap, &
135  CDZ, &
136  dt )
138  use scale_const, only: &
139  dwatr => const_dwatr
140  use scale_matrix, only: &
141  matrix_solver_tridiagonal
142  implicit none
143 
144  ! arguments
145  real(RP), intent(out) :: TEMP_t (lkmax,ia,ja)
146  real(RP), intent(out) :: WATER_t (lkmax,ia,ja)
147 
148  real(RP), intent(in) :: TEMP (lkmax,ia,ja)
149  real(RP), intent(in) :: WATER (lkmax,ia,ja)
150  real(RP), intent(in) :: WaterLimit (ia,ja)
151  real(RP), intent(in) :: ThermalCond (ia,ja)
152  real(RP), intent(in) :: HeatCapacity(ia,ja)
153  real(RP), intent(in) :: WaterDiff (ia,ja)
154  real(RP), intent(in) :: SFLX_GH (ia,ja)
155  real(RP), intent(in) :: SFLX_prec (ia,ja)
156  real(RP), intent(in) :: SFLX_evap (ia,ja)
157  real(RP), intent(in) :: CDZ (lkmax)
158  real(DP), intent(in) :: dt
159 
160  ! work
161  real(RP) :: TEMP1 (lkmax,ia,ja)
162  real(RP) :: WATER1(lkmax,ia,ja)
163 ! real(RP) :: RUNOFF(IA,JA)
164  real(RP) :: SOIL_DENSCS(ia,ja)
165 
166  real(RP) :: U (lkmax-1,ia,ja)
167  real(RP) :: M (lkmax-1,ia,ja)
168  real(RP) :: L (lkmax-1,ia,ja)
169  real(RP) :: Vin (lkmax-1,ia,ja)
170  real(RP) :: Vout(lkmax-1,ia,ja)
171 
172  integer :: k, i, j
173  !---------------------------------------------------------------------------
174 
175  if( io_l ) write(io_fid_log,*) '*** Land step: Slab'
176 
177  ! Solve diffusion of soil moisture (tridiagonal matrix)
178  do j = js, je
179  do i = is, ie
180  l(lks,i,j) = 0.0_rp
181  u(lks,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(lks) * ( cdz(lks) + cdz(lks+1) ) ) * dt
182  m(lks,i,j) = 1.0_rp - l(lks,i,j) - u(lks,i,j)
183  enddo
184  enddo
185 
186  do j = js, je
187  do i = is, ie
188  do k = lks+1, lke-1
189  l(k,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(k) * ( cdz(k) + cdz(k-1) ) ) * dt
190  u(k,i,j) = -2.0_rp * waterdiff(i,j) / ( cdz(k) * ( cdz(k) + cdz(k+1) ) ) * dt
191  m(k,i,j) = 1.0_rp - l(k,i,j) - u(k,i,j)
192  enddo
193  enddo
194  enddo
195 
196  do j = js, je
197  do i = is, ie
198  vin(lks ,i,j) = water(lks,i,j) &
199  + ( sflx_prec(i,j) - sflx_evap(i,j) ) / ( cdz(lks) * dwatr ) * dt ! input from atmosphere
200 
201  do k = lks+1, lke-2
202  vin(k,i,j) = water(k,i,j)
203  enddo
204 
205  vin(lke-1,i,j) = water(lke-1,i,j) &
206  - u(lke-1,i,j) * water(lke,i,j)
207  enddo
208  enddo
209 
210  call matrix_solver_tridiagonal( lkmax-1, & ! [IN]
211  ia, is, ie, & ! [IN]
212  ja, js, je, & ! [IN]
213  u(:,:,:), & ! [IN]
214  m(:,:,:), & ! [IN]
215  l(:,:,:), & ! [IN]
216  vin(:,:,:), & ! [IN]
217  vout(:,:,:) ) ! [OUT]
218 
219  do j = js, je
220  do i = is, ie
221  do k = lks, lke-1
222  water1(k,i,j) = vout(k,i,j)
223  enddo
224  enddo
225  enddo
226 
227  ! lowest layer treatment
228  if ( land_phy_update_bottom_water ) then
229  do j = js, je
230  do i = is, ie
231  water1(lke,i,j) = vout(lke-1,i,j)
232  enddo
233  enddo
234  else
235  do j = js, je
236  do i = is, ie
237  water1(lke,i,j) = water(lke,i,j)
238  enddo
239  enddo
240  endif
241 
242 
243 
244  ! runoff of soil moisture (vertical sum)
245  do j = js, je
246  do i = is, ie
247 ! RUNOFF(i,j) = 0.0_RP
248  do k = lks, lke
249 ! RUNOFF(i,j) = RUNOFF(i,j) &
250 ! + max( WATER1(k,i,j) - WaterLimit(i,j), 0.0_RP ) * CDZ(k) * DWATR
251 
252  water1(k,i,j) = min( water1(k,i,j), waterlimit(i,j) )
253  enddo
254  enddo
255  enddo
256 
257 
258 
259  ! Solve diffusion of soil temperature (tridiagonal matrix)
260  do j = js, je
261  do i = is, ie
262  soil_denscs(i,j) = ( 1.0_rp - waterlimit(i,j) ) * heatcapacity(i,j)
263  enddo
264  enddo
265 
266  do j = js, je
267  do i = is, ie
268  l(lks,i,j) = 0.0_rp
269  u(lks,i,j) = -2.0_rp * thermalcond(i,j) / ( soil_denscs(i,j) + water(lks,i,j)*water_denscs ) &
270  / ( cdz(lks) * ( cdz(lks) + cdz(lks+1) ) ) * dt
271  m(lks,i,j) = 1.0_rp - l(lks,i,j) - u(lks,i,j)
272  enddo
273  enddo
274 
275  do j = js, je
276  do i = is, ie
277  do k = lks+1, lke-1
278  l(k,i,j) = -2.0_rp * thermalcond(i,j) / ( soil_denscs(i,j) + water(k,i,j)*water_denscs ) &
279  / ( cdz(k) * ( cdz(k) + cdz(k-1) ) ) * dt
280  u(k,i,j) = -2.0_rp * thermalcond(i,j) / ( soil_denscs(i,j) + water(k,i,j)*water_denscs ) &
281  / ( cdz(k) * ( cdz(k) + cdz(k+1) ) ) * dt
282  m(k,i,j) = 1.0_rp - l(k,i,j) - u(k,i,j)
283  enddo
284  enddo
285  enddo
286 
287  do j = js, je
288  do i = is, ie
289  vin(lks ,i,j) = temp(lks,i,j) &
290  - sflx_gh(i,j) / ( soil_denscs(i,j) + water(lks,i,j)*water_denscs ) / cdz(lks) * dt ! input from atmosphere
291 
292  do k = lks+1, lke-2
293  vin(k,i,j) = temp(k,i,j)
294  enddo
295 
296  vin(lke-1,i,j) = temp(lke-1,i,j) &
297  - u(lke-1,i,j) * temp(lke,i,j)
298  enddo
299  enddo
300 
301  call matrix_solver_tridiagonal( lkmax-1, & ! [IN]
302  ia, is, ie, & ! [IN]
303  ja, js, je, & ! [IN]
304  u(:,:,:), & ! [IN]
305  m(:,:,:), & ! [IN]
306  l(:,:,:), & ! [IN]
307  vin(:,:,:), & ! [IN]
308  vout(:,:,:) ) ! [OUT]
309 
310  do j = js, je
311  do i = is, ie
312  do k = lks, lke-1
313  temp1(k,i,j) = vout(k,i,j)
314  enddo
315  enddo
316  enddo
317 
318  ! lowest layer treatment
319  if ( land_phy_update_bottom_temp ) then
320  do j = js, je
321  do i = is, ie
322  temp1(lke,i,j) = vout(lke-1,i,j)
323  enddo
324  enddo
325  else
326  do j = js, je
327  do i = is, ie
328  temp1(lke,i,j) = temp(lke,i,j)
329  enddo
330  enddo
331  endif
332 
333 
334  ! calculate tendency
335  if( land_phy_slab_const ) then
336  do j = js, je
337  do i = is, ie
338  do k = lks, lke
339  temp_t(k,i,j) = 0.0_rp
340  water_t(k,i,j) = 0.0_rp
341  enddo
342  enddo
343  enddo
344  else
345  do j = js, je
346  do i = is, ie
347  if( is_lnd(i,j) ) then
348  do k = lks, lke
349  temp_t(k,i,j) = ( temp1(k,i,j) - temp(k,i,j) ) / dt
350  water_t(k,i,j) = ( water1(k,i,j) - water(k,i,j) ) / dt
351  enddo
352  else
353  do k = lks, lke
354  temp_t(k,i,j) = 0.0_rp
355  water_t(k,i,j) = 0.0_rp
356  enddo
357  endif
358  enddo
359  enddo
360  endif
361 
362  return
363  end subroutine land_phy_slab
364 
365 end module scale_land_phy_slab
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:87
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
module STDIO
Definition: scale_stdio.F90:12
subroutine, public land_phy_slab_setup(LAND_TYPE)
Setup.
module grid index
integer, public ia
of x whole cells (local, with HALO)
module LANDUSE
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.
integer, public js
start point of inner domain: y, local
module MATRIX
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
module LAND / Physics Slab model
module land grid index
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)