next up previous contents
Next: SoilThermProp Up: Fortran: Module Interface SoilTemperatureMod Previous: Fortran: Module Interface SoilTemperatureMod   Contents

SoilTemperature


INTERFACE:

   subroutine SoilTemperature(lbl, ubl, lbc, ubc, num_urbanl, filter_urbanl, &
                              num_nolakec, filter_nolakec, xmf, fact)
DESCRIPTION:

Snow and soil temperatures including phase change o The volumetric heat capacity is calculated as a linear combination in terms of the volumetric fraction of the constituent phases. o The thermal conductivity of soil is computed from the algorithm of Johansen (as reported by Farouki 1981), and the conductivity of snow is from the formulation used in SNTHERM (Jordan 1991). o Boundary conditions: F = Rnet - Hg - LEg (top), F= 0 (base of the soil column). o Soil / snow temperature is predicted from heat conduction in 10 soil layers and up to 5 snow layers. The thermal conductivities at the interfaces between two neighboring layers (j, j+1) are derived from an assumption that the flux across the interface is equal to that from the node j to the interface and the flux from the interface to the node j+1. The equation is solved using the Crank-Nicholson method and results in a tridiagonal system equation.


USES:

     use shr_kind_mod  , only : r8 => shr_kind_r8
     use clmtype
     use clm_atmlnd    , only : clm_a2l
     use clm_time_manager  , only : get_step_size
     use clm_varctl    , only : iulog
     use clm_varcon    , only : sb, capr, cnfac, hvap, istice_mec, isturb, &
                                icol_roof, icol_sunwall, icol_shadewall, &
                                icol_road_perv, icol_road_imperv, istwet
     use clm_varpar    , only : nlevsno, nlevgrnd, max_pft_per_col, nlevurb
     use TridiagonalMod, only : Tridiagonal
ARGUMENTS:
     implicit none
     integer , intent(in)  :: lbc, ubc                    ! column bounds
     integer , intent(in)  :: num_nolakec                 ! number of column non-lake points in column filter
     integer , intent(in)  :: filter_nolakec(ubc-lbc+1)   ! column filter for non-lake points
     integer , intent(in)  :: lbl, ubl                    ! landunit-index bounds
     integer , intent(in)  :: num_urbanl                  ! number of urban landunits in clump
     integer , intent(in)  :: filter_urbanl(ubl-lbl+1)    ! urban landunit filter
     real(r8), intent(out) :: xmf(lbc:ubc)                ! total latent heat of phase change of ground water
     real(r8), intent(out) :: fact(lbc:ubc, -nlevsno+1:nlevgrnd) ! used in computing tridiagonal matrix
CALLED FROM:
   subroutine Biogeophysics2 in module Biogeophysics2Mod
REVISION HISTORY:
   15 September 1999: Yongjiu Dai; Initial code
   15 December 1999:  Paul Houser and Jon Radakovich; F90 Revision
   12/19/01, Peter Thornton
   Changed references for tg to t_grnd, for consistency with the
   rest of the code (tg eliminated as redundant)
   2/14/02, Peter Thornton: Migrated to new data structures. Added pft loop
   in calculation of net ground heat flux.
   3/18/08, David Lawrence: Change nlevsoi to nlevgrnd for deep soil
   03/28/08, Mark Flanner: Changes to allow solar radiative absorption in all snow layers and top soil layer
