5.13.3-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
wrapper.f90
Go to the documentation of this file.
1!> \file wrapper.f90
2!> \brief Module to wrap mHM with f2py to control it with Python.
3!> \authors Sebastian Mueller
4!> \date Nov 2022
5!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
6!! mHM is released under the LGPLv3+ license \license_note
7
8!> \brief Python wrapper module to control a mHM model.
9!> \date Nov 2022
10!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
11!! mHM is released under the LGPLv3+ license \license_note
12!> \ingroup mhm
13module model
14 use mo_kind, only : i4, dp
15 implicit none
16contains
17
18 !> \brief Configure coupling mode of mHM.
19 !> \details This will prevent reading the "coupling" namelist.
20 subroutine config_coupling( &
21 couple_case, &
22 meteo_timestep, &
23 meteo_time_ref_endpoint, &
24 meteo_expect_pre, &
25 meteo_expect_temp, &
26 meteo_expect_pet, &
27 meteo_expect_tmin, &
28 meteo_expect_tmax, &
29 meteo_expect_netrad, &
30 meteo_expect_absvappress, &
31 meteo_expect_windspeed, &
32 meteo_expect_ssrd, &
33 meteo_expect_strd, &
34 meteo_expect_tann &
35 )
37 implicit none
38 integer(i4), intent(in) :: couple_case !< coupling case
39 !f2py integer :: couple_case = 1
40 integer(i4), intent(in) :: meteo_timestep !< timestep for meteo-data from coupling
41 !f2py integer :: meteo_timestep = 24
42 logical, intent(in) :: meteo_time_ref_endpoint !< expect meteo has time reference point at end of time interval
43 !f2py logical :: meteo_time_ref_endpoint = 0
44 logical, intent(in) :: meteo_expect_pre !< expect meteo from coupling: [mm] Precipitation
45 !f2py logical :: meteo_expect_pre = 0
46 logical, intent(in) :: meteo_expect_temp !< expect meteo from coupling: [degC] Air temperature
47 !f2py logical :: meteo_expect_temp = 0
48 logical, intent(in) :: meteo_expect_pet !< expect meteo from coupling: [mm TS-1] Potential evapotranspiration
49 !f2py logical :: meteo_expect_pet = 0
50 logical, intent(in) :: meteo_expect_tmin !< expect meteo from coupling: [degC] minimum daily air temperature
51 !f2py logical :: meteo_expect_tmin = 0
52 logical, intent(in) :: meteo_expect_tmax !< expect meteo from coupling: [degC] maximum daily air temperature
53 !f2py logical :: meteo_expect_tmax = 0
54 logical, intent(in) :: meteo_expect_netrad !< expect meteo from coupling: [W m2] net radiation
55 !f2py logical :: meteo_expect_netrad = 0
56 logical, intent(in) :: meteo_expect_absvappress !< expect meteo from coupling: [Pa] absolute vapour pressure
57 !f2py logical :: meteo_expect_absvappress = 0
58 logical, intent(in) :: meteo_expect_windspeed !< expect meteo from coupling: [m s-1] windspeed
59 !f2py logical :: meteo_expect_windspeed = 0
60 logical, intent(in) :: meteo_expect_ssrd !< expect meteo from coupling: [W m2] short wave radiation
61 !f2py logical :: meteo_expect_ssrd = 0
62 logical, intent(in) :: meteo_expect_strd !< expect meteo from coupling: [W m2] long wave radiation
63 !f2py logical :: meteo_expect_strd = 0
64 logical, intent(in) :: meteo_expect_tann !< expect meteo from coupling: [degC] annual mean air temperature
65 !f2py logical :: meteo_expect_tann = 0
66 call couple_cfg%set_config( &
67 couple_case, &
68 meteo_timestep, &
69 meteo_time_ref_endpoint, &
70 meteo_expect_pre, &
71 meteo_expect_temp, &
72 meteo_expect_pet, &
73 meteo_expect_tmin, &
74 meteo_expect_tmax, &
75 meteo_expect_netrad, &
76 meteo_expect_absvappress, &
77 meteo_expect_windspeed, &
78 meteo_expect_ssrd, &
79 meteo_expect_strd, &
80 meteo_expect_tann &
81 )
82 end subroutine config_coupling
83
84 !> \brief Initialize a mHM model.
85 subroutine init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
87 implicit none
88 character(*), intent(in) :: namelist_mhm !< path to mHM configuration namelist
89 !f2py character(*) :: namelist_mhm = "mhm.nml"
90 character(*), intent(in) :: namelist_mhm_param !< path to mHM parameter namelist
91 !f2py character(*) :: namelist_mhm_param = "mhm_parameter.nml"
92 character(*), intent(in) :: namelist_mhm_output !< path to mHM output namelist
93 !f2py character(*) :: namelist_mhm_output = "mhm_outputs.nml"
94 character(*), intent(in) :: namelist_mrm_output !< path to mRM output namelist
95 !f2py character(*) :: namelist_mrm_output = "mrm_outputs.nml"
96 character(*), intent(in) :: cwd !< desired working directory
97 !f2py character(*) :: cwd = "."
98 call mhm_interface_init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
99 end subroutine init
100
101 !> \brief Execute a mHM model.
102 subroutine run()
104 implicit none
106 end subroutine run
107
108 !> \brief disable all mHM/mRM outputs during runtime.
109 subroutine disable_output()
112 implicit none
113 outputflxstate_mrm = .false.
114 outputflxstate = .false.
115 end subroutine disable_output
116
117 !> \brief Execute a mHM model with given parameters.
118 subroutine run_with_parameter(parameter, n)
119 use mo_mhm_eval, only: mhm_eval
120 implicit none
121 integer(i4) :: n !< number of parameters
122 real(dp), intent(in) :: parameter(n) !< parameter
123 call mhm_eval(parameter)
124 end subroutine run_with_parameter
125
126 !> \brief Execute a mHM model or an optimization depending on the configuration.
127 subroutine run_or_optimize()
129 use mo_mhm_interface, only: &
132 implicit none
133 ! RUN OR OPTIMIZE
134 if (optimize) then
136 else
137 ! single mhm run with current settings
139 end if
140 end subroutine run_or_optimize
141
142 !> \brief Get the mHM version.
143 subroutine version(ver_string)
144 use mo_file, only: mhm_version => version
145 implicit none
146 character(64), intent(out) :: ver_string
147 ver_string = mhm_version
148 end subroutine version
149
150 !> \brief Finalize a mHM model.
151 subroutine finalize()
153 implicit none
155 end subroutine finalize
156
157 !> \brief Set verbosity level of mHM.
158 subroutine set_verbosity(level)
160 implicit none
161 integer(i4), intent(in) :: level !< verbosity level (0, 1, 2)
162 call set_verbosity_level(level)
163 end subroutine set_verbosity
164end module model
165
166!> \brief Python wrapper module to control a mHM model run per time step.
167!> \date Nov 2022
168!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
169!! mHM is released under the LGPLv3+ license \license_note
170!> \ingroup mhm
171module run
172 use mo_kind, only : i4, dp
173 implicit none
174contains
175 !> \brief Prepare a mHM model run.
176 subroutine prepare()
178 implicit none
180 end subroutine prepare
181
182 !> \brief Get the number of domains of the current mHM model run.
183 subroutine get_ndomains(n)
185 implicit none
186 integer(i4), intent(out) :: n !< number of domains
188 end subroutine get_ndomains
189
190 !> \brief Prepare a certain domain of the current mHM model run.
191 subroutine prepare_domain(domain)
193 implicit none
194 integer(i4), intent(in) :: domain !< domain index (1 based and 1 by default)
195 !f2py integer :: domain = 1
197 end subroutine prepare_domain
198
199 !> \brief Check if the current mHM model time loop is finished.
200 subroutine finished(output)
202 implicit none
203 logical, intent(out) :: output !< whether the current time loop finished
204 call mhm_interface_run_finished(output)
205 end subroutine finished
206
207 !> \brief Do one time-step on the current domain of the current mHM model run.
212 end subroutine do_time_step
213
214 !> \brief Write output for the current domain of the current mHM model run.
219 end subroutine write_output
220
221 !> \brief Finalize the current domain of the current mHM model run.
226 end subroutine finalize_domain
227
228 !> \brief Finalize the current mHM model run.
229 subroutine finalize()
231 implicit none
233 end subroutine finalize
234
235 !> \brief Get the current time the current domain of the current mHM model run.
236 subroutine current_time(year, month, day, hour)
238 implicit none
239 integer(i4), intent(out) :: year !< current year
240 integer(i4), intent(out) :: month !< current month
241 integer(i4), intent(out) :: day !< current day
242 integer(i4), intent(out) :: hour !< current 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
247 end subroutine current_time
248end module run
249
250!> \brief Python wrapper module to get internal variables of a mHM model run.
251!> \date Nov 2022
252!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
253!! mHM is released under the LGPLv3+ license \license_note
254!> \ingroup mhm
255module get
256 use mo_kind, only : i4, dp
257 implicit none
258contains
259
260 !> \name constants
261 !> \brief access constants of mHM
262
263 !> \brief Get the number of soil horizons in mHM.
264 subroutine number_of_horizons(n)
266 implicit none
267 integer(i4), intent(out) :: n !< number of soil horizons
269 end subroutine number_of_horizons
270
271 !> \name parameter
272
273 !> \brief Get the number of parameters in mHM.
274 subroutine parameter_length(length)
276 implicit none
277 integer(i4), intent(out) :: length !< length of the parmeter array
279 end subroutine parameter_length
280
281 !> \name parameter
282
283 !> \brief Get the parameter settings of mHM.
284 subroutine parameter_config(config, n)
286
287 implicit none
288 integer(i4) :: n !< number of parameters
289 real(dp), intent(out) :: config(n, 5) !< parameter configuration
290
291 config = global_parameters
292 end subroutine parameter_config
293
294 !> \name parameter
295 !> \brief access parameters of mHM
296
297 !> \brief Get the parameter names of mHM.
298 subroutine parameter_name(para_name, n)
300
301 implicit none
302 character(256), intent(out) :: para_name !< parameter name
303 !f2py character(f2py_len=256) para_name
304 integer(i4), intent(in) :: n !< number of parameters
305
306 para_name = global_parameters_name(n)
307 end subroutine parameter_name
308
309 !> \name runoff
310
311 !> \brief Get the shape of mHM model runoff output for evaluation.
312 subroutine runoff_eval_length(gauge_id, length)
315 implicit none
316 integer(i4), intent(in) :: gauge_id !< gauge id
317 integer(i4), intent(out) :: length !< length of the runoff time series
318
319 integer(i4) :: iDomain
320
321 ! extract domain Id from gauge Id
322 idomain = gauge%domainId(gauge_id)
323 ! get length of evaluation period times TPD_obs
324 length = (evalper(idomain)%julEnd - evalper(idomain)%julStart + 1) * nmeasperday
325 end subroutine runoff_eval_length
326
327 !> \brief Get the mHM model runoff output.
328 subroutine runoff_eval(gauge_id, output, m)
331 implicit none
332 integer(i4) :: m !< number of time-steps for evaluation
333 real(dp), intent(out) :: output(m, 2) !< runoff
334 integer(i4), intent(in) :: gauge_id !< gauge id
335
336 ! aggregated simulated runoff
337 real(dp), dimension(:), allocatable :: runoff_agg, runoff_obs
338 logical, dimension(:), allocatable :: runoff_obs_mask
339
340 call extract_runoff(gauge_id, mrm_runoff, runoff_agg, runoff_obs, runoff_obs_mask)
341 output(:, 1) = runoff_agg
342 output(:, 2) = runoff_obs
343 end subroutine runoff_eval
344
345 !> \name runoff
346
347 !> \brief Get the shape of mHM model runoff output.
348 subroutine runoff_shape(shp)
350 implicit none
351 integer(i4), intent(out) :: shp(2) !< 2D shape of the runoff
352 shp = shape(mrm_runoff)
353 end subroutine runoff_shape
354
355 !> \name runoff
356 !> \brief access generated runoff of mHM
357
358 !> \brief Get the mHM model runoff output.
359 subroutine runoff(output, m, n)
361 implicit none
362 integer(i4) :: m !< number of time-steps
363 integer(i4) :: n !< number of gauges
364 real(dp), intent(out) :: output(m, n) !< runoff
365 output = mrm_runoff
366 end subroutine runoff
367
368 !> \name Level 0
369
370 !> \brief Get number of unmasked celles on Level-0 of the mHM model.
371 subroutine l0_domain_size(n, domain)
372 use mo_common_variables, only : level0
374 implicit none
375 integer(i4), intent(out) :: n !< number of unmasked celles
376 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
377 !f2py integer :: domain = 0
378 integer(i4) :: iDomain, i
379 i = domain
380 if ( i == 0 ) i = run_cfg%selected_domain
381 idomain = run_cfg%get_domain_index(i)
382 n = level0(idomain)%nCells
383 end subroutine l0_domain_size
384
385 !> \name Level 0
386
387 !> \brief Get the shape of Level-0 of the mHM model.
388 subroutine l0_domain_shape(shp, domain)
389 use mo_common_variables, only : level0
391 implicit none
392 integer(i4), intent(out) :: shp(2) !< shape of Level-0
393 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
394 !f2py integer :: domain = 0
395 integer(i4) :: iDomain, i
396 i = domain
397 if ( i == 0 ) i = run_cfg%selected_domain
398 idomain = run_cfg%get_domain_index(i)
399 shp = shape(level0(idomain)%mask)
400 end subroutine l0_domain_shape
401
402 !> \name Level 0
403
404 !> \brief Get the mask of Level-0 of the mHM model.
405 subroutine l0_domain_mask(mask, n, m, domain)
406 use mo_common_variables, only : level0
408 implicit none
409 integer(i4) :: m !< number of columns
410 integer(i4) :: n !< number of rows
411 logical, intent(out) :: mask(m, n) !< mask at Level-0
412 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
413 !f2py integer :: domain = 0
414 integer(i4) :: iDomain, i
415 i = domain
416 if ( i == 0 ) i = run_cfg%selected_domain
417 idomain = run_cfg%get_domain_index(i)
418 mask = level0(idomain)%mask
419 end subroutine l0_domain_mask
420
421 !> \name Level 0
422
423 !> \brief Get the information of Level-0 of the mHM model.
424 subroutine l0_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
425 use mo_common_variables, only : level0
428 implicit none
429 integer(i4), intent(out) :: ncols !< number of columns
430 integer(i4), intent(out) :: nrows !< number of rows
431 integer(i4), intent(out) :: ncells !< number of cells
432 real(dp), intent(out) :: xll !< x coordinate of lower left corner
433 real(dp), intent(out) :: yll !< y coordinate of lower left corner
434 real(dp), intent(out) :: cell_size !< cell-size
435 real(dp), intent(out) :: no_data !< no data value
436 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
437 !f2py integer :: domain = 0
438 integer(i4) :: iDomain, i
439 i = domain
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
448 ! no_data = level0(iDomain)%nodata_value
449 no_data = nodata_dp
450 end subroutine l0_domain_info
451
452 !> \name Level 0
453 !> \brief access Level 0 information and variables
454
455 !> \brief Get a variable on Level-0 of the mHM model.
456 subroutine l0_variable(output, name, idx, n)
457 use mo_message, only : error_message
461 implicit none
462 integer(i4) :: n !< size of the variable
463 real(dp), intent(out) :: output(n) !< the desired variable
464 character(*), intent(in) :: name !< name to select the variable
465 integer(i4), intent(in) :: idx !< optional index if the variable has multiple layer (1 by default)
466 !f2py integer :: idx = 1
467 integer(i4) :: iDomain, s0, e0
468
469 idomain = run_cfg%get_domain_index(run_cfg%selected_domain)
470 s0 = level0(domainmeta%L0DataFrom(idomain))%iStart
471 e0 = level0(domainmeta%L0DataFrom(idomain))%iEnd
472 select case(name)
473 case("L0_GRIDDED_LAI")
474 output = l0_gridded_lai(s0 : e0, idx)
475 case default
476 call error_message("mhm.get.L0_variable: unkown variable name: ", name)
477 end select
478 end subroutine l0_variable
479
480 !> \name Level 1
481
482 !> \brief Get number of unmasked celles on Level-1 of the mHM model.
483 subroutine l1_domain_size(n, domain)
484 use mo_common_variables, only : level1
486 implicit none
487 integer(i4), intent(out) :: n !< number of unmasked celles
488 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
489 !f2py integer :: domain = 0
490 integer(i4) :: iDomain, i
491 i = domain
492 if ( i == 0 ) i = run_cfg%selected_domain
493 idomain = run_cfg%get_domain_index(i)
494 n = level1(idomain)%nCells
495 end subroutine l1_domain_size
496
497 !> \name Level 1
498
499 !> \brief Get the shape of Level-1 of the mHM model.
500 subroutine l1_domain_shape(shp, domain)
501 use mo_common_variables, only : level1
503 implicit none
504 integer(i4), intent(out) :: shp(2) !< shape of Level-1
505 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
506 !f2py integer :: domain = 0
507 integer(i4) :: iDomain, i
508 i = domain
509 if ( i == 0 ) i = run_cfg%selected_domain
510 idomain = run_cfg%get_domain_index(i)
511 shp = shape(level1(idomain)%mask)
512 end subroutine l1_domain_shape
513
514 !> \name Level 1
515
516 !> \brief Get the mask of Level-1 of the mHM model.
517 subroutine l1_domain_mask(mask, n, m, domain)
518 use mo_common_variables, only : level1
520 implicit none
521 integer(i4) :: m !< number of columns
522 integer(i4) :: n !< number of rows
523 logical, intent(out) :: mask(m, n) !< mask at Level-1
524 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
525 !f2py integer :: domain = 0
526 integer(i4) :: iDomain, i
527 i = domain
528 if ( i == 0 ) i = run_cfg%selected_domain
529 idomain = run_cfg%get_domain_index(i)
530 mask = level1(idomain)%mask
531 end subroutine l1_domain_mask
532
533 !> \name Level 1
534 !> \brief access Level 1 information and variables
535
536 !> \brief Get the information of Level-1 of the mHM model.
537 subroutine l1_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
538 use mo_common_variables, only : level1
541 implicit none
542 integer(i4), intent(out) :: ncols !< number of columns
543 integer(i4), intent(out) :: nrows !< number of rows
544 integer(i4), intent(out) :: ncells !< number of cells
545 real(dp), intent(out) :: xll !< x coordinate of lower left corner
546 real(dp), intent(out) :: yll !< y coordinate of lower left corner
547 real(dp), intent(out) :: cell_size !< cell-size
548 real(dp), intent(out) :: no_data !< no data value
549 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
550 !f2py integer :: domain = 0
551 integer(i4) :: iDomain, i
552 i = domain
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
561 ! no_data = level1(iDomain)%nodata_value
562 no_data = nodata_dp
563 end subroutine l1_domain_info
564
565 !> \name Level 1
566 !> \brief access Level 1 information and variables
567
568 !> \brief Get a variable on Level-1 of the mHM model.
569 subroutine l1_variable(output, name, idx, n)
570 use mo_message, only : error_message
572 use mo_mpr_global_variables, only : &
573 l1_fsealed, & ! 1d
574 l1_soilmoistsat ! 2d
575 use mo_global_variables, only : &
576 l1_inter, & ! 1d
577 l1_snowpack, & ! 1d
578 l1_soilmoist, & ! 2d
579 l1_sealstw, & ! 1d
580 l1_unsatstw, & ! 1d
581 l1_satstw, & ! 1d
582 l1_neutrons, & ! 1d
583 l1_pet_calc, & ! 1d
584 l1_temp_calc, & ! 1d
585 l1_prec_calc, & ! 1d
586 l1_aetsoil, & ! 2d
587 l1_aetcanopy, & ! 1d
588 l1_aetsealed, & ! 1d
589 l1_total_runoff, & ! 1d
590 l1_runoffseal, & ! 1d
591 l1_fastrunoff, & ! 1d
592 l1_slowrunoff, & ! 1d
593 l1_baseflow, & ! 1d
594 l1_percol, & ! 1d
595 l1_infilsoil, & ! 2d
596 l1_preeffect, & ! 1d
597 l1_melt, & ! 1d
598 l1_rain, & ! 1d
599 l1_snow, & ! 1d
600 l1_throughfall ! 1d
601 implicit none
602 integer(i4) :: n !< size of the variable
603 real(dp), intent(out) :: output(n) !< the desired variable
604 character(*), intent(in) :: name !< name to select the variable
605 integer(i4), intent(in) :: idx !< optional index if the variable has multiple layer (1 by default)
606 !f2py integer :: idx = 1
607 select case(name)
608 case("L1_FSEALED")
609 output = l1_fsealed(run_cfg%s1 : run_cfg%e1, 1, run_cfg%domainDateTime%yId)
610 case("L1_FNOTSEALED")
611 output = run_cfg%L1_fNotSealed(run_cfg%s1 : run_cfg%e1, 1, run_cfg%domainDateTime%yId)
612 case("L1_INTER")
613 output = l1_inter(run_cfg%s1 : run_cfg%e1)
614 case("L1_SNOWPACK")
615 output = l1_snowpack(run_cfg%s1 : run_cfg%e1)
616 case("L1_SEALSTW")
617 output = l1_sealstw(run_cfg%s1 : run_cfg%e1)
618 case("L1_UNSATSTW")
619 output = l1_unsatstw(run_cfg%s1 : run_cfg%e1)
620 case("L1_SATSTW")
621 output = l1_satstw(run_cfg%s1 : run_cfg%e1)
622 case("L1_NEUTRONS")
623 output = l1_neutrons(run_cfg%s1 : run_cfg%e1)
624 case("L1_PET_CALC")
625 output = l1_pet_calc(run_cfg%s1 : run_cfg%e1)
626 case("L1_TEMP_CALC")
627 output = l1_temp_calc(run_cfg%s1 : run_cfg%e1)
628 case("L1_PREC_CALC")
629 output = l1_prec_calc(run_cfg%s1 : run_cfg%e1)
630 case("L1_AETCANOPY")
631 output = l1_aetcanopy(run_cfg%s1 : run_cfg%e1)
632 case("L1_AETSEALED")
633 output = l1_aetsealed(run_cfg%s1 : run_cfg%e1)
634 case("L1_TOTAL_RUNOFF")
635 output = l1_total_runoff(run_cfg%s1 : run_cfg%e1)
636 case("L1_RUNOFFSEAL")
637 output = l1_runoffseal(run_cfg%s1 : run_cfg%e1)
638 case("L1_FASTRUNOFF")
639 output = l1_fastrunoff(run_cfg%s1 : run_cfg%e1)
640 case("L1_SLOWRUNOFF")
641 output = l1_slowrunoff(run_cfg%s1 : run_cfg%e1)
642 case("L1_BASEFLOW")
643 output = l1_baseflow(run_cfg%s1 : run_cfg%e1)
644 case("L1_PERCOL")
645 output = l1_percol(run_cfg%s1 : run_cfg%e1)
646 case("L1_PREEFFECT")
647 output = l1_preeffect(run_cfg%s1 : run_cfg%e1)
648 case("L1_SOILMOIST")
649 output = l1_soilmoist(run_cfg%s1 : run_cfg%e1, idx)
650 case("L1_SOILMOIST_VOL")
651 output = l1_soilmoist(run_cfg%s1 : run_cfg%e1, idx) &
652 / l1_soilmoistsat(run_cfg%s1 : run_cfg%e1, idx, run_cfg%domainDateTime%yId)
653 case("L1_SOILMOIST_VOL_ALL")
654 output = sum(l1_soilmoist(run_cfg%s1 : run_cfg%e1, :), dim = 2) &
655 / sum(l1_soilmoistsat(run_cfg%s1 : run_cfg%e1, :, run_cfg%domainDateTime%yId), dim = 2)
656 case("L1_SOILMOISTSAT")
657 output = l1_soilmoistsat(run_cfg%s1 : run_cfg%e1, idx, run_cfg%domainDateTime%yId)
658 case("L1_AETSOIL")
659 output = l1_aetsoil(run_cfg%s1 : run_cfg%e1, idx)
660 case("L1_INFILSOIL")
661 output = l1_infilsoil(run_cfg%s1 : run_cfg%e1, idx)
662 case("L1_MELT")
663 output = l1_melt(run_cfg%s1 : run_cfg%e1)
664 case("L1_RAIN")
665 output = l1_rain(run_cfg%s1 : run_cfg%e1)
666 case("L1_SNOW")
667 output = l1_snow(run_cfg%s1 : run_cfg%e1)
668 case("L1_THROUGHFALL")
669 output = l1_throughfall(run_cfg%s1 : run_cfg%e1)
670 case default
671 call error_message("mhm.get.L1_variable: unkown variable name: ", name)
672 end select
673 end subroutine l1_variable
674
675 !> \name Level 11
676
677 !> \brief Get number of unmasked celles on Level-11 of the mHM model.
678 subroutine l11_domain_size(n, domain)
681 implicit none
682 integer(i4), intent(out) :: n !< number of unmasked celles
683 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
684 !f2py integer :: domain = 0
685 integer(i4) :: iDomain, i
686 i = domain
687 if ( i == 0 ) i = run_cfg%selected_domain
688 idomain = run_cfg%get_domain_index(i)
689 n = level11(idomain)%nCells
690 end subroutine l11_domain_size
691
692 !> \name Level 11
693
694 !> \brief Get the shape of Level-11 of the mHM model.
695 subroutine l11_domain_shape(shp, domain)
698 implicit none
699 integer(i4), intent(out) :: shp(2) !< shape of Level-11
700 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
701 !f2py integer :: domain = 0
702 integer(i4) :: iDomain, i
703 i = domain
704 if ( i == 0 ) i = run_cfg%selected_domain
705 idomain = run_cfg%get_domain_index(i)
706 shp = shape(level11(idomain)%mask)
707 end subroutine l11_domain_shape
708
709 !> \name Level 11
710
711 !> \brief Get the mask of Level-11 of the mHM model.
712 subroutine l11_domain_mask(mask, n, m, domain)
715 implicit none
716 integer(i4) :: m !< number of columns
717 integer(i4) :: n !< number of rows
718 logical, intent(out) :: mask(m, n) !< mask at Level-11
719 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
720 !f2py integer :: domain = 0
721 integer(i4) :: iDomain, i
722 i = domain
723 if ( i == 0 ) i = run_cfg%selected_domain
724 idomain = run_cfg%get_domain_index(i)
725 mask = level11(idomain)%mask
726 end subroutine l11_domain_mask
727
728 !> \name Level 11
729
730 !> \brief Get the information of Level-11 of the mHM model.
731 subroutine l11_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
735 implicit none
736 integer(i4), intent(out) :: ncols !< number of columns
737 integer(i4), intent(out) :: nrows !< number of rows
738 integer(i4), intent(out) :: ncells !< number of cells
739 real(dp), intent(out) :: xll !< x coordinate of lower left corner
740 real(dp), intent(out) :: yll !< y coordinate of lower left corner
741 real(dp), intent(out) :: cell_size !< cell-size
742 real(dp), intent(out) :: no_data !< no data value
743 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
744 !f2py integer :: domain = 0
745 integer(i4) :: iDomain, i
746 i = domain
747 if ( i == 0 ) i = run_cfg%selected_domain
748 idomain = run_cfg%get_domain_index(i)
749 ncols = level11(idomain)%ncols
750 nrows = level11(idomain)%nrows
751 ncells = level11(idomain)%nCells
752 xll = level11(idomain)%xllcorner
753 yll = level11(idomain)%yllcorner
754 cell_size = level11(idomain)%cellsize
755 ! no_data = level11(iDomain)%nodata_value
756 no_data = nodata_dp
757 end subroutine l11_domain_info
758
759 !> \name Level 11
760 !> \brief access Level 11 information and variables
761
762 !> \brief Get a variable on Level-11 of the mHM model.
763 subroutine l11_variable(output, name, idx, n)
764 use mo_message, only : error_message
766 use mo_mrm_global_variables, only : &
767 l11_qmod, &
768 l11_qout, &
769 l11_qtin, &
770 l11_qtr
771 implicit none
772 integer(i4) :: n !< size of the variable
773 real(dp), intent(out) :: output(n) !< the desired variable
774 character(*), intent(in) :: name !< name to select the variable
775 integer(i4), intent(in) :: idx !< optional index if the variable has multiple layer (1 by default)
776 !f2py integer :: idx = 1
777 select case(name)
778 case("L11_QMOD")
779 output = l11_qmod(run_cfg%s11 : run_cfg%e11)
780 case("L11_QOUT")
781 output = l11_qout(run_cfg%s11 : run_cfg%e11)
782 case("L11_QTIN")
783 output = l11_qtin(run_cfg%s11 : run_cfg%e11, idx)
784 case("L11_QTR")
785 output = l11_qtr(run_cfg%s11 : run_cfg%e11, idx)
786 case default
787 call error_message("mhm.get.L11_variable: unkown variable name: ", name)
788 end select
789 end subroutine l11_variable
790
791 !> \name Level 2
792
793 !> \brief Get number of unmasked celles on Level-2 of the mHM model.
794 subroutine l2_domain_size(n, domain)
797 implicit none
798 integer(i4), intent(out) :: n !< number of unmasked celles
799 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
800 !f2py integer :: domain = 0
801 integer(i4) :: iDomain, i
802 i = domain
803 if ( i == 0 ) i = run_cfg%selected_domain
804 idomain = run_cfg%get_domain_index(i)
805 n = meteo_handler%level2(idomain)%nCells
806 end subroutine l2_domain_size
807
808 !> \name Level 2
809
810 !> \brief Get the shape of Level-2 of the mHM model.
811 subroutine l2_domain_shape(shp, domain)
814 implicit none
815 integer(i4), intent(out) :: shp(2) !< shape of Level-2
816 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
817 !f2py integer :: domain = 0
818 integer(i4) :: iDomain, i
819 i = domain
820 if ( i == 0 ) i = run_cfg%selected_domain
821 idomain = run_cfg%get_domain_index(i)
822 shp = shape(meteo_handler%level2(idomain)%mask)
823 end subroutine l2_domain_shape
824
825 !> \name Level 2
826
827 !> \brief Get the mask of Level-2 of the mHM model.
828 subroutine l2_domain_mask(mask, n, m, domain)
831 implicit none
832 integer(i4) :: m !< number of columns
833 integer(i4) :: n !< number of rows
834 logical, intent(out) :: mask(m, n) !< mask at Level-2
835 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
836 !f2py integer :: domain = 0
837 integer(i4) :: iDomain, i
838 i = domain
839 if ( i == 0 ) i = run_cfg%selected_domain
840 idomain = run_cfg%get_domain_index(i)
841 mask = meteo_handler%level2(idomain)%mask
842 end subroutine l2_domain_mask
843
844 !> \name Level 2
845 !> \brief access Level 2 information
846
847 !> \brief Get the information of Level-2 of the mHM model.
848 subroutine l2_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
852 implicit none
853 integer(i4), intent(out) :: ncols !< number of columns
854 integer(i4), intent(out) :: nrows !< number of rows
855 integer(i4), intent(out) :: ncells !< number of cells
856 real(dp), intent(out) :: xll !< x coordinate of lower left corner
857 real(dp), intent(out) :: yll !< y coordinate of lower left corner
858 real(dp), intent(out) :: cell_size !< cell-size
859 real(dp), intent(out) :: no_data !< no data value
860 integer(i4), intent(in) :: domain !< selected domain (0 by default for current domain)
861 !f2py integer :: domain = 0
862 integer(i4) :: iDomain, i
863 i = domain
864 if ( i == 0 ) i = run_cfg%selected_domain
865 idomain = run_cfg%get_domain_index(i)
866 ncols = meteo_handler%level2(idomain)%ncols
867 nrows = meteo_handler%level2(idomain)%nrows
868 ncells = meteo_handler%level2(idomain)%nCells
869 xll = meteo_handler%level2(idomain)%xllcorner
870 yll = meteo_handler%level2(idomain)%yllcorner
871 cell_size = meteo_handler%level2(idomain)%cellsize
872 ! no_data = meteo_handler%level2(iDomain)%nodata_value
873 no_data = nodata_dp
874 end subroutine l2_domain_info
875
876end module get
877
878!> \brief Python wrapper module to set internal variables of a mHM model run.
879!> \date Nov 2022
880!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
881!! mHM is released under the LGPLv3+ license \license_note
882!> \ingroup mhm
883module set
884 use mo_kind, only : i4, dp
885 implicit none
886contains
887
888 !> \name Level 0
889 !> \brief alter Level 0 variables
890
891 !> \brief Set a variable on Level-0 of the mHM model.
892 subroutine l0_variable(input, name, idx, n)
893 use mo_message, only: error_message
897 implicit none
898 integer(i4) :: n !< size of the variable
899 real(dp), intent(in) :: input(n) !< the variable value
900 character(*), intent(in) :: name !< name to select the variable
901 integer(i4), intent(in) :: idx !< optional index if the variable has multiple layer (1 by default)
902 !f2py integer :: idx = 1
903 integer(i4) :: iDomain, s0, e0
904
905 idomain = run_cfg%get_domain_index(run_cfg%selected_domain)
906 s0 = level0(domainmeta%L0DataFrom(idomain))%iStart
907 e0 = level0(domainmeta%L0DataFrom(idomain))%iEnd
908 select case(name)
909 case("L0_GRIDDED_LAI")
910 l0_gridded_lai(s0 : e0, idx) = input
911 case default
912 call error_message("mhm.set.L0_variable: unkown variable name: ", name)
913 end select
914 end subroutine l0_variable
915
916 !> \name Meteo
917 !> \brief set meteo values
918
919 !> \brief Set a meteo variable on Level-1 of the mHM model.
920 subroutine meteo(input, name, year, month, day, hour, n)
921 use mo_message, only: error_message
925 implicit none
926 integer(i4) :: n !< size of the variable
927 real(dp), intent(in) :: input(n) !< the variable value
928 character(*), intent(in) :: name !< name to select the variable
929 integer(i4), intent(in) :: year !< year of the input time
930 integer(i4), intent(in) :: month !< month of the input time
931 integer(i4), intent(in) :: day !< day of the input time
932 integer(i4), intent(in) :: hour !< optional hour of the input time
933 !f2py integer :: hour = 0
934 select case(name)
935 case("PRE")
936 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, pre=input)
937 case("TEMP")
938 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, temp=input)
939 case("PET")
940 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, pet=input)
941 case("TMIN")
942 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tmin=input)
943 case("TMAX")
944 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tmax=input)
945 case("NETRAD")
946 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, netrad=input)
947 case("ABSVAPPRESS")
948 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, absvappress=input)
949 case("WINDSPEED")
950 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, windspeed=input)
951 case("SSRD")
952 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, ssrd=input)
953 case("STRD")
954 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, strd=input)
955 case("TANN")
956 call meteo_handler%set_meteo(year=year, month=month, day=day, hour=hour, tann=input)
957 case default
958 call error_message("mhm.set.meteo: unkown variable name: ", name)
959 end select
960 end subroutine meteo
961end module set
Python wrapper module to get internal variables of a mHM model run.
Definition wrapper.f90:255
subroutine runoff_shape(shp)
Get the shape of mHM model runoff output.
Definition wrapper.f90:349
subroutine runoff_eval_length(gauge_id, length)
Get the shape of mHM model runoff output for evaluation.
Definition wrapper.f90:313
subroutine l11_domain_shape(shp, domain)
Get the shape of Level-11 of the mHM model.
Definition wrapper.f90:696
subroutine runoff(output, m, n)
Get the mHM model runoff output.
Definition wrapper.f90:360
subroutine parameter_config(config, n)
Get the parameter settings of mHM.
Definition wrapper.f90:285
subroutine l2_domain_size(n, domain)
Get number of unmasked celles on Level-2 of the mHM model.
Definition wrapper.f90:795
subroutine l11_domain_mask(mask, n, m, domain)
Get the mask of Level-11 of the mHM model.
Definition wrapper.f90:713
subroutine l11_domain_size(n, domain)
Get number of unmasked celles on Level-11 of the mHM model.
Definition wrapper.f90:679
subroutine l2_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-2 of the mHM model.
Definition wrapper.f90:849
subroutine l11_variable(output, name, idx, n)
Get a variable on Level-11 of the mHM model.
Definition wrapper.f90:764
subroutine l1_domain_shape(shp, domain)
Get the shape of Level-1 of the mHM model.
Definition wrapper.f90:501
subroutine parameter_length(length)
Get the number of parameters in mHM.
Definition wrapper.f90:275
subroutine l1_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-1 of the mHM model.
Definition wrapper.f90:538
subroutine parameter_name(para_name, n)
Get the parameter names of mHM.
Definition wrapper.f90:299
subroutine l0_domain_size(n, domain)
Get number of unmasked celles on Level-0 of the mHM model.
Definition wrapper.f90:372
subroutine l0_domain_shape(shp, domain)
Get the shape of Level-0 of the mHM model.
Definition wrapper.f90:389
subroutine l2_domain_mask(mask, n, m, domain)
Get the mask of Level-2 of the mHM model.
Definition wrapper.f90:829
subroutine l1_variable(output, name, idx, n)
Get a variable on Level-1 of the mHM model.
Definition wrapper.f90:570
subroutine l0_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-0 of the mHM model.
Definition wrapper.f90:425
subroutine number_of_horizons(n)
Get the number of soil horizons in mHM.
Definition wrapper.f90:265
subroutine l11_domain_info(ncols, nrows, ncells, xll, yll, cell_size, no_data, domain)
Get the information of Level-11 of the mHM model.
Definition wrapper.f90:732
subroutine l1_domain_mask(mask, n, m, domain)
Get the mask of Level-1 of the mHM model.
Definition wrapper.f90:518
subroutine l2_domain_shape(shp, domain)
Get the shape of Level-2 of the mHM model.
Definition wrapper.f90:812
subroutine runoff_eval(gauge_id, output, m)
Get the mHM model runoff output.
Definition wrapper.f90:329
subroutine l0_variable(output, name, idx, n)
Get a variable on Level-0 of the mHM model.
Definition wrapper.f90:457
subroutine l1_domain_size(n, domain)
Get number of unmasked celles on Level-1 of the mHM model.
Definition wrapper.f90:484
subroutine l0_domain_mask(mask, n, m, domain)
Get the mask of Level-0 of the mHM model.
Definition wrapper.f90:406
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), 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.
Definition mo_file.F90:29
character(len=*), parameter version
Current mHM model version (will be set to )
Definition mo_file.F90:33
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
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.
Definition wrapper.f90:13
subroutine version(ver_string)
Get the mHM version.
Definition wrapper.f90:144
subroutine set_verbosity(level)
Set verbosity level of mHM.
Definition wrapper.f90:159
subroutine disable_output()
disable all mHM/mRM outputs during runtime.
Definition wrapper.f90:110
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.
Definition wrapper.f90:36
subroutine run_or_optimize()
Execute a mHM model or an optimization depending on the configuration.
Definition wrapper.f90:128
subroutine run_with_parameter(parameter, n)
Execute a mHM model with given parameters.
Definition wrapper.f90:119
subroutine init(namelist_mhm, namelist_mhm_param, namelist_mhm_output, namelist_mrm_output, cwd)
Initialize a mHM model.
Definition wrapper.f90:86
subroutine finalize()
Finalize a mHM model.
Definition wrapper.f90:152
Python wrapper module to control a mHM model run per time step.
Definition wrapper.f90:171
subroutine write_output()
Write output for the current domain of the current mHM model run.
Definition wrapper.f90:216
subroutine prepare()
Prepare a mHM model run.
Definition wrapper.f90:177
subroutine prepare_domain(domain)
Prepare a certain domain of the current mHM model run.
Definition wrapper.f90:192
subroutine finished(output)
Check if the current mHM model time loop is finished.
Definition wrapper.f90:201
subroutine current_time(year, month, day, hour)
Get the current time the current domain of the current mHM model run.
Definition wrapper.f90:237
subroutine finalize_domain()
Finalize the current domain of the current mHM model run.
Definition wrapper.f90:223
subroutine do_time_step()
Do one time-step on the current domain of the current mHM model run.
Definition wrapper.f90:209
subroutine get_ndomains(n)
Get the number of domains of the current mHM model run.
Definition wrapper.f90:184
subroutine finalize()
Finalize the current mHM model run.
Definition wrapper.f90:230
Python wrapper module to set internal variables of a mHM model run.
Definition wrapper.f90:883
subroutine meteo(input, name, year, month, day, hour, n)
Set a meteo variable on Level-1 of the mHM model.
Definition wrapper.f90:921
subroutine l0_variable(input, name, idx, n)
Set a variable on Level-0 of the mHM model.
Definition wrapper.f90:893