next up previous contents index
Next: Stomata Up: Module CanopyFluxesMod (File: CanopyFluxesMod.F90) Previous: Module CanopyFluxesMod (File: CanopyFluxesMod.F90)   Contents   Index


CanopyFluxes


INTERFACE:

   subroutine CanopyFluxes(lbg, ubg, lbc, ubc, lbp, ubp, &
                           num_nolakep, filter_nolakep)
DESCRIPTION:

   1. Calculates the leaf temperature:
   2. Calculates the leaf fluxes, transpiration, photosynthesis and
      updates the dew accumulation due to evaporation.
   Method:
   Use the Newton-Raphson iteration to solve for the foliage
   temperature that balances the surface energy budget:
   f(t_veg) = Net radiation - Sensible - Latent = 0
   f(t_veg) + d(f)/d(t_veg) * dt_veg = 0     (*)
   Note:
   (1) In solving for t_veg, t_grnd is given from the previous timestep.
   (2) The partial derivatives of aerodynamical resistances, which cannot
       be determined analytically, are ignored for d(H)/dT and d(LE)/dT
   (3) The weighted stomatal resistance of sunlit and shaded foliage is used
   (4) Canopy air temperature and humidity are derived from => Hc + Hg = Ha
                                                            => Ec + Eg = Ea
   (5) Energy loss is due to: numerical truncation of energy budget equation
       (*); and "ecidif" (see the code) which is dropped into the sensible
       heat
   (6) The convergence criteria: the difference, del = t_veg(n+1)-t_veg(n)
       and del2 = t_veg(n)-t_veg(n-1) less than 0.01 K, and the difference
       of water flux from the leaf between the iteration step (n+1) and (n)
       less than 0.1 W/m2; or the iterative steps over 40.
USES:
     use shr_kind_mod       , only: r8 => shr_kind_r8
     use clmtype
     use time_manager       , only : get_step_size
     use clm_varpar         , only : nlevsoi, nlevsno
     use clm_varcon         , only : sb, cpair, hvap, vkc, grav, denice, &
                                     denh2o, tfrz, csoilc
     use QSatMod            , only : QSat
     use FrictionVelocityMod, only : FrictionVelocity, MoninObukIni
ARGUMENTS:
     implicit none
     integer, intent(in) :: lbg, ubg                    ! gridcell bounds
     integer, intent(in) :: lbc, ubc                    ! column bounds
     integer, intent(in) :: lbp, ubp                    ! pft bounds
     integer, intent(in) :: num_nolakep                 ! number of column non
     integer, intent(in) :: filter_nolakep(ubp-lbp+1)   ! pft filter for non-l
CALLED FROM:
   subroutine Biogeophysics1 in module Biogeophysics1Mod
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 tg to t_grnd for consistency with other routines
   1/29/02, Peter Thornton
   Migrate to new data structures, new calling protocol. For now co2 and
   o2 partial pressures are hardwired, but they should be coming in from
   forc_co2 and forc_o2. Keeping the same hardwired values as in CLM2 to
   assure bit-for-bit results in the first comparisons.
LOCAL VARIABLES:
   local pointers to implicit in variables
    integer , pointer :: ivt(:)         ! pft vegetation type
    integer , pointer :: pcolumn(:)     ! pft's column index
    integer , pointer :: plandunit(:)   ! pft's landunit index
    integer , pointer :: pgridcell(:)   ! pft's gridcell index
    real(r8), pointer :: forc_th(:)     ! atmospheric potential temperature (K
    real(r8), pointer :: t_grnd(:)      ! ground surface temperature [K]
    real(r8), pointer :: thm(:)         ! intermediate variable (forc_t+0.0098
    real(r8), pointer :: qg(:)          ! specific humidity at ground surface 
    real(r8), pointer :: thv(:)         ! virtual potential temperature (kelvi
    real(r8), pointer :: z0mv(:)        ! roughness length over vegetation, mo
    real(r8), pointer :: z0hv(:)        ! roughness length over vegetation, se
    real(r8), pointer :: z0qv(:)        ! roughness length over vegetation, la
    real(r8), pointer :: z0mg(:)        ! roughness length of ground, momentum
    real(r8), pointer :: dqgdT(:)       ! temperature derivative of "qg"
    real(r8), pointer :: htvp(:)        ! latent heat of evaporation (/sublima
    real(r8), pointer :: emv(:)         ! ground emissivity
    real(r8), pointer :: emg(:)         ! vegetation emissivity
    real(r8), pointer :: forc_pbot(:)   ! atmospheric pressure (Pa)
    real(r8), pointer :: forc_q(:)      ! atmospheric specific humidity (kg/kg
    real(r8), pointer :: forc_u(:)      ! atmospheric wind speed in east direc
    real(r8), pointer :: forc_v(:)      ! atmospheric wind speed in north dire
    real(r8), pointer :: forc_hgt_u(:)  ! observational height of wind [m]
    real(r8), pointer :: forc_rho(:)    ! density (kg/m**3)
  ! real(r8), pointer :: forc_co2(:)    ! atmospheric CO2 concentration (Pa)
  ! real(r8), pointer :: forc_o2(:)     ! atmospheric O2 concentration (Pa)
    real(r8), pointer :: forc_lwrad(:)  ! downward infrared (longwave) radiati
    real(r8), pointer :: displa(:)      ! displacement height (m)
    real(r8), pointer :: parsun(:)      ! average absorbed PAR for sunlit leav
    real(r8), pointer :: parsha(:)      ! average absorbed PAR for shaded leav
    real(r8), pointer :: elai(:)        ! one-sided leaf area index with buryi
    real(r8), pointer :: esai(:)        ! one-sided stem area index with buryi
    real(r8), pointer :: fdry(:)        ! fraction of foliage that is green an
    real(r8), pointer :: fwet(:)        ! fraction of canopy that is wet (0 to
    real(r8), pointer :: laisun(:)      ! sunlit leaf area
    real(r8), pointer :: laisha(:)      ! shaded leaf area
    integer , pointer :: frac_veg_nosno(:)  ! frac of veg not covered by snow 
    real(r8), pointer :: sabv(:)        ! solar radiation absorbed by vegetati
    real(r8), pointer :: watsat(:,:)    ! volumetric soil water at saturation 
    real(r8), pointer :: h2osoi_ice(:,:)! ice lens (kg/m2)
    real(r8), pointer :: h2osoi_liq(:,:)! liquid water (kg/m2)
    real(r8), pointer :: dz(:,:)        ! layer depth (m)
    real(r8), pointer :: t_soisno(:,:)  ! soil temperature (Kelvin)
    real(r8), pointer :: sucsat(:,:)    ! minimum soil suction (mm)
    real(r8), pointer :: bsw(:,:)       ! Clapp and Hornberger "b"
    real(r8), pointer :: rootfr(:,:)    ! fraction of roots in each soil layer
    real(r8), pointer :: smpmax(:)      ! wilting point potential in mm
    real(r8), pointer :: dleaf(:)       ! characteristic leaf dimension (m)
   local pointers to implicit inout arguments
    real(r8), pointer :: cgrnds(:)      ! deriv. of soil sensible heat flux wr
    real(r8), pointer :: cgrndl(:)      ! deriv. of soil latent heat flux wrt 
    real(r8), pointer :: t_veg(:)       ! vegetation temperature (Kelvin)
    real(r8), pointer :: t_ref2m(:)     ! 2 m height surface air temperature (
    real(r8), pointer :: q_ref2m(:)     ! 2 m height surface specific humidity
    real(r8), pointer :: h2ocan(:)      ! canopy water (mm H2O)
 #ifdef DGVM
    real(r8), pointer :: annpsnpot(:)   ! annual potential photosynthesis (umo
    real(r8), pointer :: annpsn(:)      ! annual photosynthesis (umol CO2 /m**
 #else
    real(r8) :: annpsnpot(lbp:ubp)      ! temporary, set to zer0
    real(r8) :: annpsn(lbp:ubp)         ! temporary, set to zer0
 #endif
   local pointers to implicit out arguments
    real(r8), pointer :: cgrnd(:)           ! deriv. of soil energy flux wrt t
    real(r8), pointer :: dlrad(:)           ! downward longwave radiation belo
    real(r8), pointer :: ulrad(:)           ! upward longwave radiation above 
    real(r8), pointer :: ram1(:)            ! aerodynamical resistance (s/m)
    real(r8), pointer :: btran(:)           ! transpiration wetness factor (0 
    real(r8), pointer :: rssun(:)           ! sunlit stomatal resistance (s/m)
    real(r8), pointer :: rssha(:)           ! shaded stomatal resistance (s/m)
    real(r8), pointer :: psnsun(:)          ! sunlit leaf photosynthesis (umol
    real(r8), pointer :: psnsha(:)          ! shaded leaf photosynthesis (umol
    real(r8), pointer :: qflx_tran_veg(:)   ! vegetation transpiration (mm H2O
    real(r8), pointer :: dt_veg(:)          ! change in t_veg, last iteration 
    real(r8), pointer :: qflx_evap_veg(:)   ! vegetation evaporation (mm H2O/s
    real(r8), pointer :: eflx_sh_veg(:)     ! sensible heat flux from leaves (
    real(r8), pointer :: taux(:)            ! wind (shear) stress: e-w (kg/m/s
    real(r8), pointer :: tauy(:)            ! wind (shear) stress: n-s (kg/m/s
    real(r8), pointer :: eflx_sh_grnd(:)    ! sensible heat flux from ground (
    real(r8), pointer :: qflx_evap_soi(:)   ! soil evaporation (mm H2O/s) (+ =
    real(r8), pointer :: fpsn(:)            ! photosynthesis (umol CO2 /m**2 /
    real(r8), pointer :: rootr(:,:)         ! effective fraction of roots in e



Mariana Vertenstein 2004-06-21