LOCAL VARIABLES:
   local pointers to original implicit in arguments
     integer , pointer :: pgridcell(:)       ! pft's gridcell index
     integer , pointer :: plandunit(:)       ! pft's landunit index
     integer , pointer :: clandunit(:)       ! column's landunit
     integer , pointer :: ltype(:)           ! landunit type
     integer , pointer :: ctype(:)           ! column type
     integer , pointer :: npfts(:)           ! column's number of pfts 
     integer , pointer :: pfti(:)            ! column's beginning pft index 
     real(r8), pointer :: pwtcol(:)          ! weight of pft relative to column
     real(r8), pointer :: pwtgcell(:)        ! weight of pft relative to corresponding gridcell
     real(r8), pointer :: forc_lwrad(:)      ! downward infrared (longwave) radiation (W/m**2)
     integer , pointer :: snl(:)             ! number of snow layers
     real(r8), pointer :: htvp(:)            ! latent heat of vapor of water (or sublimation) [j/kg]
     real(r8), pointer :: emg(:)             ! ground emissivity
     real(r8), pointer :: cgrnd(:)           ! deriv. of soil energy flux wrt to soil temp [w/m2/k]
     real(r8), pointer :: dlrad(:)           ! downward longwave radiation blow the canopy [W/m2]
     real(r8), pointer :: sabg(:)            ! solar radiation absorbed by ground (W/m**2)
     integer , pointer :: frac_veg_nosno(:)  ! fraction of vegetation not covered by snow (0 OR 1 now) [-] (new)
     real(r8), pointer :: eflx_sh_grnd(:)    ! sensible heat flux from ground (W/m**2) [+ to atm]
     real(r8), pointer :: qflx_evap_soi(:)   ! soil evaporation (mm H2O/s) (+ = to atm)
     real(r8), pointer :: qflx_tran_veg(:)   ! vegetation transpiration (mm H2O/s) (+ = to atm)
     real(r8), pointer :: zi(:,:)            ! interface level below a "z" level (m)
     real(r8), pointer :: dz(:,:)            ! layer depth (m)
     real(r8), pointer :: z(:,:)             ! layer thickness (m)
     real(r8), pointer :: t_soisno(:,:)      ! soil temperature (Kelvin)
     real(r8), pointer :: eflx_lwrad_net(:)  ! net infrared (longwave) rad (W/m**2) [+ = to atm]
     real(r8), pointer :: tssbef(:,:)        ! temperature at previous time step [K]
     real(r8), pointer :: t_building(:)      ! internal building temperature (K)
     real(r8), pointer :: t_building_max(:)  ! maximum internal building temperature (K)
     real(r8), pointer :: t_building_min(:)  ! minimum internal building temperature (K)
     real(r8), pointer :: hc_soi(:)          ! soil heat content (MJ/m2)
     real(r8), pointer :: hc_soisno(:)       ! soil plus snow plus lake heat content (MJ/m2)
     real(r8), pointer :: eflx_fgr12(:)      ! heat flux between soil layer 1 and 2 (W/m2)
     real(r8), pointer :: eflx_traffic(:)    ! traffic sensible heat flux (W/m**2)
     real(r8), pointer :: eflx_wasteheat(:)  ! sensible heat flux from urban heating/cooling sources of waste heat (W/m**2)
     real(r8), pointer :: eflx_wasteheat_pft(:)  ! sensible heat flux from urban heating/cooling sources of waste heat (W/m**2)
     real(r8), pointer :: eflx_heat_from_ac(:) !sensible heat flux put back into canyon due to removal by AC (W/m**2)
     real(r8), pointer :: eflx_heat_from_ac_pft(:) !sensible heat flux put back into canyon due to removal by AC (W/m**2)
     real(r8), pointer :: eflx_traffic_pft(:)    ! traffic sensible heat flux (W/m**2)
     real(r8), pointer :: eflx_anthro(:)         ! total anthropogenic heat flux (W/m**2)
     real(r8), pointer :: canyon_hwr(:)      ! urban canyon height to width ratio
     real(r8), pointer :: wtlunit_roof(:)    ! weight of roof with respect to landunit
     real(r8), pointer :: eflx_bot(:)        ! heat flux from beneath column (W/m**2) [+ = upward]
   
   local pointers to  original implicit inout arguments
     real(r8), pointer :: t_grnd(:)          ! ground surface temperature [K]
   local pointers to original implicit out arguments
     real(r8), pointer :: eflx_gnet(:)          ! net ground heat flux into the surface (W/m**2)
     real(r8), pointer :: dgnetdT(:)            ! temperature derivative of ground net heat flux
     real(r8), pointer :: eflx_building_heat(:) ! heat flux from urban building interior to walls, roof (W/m**2)
 
   variables needed for SNICAR
     real(r8), pointer :: sabg_lyr(:,:)      ! absorbed solar radiation (pft,lyr) [W/m2]
     real(r8), pointer :: h2osno(:)          ! total snow water (col) [kg/m2]
     real(r8), pointer :: h2osoi_liq(:,:)    ! liquid water (col,lyr) [kg/m2]
     real(r8), pointer :: h2osoi_ice(:,:)    ! ice content (col,lyr) [kg/m2]
 
   Urban building HAC fluxes
     real(r8), pointer :: eflx_urban_ac(:)      ! urban air conditioning flux (W/m**2)
     real(r8), pointer :: eflx_urban_heat(:)    ! urban heating flux (W/m**2)
   !OTHER LOCAL VARIABLES:



Erik Kluzek 2011-06-15