14 use mo_kind,
only : i4, dp
23 meteo_time_ref_endpoint, &
29 meteo_expect_netrad, &
30 meteo_expect_absvappress, &
31 meteo_expect_windspeed, &
38 integer(i4),
intent(in) :: couple_case
40 integer(i4),
intent(in) :: meteo_timestep
42 logical,
intent(in) :: meteo_time_ref_endpoint
44 logical,
intent(in) :: meteo_expect_pre
46 logical,
intent(in) :: meteo_expect_temp
48 logical,
intent(in) :: meteo_expect_pet
50 logical,
intent(in) :: meteo_expect_tmin
52 logical,
intent(in) :: meteo_expect_tmax
54 logical,
intent(in) :: meteo_expect_netrad
56 logical,
intent(in) :: meteo_expect_absvappress
58 logical,
intent(in) :: meteo_expect_windspeed
60 logical,
intent(in) :: meteo_expect_ssrd
62 logical,
intent(in) :: meteo_expect_strd
64 logical,
intent(in) :: meteo_expect_tann
69 meteo_time_ref_endpoint, &
75 meteo_expect_netrad, &
76 meteo_expect_absvappress, &
77 meteo_expect_windspeed, &
85 subroutine init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
88 character(*),
intent(in) :: namelist_mhm
90 character(*),
intent(in) :: namelist_mhm_param
92 character(*),
intent(in) :: namelist_mhm_output
94 character(*),
intent(in) :: namelist_mrm_output
96 character(*),
intent(in) :: cwd
98 call mhm_interface_init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
122 real(dp),
intent(in) ::
parameter(n)
146 character(64),
intent(out) :: ver_string
147 ver_string = mhm_version
161 integer(i4),
intent(in) :: level
172 use mo_kind,
only : i4, dp
186 integer(i4),
intent(out) :: n
194 integer(i4),
intent(in) :: domain
203 logical,
intent(out) :: output
239 integer(i4),
intent(out) :: year
240 integer(i4),
intent(out) :: month
241 integer(i4),
intent(out) :: day
242 integer(i4),
intent(out) :: hour
243 year =
run_cfg%domainDateTime%year
244 month =
run_cfg%domainDateTime%month
245 day =
run_cfg%domainDateTime%day
246 hour =
run_cfg%domainDateTime%hour
256 use mo_kind,
only : i4, dp
267 integer(i4),
intent(out) :: n
277 integer(i4),
intent(out) :: length
289 real(dp),
intent(out) :: config(n, 5)
302 character(256),
intent(out) :: para_name
304 integer(i4),
intent(in) :: n
316 integer(i4),
intent(in) :: gauge_id
317 integer(i4),
intent(out) :: length
319 integer(i4) :: iDomain
322 idomain =
gauge%domainId(gauge_id)
333 real(dp),
intent(out) :: output(m, 2)
334 integer(i4),
intent(in) :: gauge_id
337 real(dp),
dimension(:),
allocatable :: runoff_agg, runoff_obs
338 logical,
dimension(:),
allocatable :: runoff_obs_mask
341 output(:, 1) = runoff_agg
342 output(:, 2) = runoff_obs
351 integer(i4),
intent(out) :: shp(2)
364 real(dp),
intent(out) :: output(m, n)
375 integer(i4),
intent(out) :: n
376 integer(i4),
intent(in) :: domain
378 integer(i4) :: iDomain, i
380 if ( i == 0 ) i =
run_cfg%selected_domain
381 idomain =
run_cfg%get_domain_index(i)
382 n =
level0(idomain)%nCells
392 integer(i4),
intent(out) :: shp(2)
393 integer(i4),
intent(in) :: domain
395 integer(i4) :: iDomain, i
397 if ( i == 0 ) i =
run_cfg%selected_domain
398 idomain =
run_cfg%get_domain_index(i)
399 shp = shape(
level0(idomain)%mask)
411 logical,
intent(out) :: mask(m, n)
412 integer(i4),
intent(in) :: domain
414 integer(i4) :: iDomain, i
416 if ( i == 0 ) i =
run_cfg%selected_domain
417 idomain =
run_cfg%get_domain_index(i)
418 mask =
level0(idomain)%mask
424 subroutine l0_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
429 integer(i4),
intent(out) :: ncols
430 integer(i4),
intent(out) :: nrows
431 integer(i4),
intent(out) :: ncells
432 real(dp),
intent(out) :: xll
433 real(dp),
intent(out) :: yll
434 real(dp),
intent(out) :: cell_size
435 real(dp),
intent(out) :: no_data
436 integer(i4),
intent(in) :: domain
438 integer(i4) :: iDomain, i
440 if ( i == 0 ) i =
run_cfg%selected_domain
441 idomain =
run_cfg%get_domain_index(i)
442 ncols =
level0(idomain)%ncols
443 nrows =
level0(idomain)%nrows
444 ncells =
level0(idomain)%nCells
445 xll =
level0(idomain)%xllcorner
446 yll =
level0(idomain)%yllcorner
447 cell_size =
level0(idomain)%cellsize
457 use mo_message,
only : error_message
463 real(dp),
intent(out) :: output(n)
464 character(*),
intent(in) :: name
465 integer(i4),
intent(in) :: idx
467 integer(i4) :: iDomain, s0, e0
473 case(
"L0_GRIDDED_LAI")
476 call error_message(
"mhm.get.L0_variable: unkown variable name: ", name)
487 integer(i4),
intent(out) :: n
488 integer(i4),
intent(in) :: domain
490 integer(i4) :: iDomain, i
492 if ( i == 0 ) i =
run_cfg%selected_domain
493 idomain =
run_cfg%get_domain_index(i)
494 n =
level1(idomain)%nCells
504 integer(i4),
intent(out) :: shp(2)
505 integer(i4),
intent(in) :: domain
507 integer(i4) :: iDomain, i
509 if ( i == 0 ) i =
run_cfg%selected_domain
510 idomain =
run_cfg%get_domain_index(i)
511 shp = shape(
level1(idomain)%mask)
523 logical,
intent(out) :: mask(m, n)
524 integer(i4),
intent(in) :: domain
526 integer(i4) :: iDomain, i
528 if ( i == 0 ) i =
run_cfg%selected_domain
529 idomain =
run_cfg%get_domain_index(i)
530 mask =
level1(idomain)%mask
537 subroutine l1_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
542 integer(i4),
intent(out) :: ncols
543 integer(i4),
intent(out) :: nrows
544 integer(i4),
intent(out) :: ncells
545 real(dp),
intent(out) :: xll
546 real(dp),
intent(out) :: yll
547 real(dp),
intent(out) :: cell_size
548 real(dp),
intent(out) :: no_data
549 integer(i4),
intent(in) :: domain
551 integer(i4) :: iDomain, i
553 if ( i == 0 ) i =
run_cfg%selected_domain
554 idomain =
run_cfg%get_domain_index(i)
555 ncols =
level1(idomain)%ncols
556 nrows =
level1(idomain)%nrows
557 ncells =
level1(idomain)%nCells
558 xll =
level1(idomain)%xllcorner
559 yll =
level1(idomain)%yllcorner
560 cell_size =
level1(idomain)%cellsize
570 use mo_message,
only : error_message
603 real(dp),
intent(out) :: output(n)
604 character(*),
intent(in) :: name
605 integer(i4),
intent(in) :: idx
610 case(
"L1_FNOTSEALED")
634 case(
"L1_TOTAL_RUNOFF")
636 case(
"L1_RUNOFFSEAL")
638 case(
"L1_FASTRUNOFF")
640 case(
"L1_SLOWRUNOFF")
650 case(
"L1_SOILMOIST_VOL")
653 case(
"L1_SOILMOIST_VOL_ALL")
656 case(
"L1_SOILMOISTSAT")
668 case(
"L1_THROUGHFALL")
671 call error_message(
"mhm.get.L1_variable: unkown variable name: ", name)
682 integer(i4),
intent(out) :: n
683 integer(i4),
intent(in) :: domain
685 integer(i4) :: iDomain, i
687 if ( i == 0 ) i =
run_cfg%selected_domain
688 idomain =
run_cfg%get_domain_index(i)
699 integer(i4),
intent(out) :: shp(2)
700 integer(i4),
intent(in) :: domain
702 integer(i4) :: iDomain, i
704 if ( i == 0 ) i =
run_cfg%selected_domain
705 idomain =
run_cfg%get_domain_index(i)
706 shp = shape(
level11(idomain)%mask)
718 logical,
intent(out) :: mask(m, n)
719 integer(i4),
intent(in) :: domain
721 integer(i4) :: iDomain, i
723 if ( i == 0 ) i =
run_cfg%selected_domain
724 idomain =
run_cfg%get_domain_index(i)
731 subroutine l11_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
736 integer(i4),
intent(out) :: ncols
737 integer(i4),
intent(out) :: nrows
738 integer(i4),
intent(out) :: ncells
739 real(dp),
intent(out) :: xll
740 real(dp),
intent(out) :: yll
741 real(dp),
intent(out) :: cell_size
742 real(dp),
intent(out) :: no_data
743 integer(i4),
intent(in) :: domain
745 integer(i4) :: iDomain, i
747 if ( i == 0 ) i =
run_cfg%selected_domain
748 idomain =
run_cfg%get_domain_index(i)
751 ncells =
level11(idomain)%nCells
752 xll =
level11(idomain)%xllcorner
753 yll =
level11(idomain)%yllcorner
754 cell_size =
level11(idomain)%cellsize
764 use mo_message,
only : error_message
773 real(dp),
intent(out) :: output(n)
774 character(*),
intent(in) :: name
775 integer(i4),
intent(in) :: idx
787 call error_message(
"mhm.get.L11_variable: unkown variable name: ", name)
798 integer(i4),
intent(out) :: n
799 integer(i4),
intent(in) :: domain
801 integer(i4) :: iDomain, i
803 if ( i == 0 ) i =
run_cfg%selected_domain
804 idomain =
run_cfg%get_domain_index(i)
815 integer(i4),
intent(out) :: shp(2)
816 integer(i4),
intent(in) :: domain
818 integer(i4) :: iDomain, i
820 if ( i == 0 ) i =
run_cfg%selected_domain
821 idomain =
run_cfg%get_domain_index(i)
834 logical,
intent(out) :: mask(m, n)
835 integer(i4),
intent(in) :: domain
837 integer(i4) :: iDomain, i
839 if ( i == 0 ) i =
run_cfg%selected_domain
840 idomain =
run_cfg%get_domain_index(i)
848 subroutine l2_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
853 integer(i4),
intent(out) :: ncols
854 integer(i4),
intent(out) :: nrows
855 integer(i4),
intent(out) :: ncells
856 real(dp),
intent(out) :: xll
857 real(dp),
intent(out) :: yll
858 real(dp),
intent(out) :: cell_size
859 real(dp),
intent(out) :: no_data
860 integer(i4),
intent(in) :: domain
862 integer(i4) :: iDomain, i
864 if ( i == 0 ) i =
run_cfg%selected_domain
865 idomain =
run_cfg%get_domain_index(i)
884 use mo_kind,
only : i4, dp
893 use mo_message,
only: error_message
899 real(dp),
intent(in) :: input(n)
900 character(*),
intent(in) :: name
901 integer(i4),
intent(in) :: idx
903 integer(i4) :: iDomain, s0, e0
909 case(
"L0_GRIDDED_LAI")
912 call error_message(
"mhm.set.L0_variable: unkown variable name: ", name)
920 subroutine meteo(input, name, year, month, day, hour, n)
921 use mo_message,
only: error_message
927 real(dp),
intent(in) :: input(n)
928 character(*),
intent(in) :: name
929 integer(i4),
intent(in) :: year
930 integer(i4),
intent(in) :: month
931 integer(i4),
intent(in) :: day
932 integer(i4),
intent(in) :: hour
936 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, pre=input)
938 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, temp=input)
940 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, pet=input)
942 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tmin=input)
944 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tmax=input)
946 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, netrad=input)
948 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, absvappress=input)
950 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, windspeed=input)
952 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, ssrd=input)
954 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, strd=input)
956 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tann=input)
958 call error_message(
"mhm.set.meteo: unkown variable name: ", name)
Python wrapper module to get internal variables of a mHM model run.
subroutine runoff_shape(shp)
Get the shape of mHM model runoff output.
subroutine runoff_eval_length(gauge_id, length)
Get the shape of mHM model runoff output for evaluation.
subroutine l11_domain_shape(shp, domain)
Get the shape of Level-11 of the mHM model.
subroutine runoff(output, m, n)
Get the mHM model runoff output.
subroutine parameter_config(config, n)
Get the parameter settings of mHM.
subroutine l2_domain_size(n, domain)
Get number of unmasked celles on Level-2 of the mHM model.
subroutine l11_domain_mask(mask, n, m, domain)
Get the mask of Level-11 of the mHM model.
subroutine l11_domain_size(n, domain)
Get number of unmasked celles on Level-11 of the mHM model.
subroutine l2_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-2 of the mHM model.
subroutine l11_variable(output, name, idx, n)
Get a variable on Level-11 of the mHM model.
subroutine l1_domain_shape(shp, domain)
Get the shape of Level-1 of the mHM model.
subroutine parameter_length(length)
Get the number of parameters in mHM.
subroutine l1_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-1 of the mHM model.
subroutine parameter_name(para_name, n)
Get the parameter names of mHM.
subroutine l0_domain_size(n, domain)
Get number of unmasked celles on Level-0 of the mHM model.
subroutine l0_domain_shape(shp, domain)
Get the shape of Level-0 of the mHM model.
subroutine l2_domain_mask(mask, n, m, domain)
Get the mask of Level-2 of the mHM model.
subroutine l1_variable(output, name, idx, n)
Get a variable on Level-1 of the mHM model.
subroutine l0_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-0 of the mHM model.
subroutine number_of_horizons(n)
Get the number of soil horizons in mHM.
subroutine l11_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-11 of the mHM model.
subroutine l1_domain_mask(mask, n, m, domain)
Get the mask of Level-1 of the mHM model.
subroutine l2_domain_shape(shp, domain)
Get the shape of Level-2 of the mHM model.
subroutine runoff_eval(gauge_id, output, m)
Get the mHM model runoff output.
subroutine l0_variable(output, name, idx, n)
Get a variable on Level-0 of the mHM model.
subroutine l1_domain_size(n, domain)
Get number of unmasked celles on Level-1 of the mHM model.
subroutine l0_domain_mask(mask, n, m, domain)
Get the mask of Level-0 of the mHM model.
Provides constants commonly used by mHM, mRM and MPR.
real(dp), parameter, public nodata_dp
Provides structures needed by mHM, mRM and/or mpr.
integer(i4), public ntstepday
integer(i4), dimension(:), allocatable, public warmingdays
type(period), dimension(:), allocatable, public evalper
Provides structures needed by mhm_eval to store current run config.
type(run_cfg_t), public run_cfg
This is a container to hold all information while running mHM.
Provides structures needed by mHM, mRM and/or mpr.
real(dp), dimension(:, :), allocatable, target, public global_parameters
character(256), dimension(:), allocatable, public global_parameters_name
type(domain_meta), public domainmeta
type(grid), dimension(:), allocatable, target, public level1
type(grid), dimension(:), allocatable, target, public level0
Provides file names and units for mHM.
character(len=*), parameter version
Current mHM model version (will be set to )
Main global variables for mHM.
type(meteo_handler_type), public meteo_handler
the meteo handler class
real(dp), dimension(:), allocatable, public l1_slowrunoff
[mm TS-1] Slow runoff component
real(dp), dimension(:, :), allocatable, public l1_aetsoil
[mm TS-1] Actual ET from soil layers
real(dp), dimension(:), allocatable, public l1_rain
[mm TS-1] Rain precipitation depth
real(dp), dimension(:), allocatable, public l1_aetsealed
[mm TS-1] Real evap.
real(dp), dimension(:), allocatable, public l1_satstw
[mm] groundwater storage
logical, dimension(noutflxstate) outputflxstate
Define model outputs see "mhm_outputs.nml" dim1 = number of output variables to be written.
real(dp), dimension(:), allocatable, public l1_aetcanopy
[mm TS-1] Real evaporation intensity from canopy
real(dp), dimension(:), allocatable, public l1_snow
[mm TS-1] Snow precipitation depth
real(dp), dimension(:), allocatable, public l1_preeffect
[mm TS-1] Effective precip.
real(dp), dimension(:), allocatable, public l1_inter
[mm] Interception
real(dp), dimension(:), allocatable, public l1_sealstw
[mm] Retention storage of impervious areas
real(dp), dimension(:), allocatable, public l1_percol
[mm TS-1] Percolation.
type(couple_cfg_type), public couple_cfg
coupling configuration class
real(dp), dimension(:, :), allocatable, public l1_soilmoist
[mm] Soil moisture of each horizon
real(dp), dimension(:), allocatable, public l1_pet_calc
[mm TS-1] estimated/corrected potential evapotranspiration
real(dp), dimension(:), allocatable, public l1_unsatstw
[mm] upper soil storage
real(dp), dimension(:), allocatable, public l1_melt
[mm TS-1] Melting snow depth.
real(dp), dimension(:), allocatable, public l1_fastrunoff
[mm TS-1] Fast runoff component
real(dp), dimension(:), allocatable, public l1_snowpack
[mm] Snowpack
real(dp), dimension(:), allocatable, public l1_neutrons
[mm] Ground Albedo Neutrons
real(dp), dimension(:), allocatable, public l1_prec_calc
[mm TS-1] precipitation for current time step
real(dp), dimension(:), allocatable, public l1_runoffseal
[mm TS-1] Direct runoff from impervious areas
real(dp), dimension(:, :), allocatable, public l1_infilsoil
[mm TS-1] Infiltration intensity each soil horizon
real(dp), dimension(:), allocatable, public l1_baseflow
[mm TS-1] Baseflow
real(dp), dimension(:), allocatable, public l1_throughfall
[mm TS-1] Throughfall.
real(dp), dimension(:), allocatable, public l1_temp_calc
[degC] temperature for current time step
real(dp), dimension(:), allocatable, public l1_total_runoff
[m3 TS-1] Generated runoff
Module to parse command line arguments of mHM.
subroutine, public set_verbosity_level(level)
Set the verbosity level of mHM.
Runs mhm with a specific parameter set and returns required variables, e.g.
subroutine, public mhm_eval(parameterset, opti_domain_indices, runoff, smoptisim, neutronsoptisim, etoptisim, twsoptisim, bfi)
Runs mhm with a specific parameter set and returns required variables, e.g.
Module providing interfaces for running preconfigured mHM.
subroutine mhm_interface_run_get_ndomains(ndomains)
get number of domains for looping
subroutine mhm_interface_run_prepare(parameterset, opti_domain_indices, runoff_present, bfi_present)
prepare single run of mHM
subroutine mhm_interface_run_write_output()
write output after current time-step
subroutine mhm_interface_run_prepare_domain(domain, etoptisim, twsoptisim, neutronsoptisim, smoptisim)
prepare single domain to run mHM on
subroutine mhm_interface_run_finished(time_loop_finished)
check if current time loop is finished
subroutine mhm_interface_run_do_time_step()
do one time-step on current domain
subroutine mhm_interface_run_finalize_domain()
finalize current domain after running
subroutine mhm_interface_run_finalize(runoff, bfi)
finalize run
Module providing interfaces for mHM.
subroutine, public mhm_interface_init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
initialize mHM from given namelist paths.
subroutine, public mhm_interface_run()
Run mHM with current settings.
subroutine, public mhm_interface_get_parameter_number(n)
Get number of current global parameter value of mHM.
subroutine, public mhm_interface_run_optimization()
Run mHM optimization with current settings.
subroutine, public mhm_interface_finalize()
Write mHM restart.
Global variables for mpr only.
integer(i4), public nsoilhorizons_mhm
real(dp), dimension(:, :), allocatable, public l0_gridded_lai
real(dp), dimension(:, :, :), allocatable, public l1_fsealed
real(dp), dimension(:, :, :), allocatable, public l1_soilmoistsat
Global variables for mRM only.
real(dp), dimension(:, :), allocatable, public l11_qtin
real(dp), dimension(:), allocatable, public l11_qout
real(dp), dimension(:, :), allocatable, public mrm_runoff
type(gaugingstation), public gauge
logical, dimension(noutflxstate) outputflxstate_mrm
Define model outputs see "mhm_outputs.nml".
real(dp), dimension(:), allocatable, public l11_qmod
type(grid), dimension(:), allocatable, target, public level11
real(dp), dimension(:, :), allocatable, public l11_qtr
integer(i4), public nmeasperday
Objective Functions for Optimization of mHM/mRM against runoff.
subroutine, public extract_runoff(gaugeid, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
extracts runoff data from global variables
Python wrapper module to control a mHM model.
subroutine version(ver_string)
Get the mHM version.
subroutine set_verbosity(level)
Set verbosity level of mHM.
subroutine disable_output()
disable all mHM/mRM outputs during runtime.
subroutine config_coupling(couple_case, meteo_timestep, meteo_time_ref_endpoint, meteo_expect_pre, meteo_expect_temp, meteo_expect_pet, meteo_expect_tmin, meteo_expect_tmax, meteo_expect_netrad, meteo_expect_absvappress, meteo_expect_windspeed, meteo_expect_ssrd, meteo_expect_strd, meteo_expect_tann)
Configure coupling mode of mHM.
subroutine run_or_optimize()
Execute a mHM model or an optimization depending on the configuration.
subroutine run_with_parameter(parameter, n)
Execute a mHM model with given parameters.
subroutine init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
Initialize a mHM model.
subroutine finalize()
Finalize a mHM model.
Python wrapper module to control a mHM model run per time step.
subroutine write_output()
Write output for the current domain of the current mHM model run.
subroutine prepare()
Prepare a mHM model run.
subroutine prepare_domain(domain)
Prepare a certain domain of the current mHM model run.
subroutine finished(output)
Check if the current mHM model time loop is finished.
subroutine current_time(year, month, day, hour)
Get the current time the current domain of the current mHM model run.
subroutine finalize_domain()
Finalize the current domain of the current mHM model run.
subroutine do_time_step()
Do one time-step on the current domain of the current mHM model run.
subroutine get_ndomains(n)
Get the number of domains of the current mHM model run.
subroutine finalize()
Finalize the current mHM model run.
Python wrapper module to set internal variables of a mHM model run.
subroutine meteo(input, name, year, month, day, hour, n)
Set a meteo variable on Level-1 of the mHM model.
subroutine l0_variable(input, name, idx, n)
Set a variable on Level-0 of the mHM model.