5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_mpr_pet Module Reference

MPR routine for PET. More...

Functions/Subroutines

subroutine, public pet_correctbylai (param, nodata, lcover0, lai0, mask0, cell_id0, upp_row_l1, low_row_l1, lef_col_l1, rig_col_l1, nl0_in_l1, l1_petlaicorfactor)
 estimate PET correction factor based on LAI at L1
 
subroutine, public pet_correctbyasp (id0, latitude_l0, asp0, param, nodata, fasp0)
 correction of PET
 
subroutine, public priestley_taylor_alpha (lai0, param, mask0, nodata, cell_id0, nl0_in_l1, upp_row_l1, low_row_l1, lef_col_l1, rig_col_l1, priestley_taylor_alpha1)
 Regionalization of priestley taylor alpha.
 
subroutine, public bulksurface_resistance (lai0, param, mask0, nodata, cell_id0, nl0_in_l1, upp_row_l1, low_row_l1, lef_col_l1, rig_col_l1, bulksurface_resistance1)
 Regionalization of bulk surface resistance.
 

Detailed Description

MPR routine for PET.

This module sets up pet correction factor at level-1 based on LAI

Authors
Mehmet Cuneyd Demirel, Simon Stisen
Date
May 2017

Function/Subroutine Documentation

◆ bulksurface_resistance()

subroutine, public mo_mpr_pet::bulksurface_resistance ( real(dp), dimension(:, :), intent(in)  lai0,
real(dp), intent(in)  param,
logical, dimension(:, :), intent(in)  mask0,
real(dp), intent(in)  nodata,
integer(i4), dimension(:), intent(in)  cell_id0,
integer(i4), dimension(:), intent(in)  nl0_in_l1,
integer(i4), dimension(:), intent(in)  upp_row_l1,
integer(i4), dimension(:), intent(in)  low_row_l1,
integer(i4), dimension(:), intent(in)  lef_col_l1,
integer(i4), dimension(:), intent(in)  rig_col_l1,
real(dp), dimension(:, :), intent(out)  bulksurface_resistance1 
)

Regionalization of bulk surface resistance.

estimation of bulk surface resistance Global parameters needed (see mhm_parameter.nml):

  • param(1) = stomatal_resistance
    Parameters
    [in]real(dp), dimension(:, :) :: LAI0LAI at level-0
    [in]real(dp) :: param- global parameter
    [in]logical, dimension(:, :) :: mask0mask at level 0
    [in]real(dp) :: nodata- nodata value
    [in]integer(i4), dimension(:) :: cell_id0Cell ids of hi res field
    [in]integer(i4), dimension(:) :: nL0_in_L1number of l0 cells within a l1 cell
    [in]integer(i4), dimension(:) :: Upp_row_L1upper row of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Low_row_L1lower row of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Lef_col_L1left col of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Rig_col_L1right col of a l1 cell in l0 grid
    [out]real(dp), dimension(:, :) :: bulksurface_resistance1bulk surface resistance
    Authors
    Matthias Zink
    Date
    Apr 2013

Definition at line 407 of file mo_mpr_pet.f90.

References mo_mpr_constants::lai_factor_surfresi, mo_mpr_constants::lai_offset_surfresi, mo_mpr_constants::max_surfresist, and mo_upscaling_operators::upscale_arithmetic_mean().

Referenced by mo_multi_param_reg::mpr().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pet_correctbyasp()

subroutine, public mo_mpr_pet::pet_correctbyasp ( integer(i4), dimension(:), intent(in)  id0,
real(dp), dimension(:), intent(in)  latitude_l0,
real(dp), dimension(:), intent(in)  asp0,
real(dp), dimension(3), intent(in)  param,
real(dp), intent(in)  nodata,
real(dp), dimension(:), intent(out)  fasp0 
)

correction of PET

Correction of PET based on L0 aspect data. Global parameters needed (see mhm_parameter.nml):

  • param(1) = minCorrectionFactorPET
  • param(2) = maxCorrectionFactorPET
  • param(3) = aspectTresholdPET
    Parameters
    [in]integer(i4), dimension(:) :: id0Level 0 cell id
    [in]real(dp), dimension(:) :: latitude_l0latitude on l0
    [in]real(dp), dimension(:) :: Asp0[degree] Aspect at Level 0
    [in]real(dp), dimension(3) :: paramprocess parameters
    [in]real(dp) :: nodata- no data value
    [out]real(dp), dimension(:) :: fAsp0PET correction for Aspect
    Authors
    Stephan Thober, Rohini Kumar
    Date
    Dec 2012

Definition at line 214 of file mo_mpr_pet.f90.

Referenced by mo_multi_param_reg::mpr().

Here is the caller graph for this function:

◆ pet_correctbylai()

