43   logical, 
private :: land_phy_slab_const = .false. 
    45   logical, 
private :: land_phy_update_bottom_temp  = .false. 
    46   logical, 
private :: land_phy_update_bottom_water = .false. 
    48   real(RP), 
private :: water_denscs
    50   logical, 
allocatable, 
private :: is_lnd(:,:)
    66     character(len=*), 
intent(in) :: LAND_TYPE
    68     namelist / param_land_phy_slab / &
    69        land_phy_update_bottom_temp,   &
    70        land_phy_update_bottom_water
    77     if( 
io_l ) 
write(
io_fid_log,*) 
'++++++ Module[SLAB] / Categ[LAND PHY] / Origin[SCALElib]'    81     read(
io_fid_conf,nml=param_land_phy_slab,iostat=ierr)
    83        if( 
io_l ) 
write(
io_fid_log,*) 
'*** Not found namelist. Default used.'    84     elseif( ierr > 0 ) 
then     85        write(*,*) 
'xxx Not appropriate names in namelist PARAM_LAND_PHY_SLAB. Check!'    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.
    95        write(*,*) 
'xxx wrong LAND_TYPE. Check!'    99     water_denscs = dwatr * cl
   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
   106     allocate( is_lnd(
ia,
ja) )
   113         is_lnd(i,j) = .false.
   141        matrix_solver_tridiagonal
   145     real(RP), 
intent(out) :: TEMP_t      (
lkmax,
ia,
ja)
   146     real(RP), 
intent(out) :: WATER_t     (
lkmax,
ia,
ja)
   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
   164     real(RP) :: SOIL_DENSCS(
ia,
ja)
   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)
   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)
   198        vin(
lks  ,i,j) = water(
lks,i,j) &
   199                       + ( sflx_prec(i,j) - sflx_evap(i,j) ) / ( cdz(
lks) * dwatr ) * dt 
   202           vin(k,i,j) = water(k,i,j)
   205        vin(
lke-1,i,j) = water(
lke-1,i,j) &
   206                       - u(
lke-1,i,j) * water(
lke,i,j)
   210     call matrix_solver_tridiagonal( 
lkmax-1,     & 
   222        water1(k,i,j) = vout(k,i,j)
   228     if ( land_phy_update_bottom_water ) 
then   231           water1(
lke,i,j) = vout(
lke-1,i,j)
   237           water1(
lke,i,j) = water(
lke,i,j)
   252           water1(k,i,j) = min( water1(k,i,j), waterlimit(i,j) )
   262        soil_denscs(i,j) = ( 1.0_rp - waterlimit(i,j) ) * heatcapacity(i,j)
   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)
   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)
   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 
   293           vin(k,i,j) = temp(k,i,j)
   296        vin(
lke-1,i,j) = temp(
lke-1,i,j) &
   297                       - u(
lke-1,i,j) * temp(
lke,i,j)
   301     call matrix_solver_tridiagonal( 
lkmax-1,     & 
   313        temp1(k,i,j) = vout(k,i,j)
   319     if ( land_phy_update_bottom_temp ) 
then   322           temp1(
lke,i,j) = vout(
lke-1,i,j)
   328           temp1(
lke,i,j) = temp(
lke,i,j)
   335     if( land_phy_slab_const ) 
then   339           temp_t(k,i,j) = 0.0_rp
   340           water_t(k,i,j) = 0.0_rp
   347           if( is_lnd(i,j) ) 
then   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
   354                 temp_t(k,i,j) = 0.0_rp
   355                 water_t(k,i,j) = 0.0_rp
 integer, public is
start point of inner domain: x, local 
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] 
logical, public io_l
output log or not? (this process) 
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K] 
subroutine, public land_phy_slab_setup(LAND_TYPE)
Setup. 
integer, public ia
of x whole cells (local, with HALO)
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 
integer, public ie
end point of inner domain: x, local 
logical, public io_lnml
output log or not? (for namelist, this process) 
module LAND / Physics Slab model 
integer, public io_fid_conf
Config file ID. 
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor 
integer, public io_fid_log
Log file ID. 
integer, public ja
of y whole cells (local, with HALO)