0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_optimization_utils.f90
Go to the documentation of this file.
1!> \file mo_optimization_utils.f90
2!> \brief \copybrief mo_optimization_utils
3!> \details \copydetails mo_optimization_utils
4
5!> \brief Utility functions, such as interface definitions, for optimization routines.
6!> \copyright Copyright 2005-\today, the CHS Developers, Sabine Attinger: All rights reserved.
7!! FORCES is released under the LGPLv3+ license \license_note
9
10 use mo_kind, only : dp
11
12 implicit none
13
14 !> \brief Interface for evaluation routine.
15 abstract interface
16 subroutine eval_interface(parameterset, opti_domain_indices, runoff, smOptiSim, neutronsOptiSim, etOptiSim, twsOptiSim, &
17 lake_level, lake_volume, lake_area, lake_spill, lake_outflow, BFI)
18 use mo_kind, only : dp, i4
20 real(dp), dimension(:), intent(in) :: parameterset
21 integer(i4), dimension(:), optional, intent(in) :: opti_domain_indices
22 real(dp), dimension(:, :), allocatable, optional, intent(out) :: runoff !< dim1=time dim2=gauge
23 type(optidata_sim), dimension(:), optional, intent(inout) :: smOptiSim !< dim1=ncells, dim2=time
24 type(optidata_sim), dimension(:), optional, intent(inout) :: neutronsOptiSim !< dim1=ncells, dim2=time
25 type(optidata_sim), dimension(:), optional, intent(inout) :: etOptiSim !< dim1=ncells, dim2=time
26 type(optidata_sim), dimension(:), optional, intent(inout) :: twsOptiSim !< dim1=ncells, dim2=time
27 real(dp), dimension(:, :), allocatable, optional, intent(out) :: lake_level !< dim1=time dim2=lake
28 real(dp), dimension(:, :), allocatable, optional, intent(out) :: lake_volume !< dim1=time dim2=lake
29 real(dp), dimension(:, :), allocatable, optional, intent(out) :: lake_area !< dim1=time dim2=lake
30 real(dp), dimension(:, :), allocatable, optional, intent(out) :: lake_spill !< dim1=time dim2=lake
31 real(dp), dimension(:, :), allocatable, optional, intent(out) :: lake_outflow !< dim1=time dim2=lake
32 real(dp), dimension(:), allocatable, optional, intent(out) :: BFI !< baseflow index, dim1=domainID
33 end subroutine
34 end interface
35
36 !> \brief Interface for objective function.
37 interface
38 function objective_interface (parameterset, eval, arg1, arg2, arg3)
39 use mo_kind, only : dp
40 import eval_interface
41 real(dp), intent(in), dimension(:) :: parameterset !< parameter set
42 procedure(eval_interface), INTENT(IN), pointer :: eval !< evaluation routine
43 real(dp), optional, intent(in) :: arg1 !< optional argument 1
44 real(dp), optional, intent(out) :: arg2 !< optional argument 2
45 real(dp), optional, intent(out) :: arg3 !< optional argument 3
46
47 real(dp) :: objective_interface
48 end function objective_interface
49 end interface
50
51end module mo_optimization_utils
Interface for evaluation routine.
Define number representations.
Definition mo_kind.F90:17
integer, parameter i4
4 Byte Integer Kind
Definition mo_kind.F90:40
integer, parameter dp
Double Precision Real Kind.
Definition mo_kind.F90:46
Type definitions for optimization routines.
Utility functions, such as interface definitions, for optimization routines.
type for simulated optional data