mHM
The mesoscale Hydrological Model
|
Temporal disaggregation of daily input values. More...
Functions/Subroutines | |
elemental subroutine, public | temporal_disagg_meteo_weights (meteo_val_day, meteo_val_weights, meteo_val, weights_correction) |
Temporally distribute daily mean forcings onto time step. | |
elemental subroutine, public | temporal_disagg_flux_daynight (isday, ntimesteps_day, meteo_val_day, fday_meteo_val, fnight_meteo_val, meteo_val) |
Temporally distribute daily mean forcings onto time step. | |
elemental subroutine, public | temporal_disagg_state_daynight (isday, ntimesteps_day, meteo_val_day, fday_meteo_val, fnight_meteo_val, meteo_val, add_correction) |
Temporally distribute daily mean state forcings onto time step. | |
Temporal disaggregation of daily input values.
Calculate actual values for precipitation, PET and temperature from daily mean inputs
COPYING
and COPYING.LESSER
provided with this software. The complete GNU license text can also be found at http://www.gnu.org/licenses/. elemental subroutine, public mo_meteo_temporal_tools::temporal_disagg_flux_daynight | ( | logical, intent(in) | isday, |
real(dp), intent(in) | ntimesteps_day, | ||
real(dp), intent(in) | meteo_val_day, | ||
real(dp), intent(in) | fday_meteo_val, | ||
real(dp), intent(in) | fnight_meteo_val, | ||
real(dp), intent(out) | meteo_val | ||
) |
Temporally distribute daily mean forcings onto time step.
Calculates actual meteo forcings from daily mean inputs. They are distributed with predefined factors/summands for day and night.
[in] | isday | is day (False for night) |
[in] | ntimesteps_day | number of time steps per day |
[in] | meteo_val_day | Daily meteo_val |
[in] | fday_meteo_val | Daytime fraction of meteo_val |
[in] | fnight_meteo_val | Nighttime fraction of meteo_val |
[out] | meteo_val | Actual meteo_val |
Definition at line 66 of file mo_meteo_temporal_tools.f90.
Referenced by mo_meteo_handler::get_corrected_pet(), and mo_meteo_handler::get_prec().
elemental subroutine, public mo_meteo_temporal_tools::temporal_disagg_meteo_weights | ( | real(dp), intent(in) | meteo_val_day, |
real(dp), intent(in) | meteo_val_weights, | ||
real(dp), intent(out) | meteo_val, | ||
real(dp), intent(in), optional | weights_correction | ||
) |
Temporally distribute daily mean forcings onto time step.
Calculates actual meteo forcings from daily mean inputs. They are distributed with predefined weights onto the day.
[in] | meteo_val_day | Daily meteo_val |
[in] | meteo_val_weights | weights for meteo_val |
[out] | meteo_val | Actual meteo_val |
[in] | weights_correction | Additive correction before applying weights (e.g. for temperature conversion) |
Definition at line 32 of file mo_meteo_temporal_tools.f90.
Referenced by mo_meteo_handler::get_corrected_pet(), mo_meteo_handler::get_prec(), and mo_meteo_handler::get_temp().
elemental subroutine, public mo_meteo_temporal_tools::temporal_disagg_state_daynight | ( | logical, intent(in) | isday, |
real(dp), intent(in) | ntimesteps_day, | ||
real(dp), intent(in) | meteo_val_day, | ||
real(dp), intent(in) | fday_meteo_val, | ||
real(dp), intent(in) | fnight_meteo_val, | ||
real(dp), intent(out) | meteo_val, | ||
logical, intent(in), optional | add_correction | ||
) |
Temporally distribute daily mean state forcings onto time step.
Calculates meteo forcings from daily mean inputs of a state variable. They are distributed with predefined factors/summands for day and night.
[in] | isday | is day (False for night) |
[in] | ntimesteps_day | number of time steps per day |
[in] | meteo_val_day | Daily meteo_val |
[in] | fday_meteo_val | Daytime fraction of meteo_val |
[in] | fnight_meteo_val | Nighttime fraction of meteo_val |
[out] | meteo_val | Actual meteo_val |
[in] | add_correction | if True, correcting values will be added (e.g. for temperature), otherwise used as percentage |
Definition at line 108 of file mo_meteo_temporal_tools.f90.
Referenced by mo_meteo_handler::get_ssrd(), mo_meteo_handler::get_strd(), and mo_meteo_handler::get_temp().