Add optimizee type
Overview
This MR adds an abstract optimizee
type that serves as input for the optimizers in FORCES (anneal, mcmc, dds, sce).
This type has only a single deferred type-bound procedure evaluate
which has the signature of a likelihood function:
type, abstract :: optimizee
contains
procedure(evaluate_interface), deferred :: evaluate
end type optimizee
abstract interface
function evaluate_interface(self, parameters, sigma, stddev_new, likeli_new) result(value)
use mo_kind, only : dp
import :: optimizee
class(optimizee), intent(inout) :: self
real(dp), intent(in), dimension(:) :: parameters !< parameter set
real(dp), intent(in), optional :: sigma !< likelihood: standard deviation of data
real(dp), intent(out), optional :: stddev_new !< likelihood: standard deviation of errors using parameters
real(dp), intent(out), optional :: likeli_new !< likelihood: likelihood using stddev_new
real(dp) :: value !< output value
end function evaluate_interface
end interface
All likelihood related arguments are optional so that in case of a simple cost function optimization only the parameter set is required.
For convenience, three implementations of an optimizee are provided with FORCES:
-
function_optimizee
: Holding a procedure pointer which can be set to a simple cost function taking only a 1D parameter vector:f(x)
-
likelihood_optimizee
: Holding a procedure pointer which can be set to a likelihood functionf(x, sigma, stddev_new, likeli_new)
-
eval_optimizee
: Holding two pointers: (i) for an eval function and (ii) for an objective function. This is implemented to recreate to previous behavior of the optimization implementations
Usage
In order to optimize a simple cost function, that takes a parameter array and returns a single real value, one can use the provided function_optimizee
type with any optimizer. Here an example for SCE optimizing the griewank
test function:
use mo_kind, only: dp
use mo_opt_functions, only: griewank
use mo_optimization_utils, only: function_optimizee
use mo_sce, only: sce
real(dp), dimension(10) :: pini
real(dp), dimension(10,2) :: prange
type(function_optimizee) :: objective
prange(:,1) = (/ -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0, -600.0 /)
prange(:,2) = (/ 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0, 600.0 /)
pini = (/ -.226265e+01, -.130187e+01, -.151219e+01, 0.133983e+00, 0.988159e+00, &
-.495074e+01, -.126574e+02, 0.572684e+00, 0.303864e+01, 0.343031e+01 /)
objective%func_pointer => griewank
opt = sce(objective, pini, prange)
In order to create your own optimizee
you can extend the provided abstract optimizee
type:
use mo_optimization_utils, only: optimizee
type, extends(optimizee) :: your_optimizee
contains
! to be implemented with signature: your_evaluate(self, parameters, sigma, stddev_new, likeli_new)
procedure :: evaluate => your_evaluate
end type function_optimizee
You need to provide an implementation for the evaluate
procedure (your_evaluate
in the above example) justifying the evaluate_interface
from above. You can add other attributes or procedures as you like.
Backward compatibility
Before this MR, an optimizer would need an eval function and an objective function where internally, the objective function called the given eval function. In case of SCE this looked like this using the tailor made griewank_objective
test function which takes a dummy eval:
use mo_optimization_utils, only: eval_interface, objective_interface
use mo_opt_functions, only: griewank_objective, eval_dummy
! ...
procedure(eval_interface), pointer :: eval
procedure(objective_interface), pointer :: obj_func
! ...
eval => eval_dummy
obj_func => griewank_objective
opt = sce(eval, obj_func, pini, prange)
To recreate this behavior, you now can use the provided eval_optimizee
:
use mo_optimization_utils, only: eval_optimizee
use mo_opt_functions, only: griewank_objective, eval_dummy
! ...
type(eval_optimizee) :: objective
! ...
eval_optimizee%eval_pointer => eval_dummy
eval_optimizee%obj_pointer => griewank_objective
opt = sce(eval_optimizee, pini, prange)