subroutine, public mo_mpr_pet::pet_correctbylai ( real(dp), dimension(5), intent(in)  param,
real(dp), intent(in)  nodata,
integer(i4), dimension(:), intent(in)  lcover0,
real(dp), dimension(:, :), intent(in)  lai0,
logical, dimension(:, :), intent(in)  mask0,
integer(i4), dimension(:), intent(in)  cell_id0,
integer(i4), dimension(:), intent(in)  upp_row_l1,
integer(i4), dimension(:), intent(in)  low_row_l1,
integer(i4), dimension(:), intent(in)  lef_col_l1,
integer(i4), dimension(:), intent(in)  rig_col_l1,
integer(i4), dimension(:), intent(in)  nl0_in_l1,
real(dp), dimension(:, :), intent(inout)  l1_petlaicorfactor 
)

estimate PET correction factor based on LAI at L1

estimate PET correction factor based on LAI at L1 for a given Leaf Area Index field. Global parameters needed (see mhm_parameter.nml): Process Case 5:

  • param(1) = PET_a_forest
  • param(2) = PET_a_impervious
  • param(3) = PET_a_pervious
  • param(4) = PET_b
  • param(5) = PET_c Example DSF=PET_a+PET_b*(1-exp(PET_c*LAI)) Similar to the crop coefficient concept Kc=a+b*(1-exp(c*LAI)) by Allen, R. G., L. S. Pereira, D. Raes, and M. Smith (1998), Crop evapotranspiration - Guidelines for computing crop water requirements., FAO Irrigation and drainage paper 56. See Chapter 9, Equation 97 http://www.fao.org/docrep/X0490E/x0490e0f.htm Date: 17/5/2017
    Parameters
    [in]real(dp), dimension(5) :: paramparameters
    [in]real(dp) :: nodata- nodata value
    [in]integer(i4), dimension(:) :: LCOVER0Land cover at level 0
    [in]real(dp), dimension(:, :) :: LAI0LAI at level-0
    [in]logical, dimension(:, :) :: mask0mask at L0
    [in]integer(i4), dimension(:) :: cell_id0Cell ids of hi res field
    [in]integer(i4), dimension(:) :: upp_row_L1Upper row of hi res block
    [in]integer(i4), dimension(:) :: low_row_L1Lower row of hi res block
    [in]integer(i4), dimension(:) :: lef_col_L1Left column of hi res block
    [in]integer(i4), dimension(:) :: rig_col_L1Right column of hi res block
    [in]integer(i4), dimension(:) :: nL0_in_L1Number of L0 cells within a L1 cel
    [in,out]real(dp), dimension(:, :) :: L1_petLAIcorFactorpet cor factor at level-1
    Authors
    M. Cuneyd Demirel and Simon Stisen from GEUS.dk
    Date
    May. 2017

Definition at line 80 of file mo_mpr_pet.f90.

References mo_upscaling_operators::upscale_harmonic_mean().

Referenced by mo_multi_param_reg::mpr().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ priestley_taylor_alpha()

subroutine, public mo_mpr_pet::priestley_taylor_alpha ( real(dp), dimension(:, :), intent(in)  lai0,
real(dp), dimension(:), intent(in)  param,
logical, dimension(:, :), intent(in)  mask0,
real(dp), intent(in)  nodata,
integer(i4), dimension(:), intent(in)  cell_id0,
integer(i4), dimension(:), intent(in)  nl0_in_l1,
integer(i4), dimension(:), intent(in)  upp_row_l1,
integer(i4), dimension(:), intent(in)  low_row_l1,
integer(i4), dimension(:), intent(in)  lef_col_l1,
integer(i4), dimension(:), intent(in)  rig_col_l1,
real(dp), dimension(:, :), intent(out)  priestley_taylor_alpha1 
)

Regionalization of priestley taylor alpha.

estimation of priestley taylor alpha Global parameters needed (see mhm_parameter.nml):

  • param(1) = PriestleyTaylorCoeff
  • param(2) = PriestleyTaylorLAIcorr
    Parameters
    [in]real(dp), dimension(:, :) :: LAI0LAI at level-0
    [in]real(dp), dimension(:) :: paraminput parameter
    [in]logical, dimension(:, :) :: mask0mask at level 0
    [in]real(dp) :: nodata- nodata value
    [in]integer(i4), dimension(:) :: cell_id0Cell ids of hi res field
    [in]integer(i4), dimension(:) :: nL0_in_L1number of l0 cells within a l1 cell
    [in]integer(i4), dimension(:) :: Upp_row_L1upper row of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Low_row_L1lower row of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Lef_col_L1left col of a l1 cell in l0 grid
    [in]integer(i4), dimension(:) :: Rig_col_L1right col of a l1 cell in l0 grid
    [out]real(dp), dimension(:, :) :: priestley_taylor_alpha1bulk surface resistance
    Authors
    Matthias Zink
    Date
    Apr 2013

Definition at line 311 of file mo_mpr_pet.f90.

References mo_upscaling_operators::upscale_arithmetic_mean().

Referenced by mo_multi_param_reg::mpr().

Here is the call graph for this function:
Here is the caller graph for this function: