5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_optimization_utils.f90
Go to the documentation of this file.
1!> \file mo_optimization_utils.f90
2!> \copydoc mo_optimization_utils
3
4!> \brief Utility functions, such as interface definitions, for optimization routines.
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!> \ingroup f_common
9
10 use mo_kind, only : dp
11 use mo_optimizee, only: optimizee
12 use mo_message, only : error_message
13
14 implicit none
15
16 !> \brief Interface for evaluation routine.
17 abstract interface
18 subroutine eval_interface(parameterset, opti_domain_indices, runoff, smOptiSim, neutronsOptiSim, etOptiSim, twsOptiSim, BFI)
19 use mo_kind, only : dp, i4
21 real(dp), dimension(:), intent(in) :: parameterset
22 integer(i4), dimension(:), optional, intent(in) :: opti_domain_indices
23 real(dp), dimension(:, :), allocatable, optional, intent(out) :: runoff !< dim1=time dim2=gauge
24 type(optidata_sim), dimension(:), optional, intent(inout) :: smOptiSim !< dim1=ncells, dim2=time
25 type(optidata_sim), dimension(:), optional, intent(inout) :: neutronsOptiSim !< dim1=ncells, dim2=time
26 type(optidata_sim), dimension(:), optional, intent(inout) :: etOptiSim !< dim1=ncells, dim2=time
27 type(optidata_sim), dimension(:), optional, intent(inout) :: twsOptiSim !< dim1=ncells, dim2=time
28 real(dp), dimension(:), allocatable, optional, intent(out) :: BFI !< baseflow index, dim1=domainID
29 end subroutine
30 end interface
31
32 !> \brief Interface for objective function.
33 interface
34 function objective_interface (parameterset, eval, arg1, arg2, arg3)
35 use mo_kind, only : dp
36 import eval_interface
37 real(dp), intent(in), dimension(:) :: parameterset !< parameter set
38 procedure(eval_interface), INTENT(IN), pointer :: eval !< evaluation routine
39 real(dp), optional, intent(in) :: arg1 !< optional argument 1
40 real(dp), optional, intent(out) :: arg2 !< optional argument 2
41 real(dp), optional, intent(out) :: arg3 !< optional argument 3
42
43 real(dp) :: objective_interface
44 end function objective_interface
45 end interface
46
47 !> \brief Optimizee for a eval-objective pair
48 type, extends(optimizee) :: mhm_optimizee
49 procedure(eval_interface), pointer, nopass :: eval_pointer => null() !< Pointer to the eval
50 procedure(objective_interface), pointer, nopass :: obj_pointer => null() !< Pointer to the objective
51 contains
52 procedure :: evaluate => evaluate_obj_eval
53 end type mhm_optimizee
54
55 contains
56
57 !> \brief Implementation of the evaluate procedure for a eval-objective pair
58 function evaluate_obj_eval(self, parameters, sigma, stddev_new, likeli_new) result(value)
59 class(mhm_optimizee), intent(inout) :: self
60 real(dp), dimension(:), intent(in) :: parameters
61 real(dp), intent(in), optional :: sigma
62 real(dp), intent(out), optional :: stddev_new
63 real(dp), intent(out), optional :: likeli_new
64 real(dp) :: value
65
66 ! Ensure the eval function pointer is set
67 if (.not. associated(self%eval_pointer)) then
68 call error_message("Eval function pointer is not set in mhm_optimizee!")
69 end if
70 ! Ensure the objective function pointer is set
71 if (.not. associated(self%obj_pointer)) then
72 call error_message("Objective function pointer is not set in mhm_optimizee!")
73 end if
74
75 ! Call the objective function pointer
76 value = self%obj_pointer(parameters, self%eval_pointer, sigma, stddev_new, likeli_new)
77 end function evaluate_obj_eval
78
79end module mo_optimization_utils
Interface for evaluation routine.
Type definitions for optimization routines.
Utility functions, such as interface definitions, for optimization routines.
real(dp) function evaluate_obj_eval(self, parameters, sigma, stddev_new, likeli_new)
Implementation of the evaluate procedure for a eval-objective pair.
type for simulated optional data
Optimizee for a eval-objective pair.