mHM
The mesoscale Hydrological Model
|
Multiscale parameter regionalization (MPR). More...
Functions/Subroutines | |
subroutine, public | mpr (mask0, geounit0, soilid0, asp0, gridded_lai0, lcover0, slope_emp0, y0, id0, upper_bound1, lower_bound1, left_bound1, right_bound1, n_subcells1, fsealed1, alpha1, degdayinc1, degdaymax1, degdaynopre1, fasp1, harsamcoeff1, prietayalpha1, aeroresist1, surfresist1, froots1, kfastflow1, kslowflow1, kbaseflow1, kperco1, karstloss1, soilmoistfc1, soilmoistsat1, soilmoistexp1, jarvis_thresh_c1, tempthresh1, unsatthresh1, sealedthresh1, wiltingpoint1, maxinter1, petlaicorfactor, no_count1, bulkdens1, latticewater1, cosmicl31, parameterset) |
Regionalizing and Upscaling process parameters. | |
subroutine | baseflow_param (param, geounit0, k2_0) |
baseflow recession parameter | |
subroutine | snow_acc_melt_param (param, fforest1, fiperm1, fperm1, tempthresh1, degdaynopre1, degdayinc1, degdaymax1) |
Calculates the snow parameters. | |
subroutine | iper_thres_runoff (param, sealedthresh1) |
sets the impervious layer threshold parameter for runoff generation | |
subroutine | karstic_layer (param, geounit0, mask0, sms_fc0, ksvar_v0, id0, n_subcells1, upper_bound1, lower_bound1, left_bound1, right_bound1, karstloss1, l1_kp) |
calculates the Karstic percolation loss | |
subroutine, public | canopy_intercept_param (processmatrix, param, lai0, n_subcells1, upper_bound1, lower_bound1, left_bound1, right_bound1, id0, mask0, nodata, max_intercept1) |
estimate effective maximum interception capacity at L1 | |
subroutine | aerodynamical_resistance (lai0, lcover0, param, mask0, id0, n_subcells1, upper_bound1, lower_bound1, left_bound1, right_bound1, aerodyn_resistance1) |
Regionalization of aerodynamic resistance. | |
Multiscale parameter regionalization (MPR).
This module provides the routines for multiscale parameter regionalization (MPR).
COPYING
and COPYING.LESSER
provided with this software. The complete GNU license text can also be found at http://www.gnu.org/licenses/.
|
private |
Regionalization of aerodynamic resistance.
estimation of aerodynamical resistance Global parameters needed (see mhm_parameter.nml):
[in] | real(dp), dimension(:, :) :: LAI0 | LAI at level-0 |
[in] | integer(i4), dimension(:) :: LCover0 | land cover field |
[in] | real(dp), dimension(6) :: param | input parameter |
[in] | logical, dimension(:, :) :: mask0 | mask at level 0 |
[in] | integer(i4), dimension(:) :: Id0 | Cell ids of hi res field |
[in] | integer(i4), dimension(:) :: n_subcells1 | number of l0 cells within a l1 cell |
[in] | integer(i4), dimension(:) :: upper_bound1 | upper row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: lower_bound1 | lower row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: left_bound1 | left col of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: right_bound1 | right col of a l1 cell in l0 grid |
[out] | real(dp), dimension(:, :) :: aerodyn_resistance1 | aerodynmaical resistance |
Definition at line 1203 of file mo_multi_param_reg.f90.
References mo_common_constants::eps_dp, mo_mpr_constants::karman, mo_common_constants::nodata_dp, mo_upscaling_operators::upscale_arithmetic_mean(), and mo_mpr_constants::windmeasheight.
Referenced by mpr().
|
private |
baseflow recession parameter
This subroutine calculates the baseflow recession parameter based on the geological units at the Level 0 scale. For each level 0 cell, it assigns the value specified in the parameter array param for the geological unit in this cell. Global parameters needed (see mhm_parameter.nml):
[in] | real(dp), dimension(:) :: param | list of required parameters |
[in] | integer(i4), dimension(:) :: geoUnit0 | ids of geological units at L0 |
[out] | real(dp), dimension(:) :: k2_0 | - baseflow recession parameter at Level 0 |
Definition at line 689 of file mo_multi_param_reg.f90.
References mo_mpr_global_variables::geounitlist, and mo_common_constants::nodata_dp.
Referenced by mpr().
subroutine, public mo_multi_param_reg::canopy_intercept_param | ( | integer(i4), dimension(:, :), intent(in) | processmatrix, |
real(dp), dimension(:), intent(in) | param, | ||
real(dp), dimension(:, :), intent(in) | lai0, | ||
integer(i4), dimension(:), intent(in) | n_subcells1, | ||
integer(i4), dimension(:), intent(in) | upper_bound1, | ||
integer(i4), dimension(:), intent(in) | lower_bound1, | ||
integer(i4), dimension(:), intent(in) | left_bound1, | ||
integer(i4), dimension(:), intent(in) | right_bound1, | ||
integer(i4), dimension(:), intent(in) | id0, | ||
logical, dimension(:, :), intent(in) | mask0, | ||
real(dp), intent(in) | nodata, | ||
real(dp), dimension(:, :), intent(out) | max_intercept1 | ||
) |
estimate effective maximum interception capacity at L1
estimate effective maximum interception capacity at L1 for a given Leaf Area Index field. Global parameters needed (see mhm_parameter.nml): Process Case 1:
[in] | integer(i4), dimension(:, :) :: processMatrix | indicate processes |
[in] | real(dp), dimension(:) :: param | array of global parameters |
[in] | real(dp), dimension(:, :) :: LAI0 | LAI at level-0(nCells0, time) |
[in] | integer(i4), dimension(:) :: n_subcells1 | Number of L0 cells within a L1 cell |
[in] | integer(i4), dimension(:) :: upper_bound1 | Upper row of high resolution block |
[in] | integer(i4), dimension(:) :: lower_bound1 | Lower row of high resolution block |
[in] | integer(i4), dimension(:) :: left_bound1 | Left column of high resolution block |
[in] | integer(i4), dimension(:) :: right_bound1 | Right column of high resolution block |
[in] | integer(i4), dimension(:) :: Id0 | Cell ids at level 0 |
[in] | logical, dimension(:, :) :: mask0 | mask at level 0 field |
[in] | real(dp) :: nodata | - nodata value |
[out] | real(dp), dimension(:, :) :: max_intercept1 | max interception at level-1(nCells1, time) |
Definition at line 1076 of file mo_multi_param_reg.f90.
References mo_upscaling_operators::upscale_arithmetic_mean().
Referenced by mpr().
|
private |
sets the impervious layer threshold parameter for runoff generation
to be done by Kumar .... Global parameters needed (see mhm_parameter.nml):
[in] | real(dp), dimension(1) :: param | - given threshold parameter |
[out] | real(dp), dimension(:, :, :) :: sealedThresh1 |
Definition at line 883 of file mo_multi_param_reg.f90.
Referenced by mpr().
|
private |
calculates the Karstic percolation loss
This subroutine calls first the karstic_fraction upscaling routine for determine the karstic fraction area for every Level 1 cell. Then, the karstic percolation loss is estimated given two shape parameters by
\[ karstLoss1 = 1 + ( fKarArea * param(1)) *( (-1)**INT(param(2),i4) ) \]
where \( karstLoss1 \) is the karstic percolation loss and \( fKarArea \) is the fraction of karstic area at level 1 Global parameters needed (see mhm_parameter.nml):
[in] | real(dp), dimension(3) :: param | parameters |
[in] | integer(i4), dimension(:) :: geoUnit0 | id of the Karstic formation |
[in] | logical, dimension(:, :) :: mask0 | mask at level 0 |
[in] | real(dp), dimension(:) :: SMs_FC0 | [-] soil mositure deficit from field capacity w.r.t to saturation |
[in] | real(dp), dimension(:) :: KsVar_V0 | [-] relative variability of saturated |
[in] | integer(i4), dimension(:) :: Id0 | Cell ids of hi res field |
[in] | integer(i4), dimension(:) :: n_subcells1 | number of l0 cells within a l1 cell |
[in] | integer(i4), dimension(:) :: upper_bound1 | upper row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: lower_bound1 | lower row of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: left_bound1 | left col of a l1 cell in l0 grid |
[in] | integer(i4), dimension(:) :: right_bound1 | right col of a l1 cell in l0 grid |
[out] | real(dp), dimension(:) :: karstLoss1 | [-] Karstic percolation loss |
[out] | real(dp), dimension(:) :: L1_Kp | [d-1] percolation coefficient |
Definition at line 944 of file mo_multi_param_reg.f90.
References mo_mpr_global_variables::geounitkar, mo_mpr_global_variables::geounitlist, mo_upscaling_operators::l0_fractionalcover_in_lx(), mo_common_constants::nodata_dp, mo_common_constants::nodata_i4, and mo_upscaling_operators::upscale_arithmetic_mean().
Referenced by mpr().
subroutine, public mo_multi_param_reg::mpr | ( | logical, dimension(:, :), intent(in) | mask0, |
integer(i4), dimension(:), intent(in) | geounit0, | ||
integer(i4), dimension(:, :), intent(in) | soilid0, | ||
real(dp), dimension(:), intent(in) | asp0, | ||
real(dp), dimension(:, :), intent(in) | gridded_lai0, | ||
integer(i4), dimension(:, :), intent(in) | lcover0, | ||
real(dp), dimension(:), intent(in) | slope_emp0, | ||
real(dp), dimension(:), intent(in) | y0, | ||
integer(i4), dimension(:), intent(in) | id0, | ||
integer(i4), dimension(:), intent(in) | upper_bound1, | ||
integer(i4), dimension(:), intent(in) | lower_bound1, | ||
integer(i4), dimension(:), intent(in) | left_bound1, | ||
integer(i4), dimension(:), intent(in) | right_bound1, | ||
integer(i4), dimension(:), intent(in) | n_subcells1, | ||
real(dp), dimension(:, :, :), intent(inout) | fsealed1, | ||
real(dp), dimension(:, :, :), intent(inout) | alpha1, | ||
real(dp), dimension(:, :, :), intent(inout) | degdayinc1, | ||
real(dp), dimension(:, :, :), intent(inout) | degdaymax1, | ||
real(dp), dimension(:, :, :), intent(inout) | degdaynopre1, | ||
real(dp), dimension(:, :, :), intent(inout) | fasp1, | ||
real(dp), dimension(:, :, :), intent(inout) | harsamcoeff1, | ||
real(dp), dimension(:, :, :), intent(inout) | prietayalpha1, | ||
real(dp), dimension(:, :, :), intent(inout) | aeroresist1, | ||
real(dp), dimension(:, :, :), intent(inout) | surfresist1, | ||
real(dp), dimension(:, :, :), intent(inout) | froots1, | ||
real(dp), dimension(:, :, :), intent(inout) | kfastflow1, | ||
real(dp), dimension(:, :, :), intent(inout) | kslowflow1, | ||
real(dp), dimension(:, :, :), intent(inout) | kbaseflow1, | ||
real(dp), dimension(:, :, :), intent(inout) | kperco1, | ||
real(dp), dimension(:, :, :), intent(inout) | karstloss1, | ||
real(dp), dimension(:, :, :), intent(inout) | soilmoistfc1, | ||
real(dp), dimension(:, :, :), intent(inout) | soilmoistsat1, | ||
real(dp), dimension(:, :, :), intent(inout) | soilmoistexp1, | ||
real(dp), dimension(:, :, :), intent(inout) | jarvis_thresh_c1, | ||
real(dp), dimension(:, :, :), intent(inout) | tempthresh1, | ||
real(dp), dimension(:, :, :), intent(inout) | unsatthresh1, | ||
real(dp), dimension(:, :, :), intent(inout) | sealedthresh1, | ||
real(dp), dimension(:, :, :), intent(inout) | wiltingpoint1, | ||
real(dp), dimension(:, :, :), intent(inout) | maxinter1, | ||
real(dp), dimension(:, :, :), intent(inout) | petlaicorfactor, | ||
real(dp), dimension(:, :, :), intent(inout) | no_count1, | ||
real(dp), dimension(:, :, :), intent(inout) | bulkdens1, | ||
real(dp), dimension(:, :, :), intent(inout) | latticewater1, | ||
real(dp), dimension(:, :, :), intent(inout) | cosmicl31, | ||
real(dp), dimension(:), intent(in), optional, target | parameterset | ||
) |
Regionalizing and Upscaling process parameters.
calculating process parameters at L0 scale (Regionalization), like:
and upscale these parameters to retrieve effective parameters at scale L1. Further parameter regionalizations are done for:
[in] | mask0 | mask at level 0 field |
[in] | geounit0 | L0 geological units |
[in] | soilid0 | soil Ids at level 0 |
[in] | asp0 | [degree] Aspect at Level 0 |
[in] | gridded_lai0 | LAI grid at level 0, with dim2 = time |
[in] | lcover0 | land cover at level 0 |
[in] | slope_emp0 | Empirical quantiles of slope |
[in] | id0 | Cell ids at level 0 |
[in] | upper_bound1 | Upper row of hi res block |
[in] | lower_bound1 | Lower row of hi res block |
[in] | left_bound1 | Left column of hi res block |
[in] | right_bound1 | Right column of hi res block |
[in] | n_subcells1 | Number of L0 cells within a L1 cell |
[in] | y0 | y0 at level 0 |
[in,out] | fsealed1 | [1] fraction of sealed area |
[in,out] | soilmoistexp1 | Parameter that determines the rel. contribution to SM, upscal. Bulk den. |
[in,out] | jarvis_thresh_c1 | [1] jarvis critical value for norm SWC |
[in,out] | soilmoistsat1 | [10^-3 m] depth of saturated SM |
[in,out] | soilmoistfc1 | [10^-3 m] field capacity |
[in,out] | wiltingpoint1 | [10^-3 m] permanent wilting point |
[in,out] | froots1 | fraction of roots in soil horizon |
[in,out] | tempthresh1 | [degreeC] threshold temperature for snow rain |
[in,out] | degdaynopre1 | [mm-1 degreeC-1] Degree-day factor with no precipitation |
[in,out] | degdaymax1 | [mm-1 degreeC-1] Maximum Degree-day factor |
[in,out] | degdayinc1 | [d-1 degreeC-1] Increase of the Degree-day factor per mm of increase in precipitation |
[in,out] | fasp1 | [1] PET correction for Aspect at level 1 |
[in,out] | harsamcoeff1 | [1] PET Hargreaves Samani coeff. at level 1 |
[in,out] | prietayalpha1 | [1] PET Priestley Taylor coeff. at level 1 |
[in,out] | aeroresist1 | [s m-1] PET aerodynamical resitance at level 1 |
[in,out] | surfresist1 | [s m-1] PET bulk surface resitance at level 1 |
[in,out] | sealedthresh1 | threshold parameter |
[in,out] | unsatthresh1 | [10^-3 m] Threshhold water depth in upper reservoir (for Runoff contribution) |
[in,out] | kfastflow1 | [10^-3 m] Recession coefficient of the upper reservoir, upper outlet |
[in,out] | kslowflow1 | [10^-3 m] Recession coefficient of the upper reservoir, lower outlet |
[in,out] | kbaseflow1 | Level 1 baseflow recession |
[in,out] | alpha1 | [1] Exponent for the upper reservoir |
[in,out] | kperco1 | [d-1] percolation coefficient |
[in,out] | karstloss1 | karstic percolation loss |
[in,out] | maxinter1 | max interception |
[in,out] | petlaicorfactor | pet cor factor at level-1 |
[in,out] | no_count1 | [-] inital neutron count |
[in,out] | bulkdens1 | [gcm-3] bulk density |
[in,out] | latticewater1 | [mm/mm] lattice water content |
[in,out] | cosmicl31 | [-] cosmic L3 parameter |
[in] | parameterset | array of global parameters |
Definition at line 67 of file mo_multi_param_reg.f90.
References aerodynamical_resistance(), baseflow_param(), mo_mpr_pet::bulksurface_resistance(), canopy_intercept_param(), mo_mpr_global_variables::fracsealed_cityarea, mo_common_variables::global_parameters, mo_mpr_global_variables::horizondepth_mhm, mo_mpr_global_variables::iflag_soildb, iper_thres_runoff(), karstic_layer(), mo_upscaling_operators::l0_fractionalcover_in_lx(), mo_mpr_neutrons::mpr_neutrons(), mo_mpr_runoff::mpr_runoff(), mo_mpr_soilmoist::mpr_sm(), mo_mpr_smhorizons::mpr_smhorizons(), mo_common_constants::nodata_dp, mo_common_constants::nodata_i4, mo_mpr_global_variables::nsoilhorizons_mhm, mo_mpr_pet::pet_correctbyasp(), mo_mpr_pet::pet_correctbylai(), mo_mpr_pet::priestley_taylor_alpha(), mo_common_variables::processmatrix, snow_acc_melt_param(), mo_mpr_global_variables::soildb, and mo_upscaling_operators::upscale_arithmetic_mean().
Referenced by mo_mpr_eval::mpr_eval().
|
private |
Calculates the snow parameters.
This subroutine calculates the snow parameters threshold temperature (TT), degree-day factor without precipitation (DD) and maximum degree-day factor (DDmax) as well as increase of degree-day factor per mm of increase in precipitation (IDDP). Global parameters needed (see mhm_parameter.nml):
[in] | real(dp), dimension(8) :: param | eight global parameters |
[in] | real(dp), dimension(:) :: fForest1 | [1] fraction of forest cover |
[in] | real(dp), dimension(:) :: fIperm1 | [1] fraction of sealed area |
[in] | real(dp), dimension(:) :: fPerm1 | [1] fraction of permeable area |
[out] | real(dp), dimension(:) :: tempThresh1 | [degreeC] threshold temperature for snow rain |
[out] | real(dp), dimension(:) :: degDayNoPre1 | [mm-1 degreeC-1] Degree-day factor with |
[out] | real(dp), dimension(:) :: degDayInc1 | [d-1 degreeC-1] Increase of the Degree-day |
[out] | real(dp), dimension(:) :: degDayMax1 | [mm-1 degreeC-1] Maximum Degree-day factor |
Definition at line 800 of file mo_multi_param_reg.f90.
Referenced by mpr().