mHM
The mesoscale Hydrological Model
|
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. | |
MPR routine for PET.
This module sets up pet correction factor at level-1 based on LAI
COPYING
and COPYING.LESSER
provided with this software. The complete GNU license text can also be found at http://www.gnu.org/licenses/. 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):
[in] | real(dp), dimension(:, :) :: LAI0 | LAI at level-0 |
[in] | real(dp) :: param | - global parameter |
[in] | logical, dimension(:, :) :: mask0 | mask at level 0 |
[in] | real(dp) :: nodata | - nodata value |
[in] | integer(i4), dimension(:) :: cell_id0 | Cell ids of hi res field |
[in] | integer(i4), dimension(:) :: nL0_in_L1 | number of l0 cells within a l1 cell |
[in] | integer(i4), dimension(:) :: Upp_row_L1 | upper row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Low_row_L1 | lower row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Lef_col_L1 | left col of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Rig_col_L1 | right col of a l1 cell in l0 grid |
[out] | real(dp), dimension(:, :) :: bulksurface_resistance1 | bulk surface resistance |
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().
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):
[in] | integer(i4), dimension(:) :: id0 | Level 0 cell id |
[in] | real(dp), dimension(:) :: latitude_l0 | latitude on l0 |
[in] | real(dp), dimension(:) :: Asp0 | [degree] Aspect at Level 0 |
[in] | real(dp), dimension(3) :: param | process parameters |
[in] | real(dp) :: nodata | - no data value |
[out] | real(dp), dimension(:) :: fAsp0 | PET correction for Aspect |
Definition at line 214 of file mo_mpr_pet.f90.
Referenced by mo_multi_param_reg::mpr().
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:
[in] | real(dp), dimension(5) :: param | parameters |
[in] | real(dp) :: nodata | - nodata value |
[in] | integer(i4), dimension(:) :: LCOVER0 | Land cover at level 0 |
[in] | real(dp), dimension(:, :) :: LAI0 | LAI at level-0 |
[in] | logical, dimension(:, :) :: mask0 | mask at L0 |
[in] | integer(i4), dimension(:) :: cell_id0 | Cell ids of hi res field |
[in] | integer(i4), dimension(:) :: upp_row_L1 | Upper row of hi res block |
[in] | integer(i4), dimension(:) :: low_row_L1 | Lower row of hi res block |
[in] | integer(i4), dimension(:) :: lef_col_L1 | Left column of hi res block |
[in] | integer(i4), dimension(:) :: rig_col_L1 | Right column of hi res block |
[in] | integer(i4), dimension(:) :: nL0_in_L1 | Number of L0 cells within a L1 cel |
[in,out] | real(dp), dimension(:, :) :: L1_petLAIcorFactor | pet cor factor at level-1 |
Definition at line 80 of file mo_mpr_pet.f90.
References mo_upscaling_operators::upscale_harmonic_mean().
Referenced by mo_multi_param_reg::mpr().
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):
[in] | real(dp), dimension(:, :) :: LAI0 | LAI at level-0 |
[in] | real(dp), dimension(:) :: param | input parameter |
[in] | logical, dimension(:, :) :: mask0 | mask at level 0 |
[in] | real(dp) :: nodata | - nodata value |
[in] | integer(i4), dimension(:) :: cell_id0 | Cell ids of hi res field |
[in] | integer(i4), dimension(:) :: nL0_in_L1 | number of l0 cells within a l1 cell |
[in] | integer(i4), dimension(:) :: Upp_row_L1 | upper row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Low_row_L1 | lower row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Lef_col_L1 | left col of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: Rig_col_L1 | right col of a l1 cell in l0 grid |
[out] | real(dp), dimension(:, :) :: priestley_taylor_alpha1 | bulk surface resistance |
Definition at line 311 of file mo_mpr_pet.f90.
References mo_upscaling_operators::upscale_arithmetic_mean().
Referenced by mo_multi_param_reg::mpr().