5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_mhm_eval.f90
Go to the documentation of this file.
1!> \file mo_mhm_eval.f90
2!> \brief \copybrief mo_mhm_eval
3!> \details \copydetails mo_mhm_eval
4
5!> \brief Runs mhm with a specific parameter set and returns required variables, e.g. runoff.
6!> \details Runs mhm with a specific parameter set and returns required variables, e.g. runoff.
7!> \authors Juliane Mai, Rohini Kumar
8!> \date Feb 2013
9!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
10!! mHM is released under the LGPLv3+ license \license_note
11!> \ingroup f_mhm
13
14 use mo_kind, only : i4, dp
15 use mo_optimization_types, only : optidata_sim
16 use mo_mhm_interface_run, only : &
26
27 IMPLICIT NONE
28
29 PRIVATE
30
31 PUBLIC :: mhm_eval
32
33 ! ------------------------------------------------------------------
34
35CONTAINS
36
37 ! ------------------------------------------------------------------
38
39 ! NAME
40 ! mhm_eval
41
42 ! PURPOSE
43 !> \brief Runs mhm with a specific parameter set and returns required variables, e.g. runoff.
44
45 !> \details Runs mhm with a specific parameter set and returns required variables, e.g. runoff.
46
47 ! INTENT(IN)
48 !> \param[in] "real(dp), dimension(:) :: parameterset" a set of global parameter (gamma) to run mHM, DIMENSION
49 !> [no. of global_Parameters]
50
51 ! INTENT(OUT), OPTIONAL
52 !> \param[out] "real(dp), dimension(:, :), optional :: runoff" returns runoff time series, DIMENSION
53 !> [nTimeSteps, nGaugesTotal]
54 !> \param[out] "real(dp), dimension(:, :), optional :: sm_opti" returns soil moisture time series for all
55 !> grid cells (of multiple Domains concatenated),DIMENSION [nCells, nTimeSteps]
56 !> time series, DIMENSION [nTimeSteps, nDomains]
57 !> \param[out] "real(dp), dimension(:, :), optional :: neutrons_opti" dim1=ncells, dim2=time
58 !> \param[out] "real(dp), dimension(:, :), optional :: et_opti" returns evapotranspiration time series for
59 !> \param[out] "real(dp), dimension(:, :), optional :: tws_opti" returns tws time series
60 !> all grid cells (of multiple Domains concatenated),DIMENSION [nCells, nTimeSteps]
61
62 ! HISTORY
63 !> \authors Juliane Mai, Rohini Kumar
64
65 !> \date Feb 2013
66
67 ! Modifications:
68 ! R. Kumar Jun 2013 - restart_flag_states_read is passed to mhm call for the soil moisture initalisation
69 ! R. Kumar Jun 2013 - frac_sealed_city_area is added
70 ! R. Kumar & S. Thober Aug 2013 - code change to incorporate output timestep during writing of the netcdf file
71 ! R. Kumar Aug 2013 - added iFlag_LAI_data_format to handle LAI options, and changed within the code made accordingly
72 ! R. Kumar, J. Mai Sep 2013 - Splitting allocation and initialization of arrays
73 ! R. Kumar Nov 2013 - update intent variables in documentation
74 ! L. Samaniego Nov 2013 - relational statements == to .eq., etc.
75 ! M. Zink Feb 2014 - added PET calculation: Hargreaves-Samani (Process 5)
76 ! M. Zink Mar 2014 - added inflow from upstream areas
77 ! Stephan Thober Jun 2014 - added chunk read for meteorological input
78 ! Stephan Thober Jun 2014 - updated flag for read_restart
79 ! M. Cuntz & J. Mai Nov 2014 - LAI input from daily, monthly or yearly files
80 ! Matthias Zink Dec 2014 - adopted inflow gauges to ignore headwater cells
81 ! Stephan Thober Aug 2015 - moved writing of daily discharge to mo_write_routing, included routing related variables from mRM
82 ! David Schaefer Aug 2015 - changed to new netcdf-writing scheme
83 ! Stephan Thober Sep 2015 - updated mrm_routing call
84 ! O. Rakovec, R. Kumar Oct 2015 - added optional output for Domain averaged TWS
85 ! Rohini Kumar Mar 2016 - changes for handling multiple soil database options
86 ! Stephan Thober Nov 2016 - added two options for routing
87 ! Rohini Kuamr Dec 2016 - option to handle monthly mean gridded fields of LAI
88 ! Stephan Thober Jan 2017 - added prescribed weights for tavg and pet
89 ! Zink M. Demirel C. Mar 2017 - Added Jarvis soil water stress function at SM process(3)
90 ! Robert Schweppe Dec 2017 - extracted call to mpr from inside mhm
91 ! Robert Schweppe Jun 2018 - refactoring and reformatting
92 ! Stephan Thober Jan 2022 - added nTstepForcingDay and is_hourly_forcing flag
93
94 SUBROUTINE mhm_eval(parameterset, opti_domain_indices, runoff, smOptiSim, neutronsOptiSim, etOptiSim, twsOptiSim, BFI)
95 implicit none
96
97 !> a set of global parameter (gamma) to run mHM, DIMENSION [no. of global_Parameters]
98 real(dp), dimension(:), intent(in) :: parameterset
99 !> selected domains for optimization
100 integer(i4), dimension(:), optional, intent(in) :: opti_domain_indices
101 !> returns runoff time series, DIMENSION [nTimeSteps, nGaugesTotal]
102 real(dp), dimension(:, :), allocatable, optional, intent(out) :: runoff
103 !> returns soil moisture time series for all grid cells (of multiple Domains concatenated),DIMENSION [nCells, nTimeSteps]
104 type(optidata_sim), dimension(:), optional, intent(inout) :: smoptisim
105 !> dim1=ncells, dim2=time
106 type(optidata_sim), dimension(:), optional, intent(inout) :: neutronsoptisim
107 !> returns evapotranspiration time series for all grid cells (of multiple Domains concatenated),DIMENSION [nCells, nTimeSteps]
108 type(optidata_sim), dimension(:), optional, intent(inout) :: etoptisim
109 !> returns tws time series for all grid cells (of multiple Domains concatenated),DIMENSION [nCells, nTimeSteps]
110 type(optidata_sim), dimension(:), optional, intent(inout) :: twsoptisim
111 !> baseflow index, dim1=domainID
112 real(dp), dimension(:), allocatable, optional, intent(out) :: bfi
113
114 ! number of domains simulated in this mhm_eval run. Depends on opti_function
115 integer(i4) :: ndomains, ii
116
117 ! flag to check if time loop is finished
118 logical :: time_loop_finished
119
120 ! prepare the mhm run
121 call mhm_interface_run_prepare(parameterset, opti_domain_indices, present(runoff), present(bfi))
122
123 ! get number of domains to loop over
124 call mhm_interface_run_get_ndomains(ndomains)
125
126 ! loop over Domains
127 domainloop: do ii = 1, ndomains
128
129 ! prepare current domain
130 call mhm_interface_run_prepare_domain(ii, etoptisim, twsoptisim, neutronsoptisim, smoptisim)
131
132 ! run time-loop at least once
133 time_loop_finished = .false.
134
135 ! Loop over time
136 timeloop: do while(.not. time_loop_finished)
137
138 ! do one time-step on current domain
140
141 ! write output
143
144 ! update optisim data
145 call mhm_interface_run_update_optisim(etoptisim, twsoptisim, neutronsoptisim, smoptisim)
146
147 ! check whether to run the time-loop further
148 call mhm_interface_run_finished(time_loop_finished)
149
150 end do timeloop !<< TIME STEPS LOOP
151
152 ! finalize domain
154
155 end do domainloop !<< Domain LOOP
156
157 ! SET RUNOFF OUTPUT VARIABLE; reset init-flag for MPR
158 call mhm_interface_run_finalize(runoff, bfi)
159
160 end SUBROUTINE mhm_eval
161
162END MODULE mo_mhm_eval
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_update_optisim(etoptisim, twsoptisim, neutronsoptisim, smoptisim)
add simulation data to optimization data types
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