next up previous contents
Next: Fortran: Module Interface BareGroundFluxesMod Up: Fortran: Module Interface BalanceCheckMod Previous: BeginWaterBalance   Contents

BalanceCheck


INTERFACE:

   subroutine BalanceCheck(lbp, ubp, lbc, ubc, lbl, ubl, lbg, ubg)
DESCRIPTION:

This subroutine accumulates the numerical truncation errors of the water and energy balance calculation. It is helpful to see the performance of the process of integration.

The error for energy balance:

error = abs(Net radiation - change of internal energy - Sensible heat - Latent heat)

The error for water balance:

error = abs(precipitation - change of water storage - evaporation - runoff)


USES:

     use clmtype
     use clm_atmlnd   , only : clm_a2l
     use subgridAveMod
     use clm_time_manager , only : get_step_size, get_nstep
     use clm_varcon   , only : isturb, icol_roof, icol_sunwall, icol_shadewall, &
                               spval, icol_road_perv, icol_road_imperv, istice_mec
     use clm_varctl   , only : glc_dyntopo
ARGUMENTS:
     implicit none
     integer :: lbp, ubp ! pft-index bounds
     integer :: lbc, ubc ! column-index bounds
     integer :: lbl, ubl ! landunit-index bounds
     integer :: lbg, ubg ! grid-index bounds
CALLED FROM:
   subroutine clm_driver
REVISION HISTORY:
   15 September 1999: Yongjiu Dai; Initial code
   15 December 1999:  Paul Houser and Jon Radakovich; F90 Revision
   10 November 2000: Mariana Vertenstein
   Migrated to new data structures by Mariana Vertenstein and
   Peter Thornton
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 :: cgridcell(:)       ! column's gridcell index
     integer , pointer :: ltype(:)           ! landunit type 
     integer , pointer :: ctype(:)           ! column type 
     real(r8), pointer :: pwtgcell(:)        ! pft's weight relative to corresponding gridcell
     real(r8), pointer :: cwtgcell(:)        ! column's weight relative to corresponding gridcell
     real(r8), pointer :: forc_rain(:)       ! rain rate [mm/s]
     real(r8), pointer :: forc_snow(:)       ! snow rate [mm/s]
     real(r8), pointer :: forc_lwrad(:)      ! downward infrared (longwave) radiation (W/m**2)
     real(r8), pointer :: endwb(:)           ! water mass end of the time step
     real(r8), pointer :: begwb(:)           ! water mass begining of the time step
     real(r8), pointer :: fsa(:)             ! solar radiation absorbed (total) (W/m**2)
     real(r8), pointer :: fsr(:)             ! solar radiation reflected (W/m**2)
     real(r8), pointer :: eflx_lwrad_out(:)  ! emitted infrared (longwave) radiation (W/m**2)
     real(r8), pointer :: eflx_lwrad_net(:)  ! net infrared (longwave) rad (W/m**2) [+ = to atm]
     real(r8), pointer :: sabv(:)            ! solar radiation absorbed by vegetation (W/m**2)
     real(r8), pointer :: sabg(:)            ! solar radiation absorbed by ground (W/m**2)
     real(r8), pointer :: eflx_sh_tot(:)     ! total sensible heat flux (W/m**2) [+ to atm]
     real(r8), pointer :: eflx_sh_totg(:)    ! total sensible heat flux at grid level (W/m**2) [+ to atm]
     real(r8), pointer :: eflx_dynbal(:)     ! energy conversion flux due to dynamic land cover change(W/m**2) [+ to atm]
     real(r8), pointer :: eflx_lh_tot(:)     ! total latent heat flux (W/m8*2)  [+ to atm]
     real(r8), pointer :: eflx_soil_grnd(:)  ! soil heat flux (W/m**2) [+ = into soil]
     real(r8), pointer :: qflx_evap_tot(:)   ! qflx_evap_soi + qflx_evap_can + qflx_tran_veg
     real(r8), pointer :: qflx_irrig(:)      ! irrigation flux (mm H2O /s)
     real(r8), pointer :: qflx_surf(:)       ! surface runoff (mm H2O /s)
     real(r8), pointer :: qflx_qrgwl(:)      ! qflx_surf at glaciers, wetlands, lakes
     real(r8), pointer :: qflx_drain(:)      ! sub-surface runoff (mm H2O /s)
     real(r8), pointer :: qflx_runoff(:)     ! total runoff (mm H2O /s)
     real(r8), pointer :: qflx_runoffg(:)    ! total runoff at gridcell level inc land cover change flux (mm H2O /s)
     real(r8), pointer :: qflx_liq_dynbal(:) ! liq runoff due to dynamic land cover change (mm H2O /s)
     real(r8), pointer :: qflx_snwcp_ice(:)  ! excess snowfall due to snow capping (mm H2O /s) [+]`
     real(r8), pointer :: qflx_glcice(:)     ! flux of new glacier ice (mm H2O /s) [+ if ice grows]
     real(r8), pointer :: qflx_snwcp_iceg(:) ! excess snowfall due to snow cap inc land cover change flux (mm H20/s)
     real(r8), pointer :: qflx_ice_dynbal(:) ! ice runoff due to dynamic land cover change (mm H2O /s)
     real(r8), pointer :: forc_solad(:,:)    ! direct beam radiation (vis=forc_sols , nir=forc_soll )
     real(r8), pointer :: forc_solai(:,:)    ! diffuse radiation     (vis=forc_solsd, nir=forc_solld)
     real(r8), pointer :: eflx_traffic_pft(:)    ! traffic sensible heat flux (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 :: canyon_hwr(:)      ! ratio of building height to street width
     real(r8), pointer :: eflx_heat_from_ac_pft(:) !sensible heat flux put back into canyon due to removal by AC (W/m**2)
   local pointers to original implicit out arguments
     real(r8), pointer :: errh2o(:)          ! water conservation error (mm H2O)
     real(r8), pointer :: errsol(:)          ! solar radiation conservation error (W/m**2)
     real(r8), pointer :: errlon(:)          ! longwave radiation conservation error (W/m**2)
     real(r8), pointer :: errseb(:)          ! surface energy conservation error (W/m**2)
     real(r8), pointer :: netrad(:)          ! net radiation (positive downward) (W/m**2)
     real(r8), pointer :: errsoi_col(:)      ! column-level soil/lake energy conservation error (W/m**2)



Erik Kluzek 2011-06-15