cashocs.OptimalControlProblem#
- class cashocs.OptimalControlProblem(state_forms: list[ufl.Form] | ufl.Form, bcs_list: list[list[fenics.DirichletBC]] | list[fenics.DirichletBC] | fenics.DirichletBC, cost_functional_form: list[_typing.CostFunctional] | _typing.CostFunctional, states: list[fenics.Function] | fenics.Function, controls: list[fenics.Function] | fenics.Function, adjoints: list[fenics.Function] | fenics.Function, config: io.Config | None = None, riesz_scalar_products: list[ufl.Form] | ufl.Form | None = None, control_constraints: list[list[float | fenics.Function]] | None = None, initial_guess: list[fenics.Function] | None = None, ksp_options: _typing.KspOption | list[_typing.KspOption] | None = None, adjoint_ksp_options: _typing.KspOption | list[_typing.KspOption] | None = None, gradient_ksp_options: _typing.KspOption | list[_typing.KspOption] | None = None, desired_weights: list[float] | None = None, control_bcs_list: list[list[fenics.DirichletBC]] | list[fenics.DirichletBC] | fenics.DirichletBC | None = None, preconditioner_forms: list[ufl.Form] | ufl.Form | None = None, pre_callback: Callable | None = None, post_callback: Callable | None = None, linear_solver: _utils.linalg.LinearSolver | None = None, adjoint_linear_solver: _utils.linalg.LinearSolver | None = None, newton_linearizations: ufl.Form | list[ufl.Form] | None = None, excluded_from_time_derivative: list[int] | list[list[int]] | list[None] | None = None)[source]#
Bases:
OptimizationProblem
Implements an optimal control problem.
This class is used to define an optimal control problem, and also to solve it subsequently. For a detailed documentation, see the examples in the tutorial. For easier input, when considering single (state or control) variables, these do not have to be wrapped into a list. Note, that in the case of multiple variables these have to be grouped into ordered lists, where state_forms, bcs_list, states, adjoints have to have the same order (i.e.
[y1, y2]
and[p1, p2]
, wherep1
is the adjoint ofy1
and so on).Initializes self.
- Parameters:
state_forms (list[ufl.Form] | ufl.Form) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.
bcs_list (list[list[fenics.DirichletBC]] | list[fenics.DirichletBC] | fenics.DirichletBC) – The list of
fenics.DirichletBC
objects describing Dirichlet (essential) boundary conditions. If this isNone
, then no Dirichlet boundary conditions are imposed.cost_functional_form (list[_typing.CostFunctional] | _typing.CostFunctional) – UFL form of the cost functional. Can also be a list of summands of the cost functional
states (list[fenics.Function] | fenics.Function) – The state variable(s), can either be a
fenics.Function
, or a list of these.controls (list[fenics.Function] | fenics.Function) – The control variable(s), can either be a
fenics.Function
, or a list of these.adjoints (list[fenics.Function] | fenics.Function) – The adjoint variable(s), can either be a
fenics.Function
, or a (ordered) list of these.config (io.Config | None) – The config file for the problem, generated via
cashocs.load_config()
. Alternatively, this can also beNone
, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in thesolve
method. The default isNone
.riesz_scalar_products (list[ufl.Form] | ufl.Form | None) – The scalar products of the control space. Can either be None, a single UFL form, or a (ordered) list of UFL forms. If
None
, the \(L^2(\Omega)\) product is used (default isNone
).control_constraints (list[list[float | fenics.Function]] | None) – Box constraints posed on the control,
None
means that there are none (default isNone
). The (inner) lists should contain two elements of the form[u_a, u_b]
, whereu_a
is the lower, andu_b
the upper bound.initial_guess (list[fenics.Function] | None) – list of functions that act as initial guess for the state variables, should be valid input for
fenics.assign()
. Defaults toNone
, which means a zero initial guess.ksp_options (_typing.KspOption | list[_typing.KspOption] | None) – A list of dicts corresponding to command line options for PETSc, used to solve the state systems. If this is
None
, then the direct solver mumps is used (default isNone
).adjoint_ksp_options (_typing.KspOption | list[_typing.KspOption] | None) – A list of dicts corresponding to command line options for PETSc, used to solve the adjoint systems. If this is
None
, then the same options as for the state systems are used (default isNone
).gradient_ksp_options (_typing.KspOption | list[_typing.KspOption] | None) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is
None
, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).desired_weights (list[float] | None) – A list of values for scaling the cost functional terms. If this is supplied, the cost functional has to be given as list of summands. The individual terms are then scaled, so that term i has the magnitude of desired_weights[i] for the initial iteration. In case that desired_weights is None, no scaling is performed. Default is None.
control_bcs_list (list[list[fenics.DirichletBC]] | list[fenics.DirichletBC] | fenics.DirichletBC | None) – A list of boundary conditions for the control variables. This is passed analogously to
bcs_list
. Default isNone
.preconditioner_forms (list[ufl.Form] | ufl.Form | None) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.
pre_callback (Callable | None) – A function (without arguments) that will be called before each solve of the state system
post_callback (Callable | None) – A function (without arguments) that will be called after the computation of the gradient.
linear_solver (_utils.linalg.LinearSolver | None) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.
adjoint_linear_solver (_utils.linalg.LinearSolver | None) – The linear solver (KSP) which is used to solve the (linear) adjoint system.
newton_linearizations (ufl.Form | list[ufl.Form] | None) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.
excluded_from_time_derivative (list[int] | list[list[int]] | list[None] | None) – For each state equation, a list of indices which are not part of the first order time derivative for pseudo time stepping. Example: Pressure for incompressible flow. Default is None.
Examples
Examples how to use this class can be found in the tutorial.
Methods Summary
Solves the adjoint system.
Solves the Riesz problem to determine the gradient.
Solves the state system.
gradient_test
([u, h, rng])Performs a Taylor test to verify correctness of the computed gradient.
initialize_solve_parameters
([algorithm, ...])Solves the optimization problem by the method specified in the config file.
inject_post_callback
(function)Changes the a-posteriori callback of the OptimizationProblem to function.
inject_pre_callback
(function)Changes the a-priori callback of the OptimizationProblem to function.
inject_pre_post_callback
(pre_function, ...)Changes the a-priori (pre) and a-posteriori (post) callbacks.
shift_cost_functional
([shift])Shifts the cost functional by a constant.
solve
([algorithm, rtol, atol, max_iter])Solves the optimization problem by the method specified in the config file.
supply_adjoint_forms
(adjoint_forms, ...)Overwrites the computed weak forms of the adjoint system.
supply_custom_forms
(derivatives, ...)Overrides both adjoint system and derivatives with user input.
supply_derivatives
(derivatives)Overwrites the derivatives of the reduced cost functional w.r.t.
Attributes Summary
Methods Documentation
- compute_adjoint_variables() None #
Solves the adjoint system.
This can be used for debugging purposes and solver validation. Updates / overwrites the user input for the adjoint variables. The solution of the corresponding state system needed to determine the adjoints is carried out automatically.
- Return type:
None
- compute_gradient() list[fenics.Function] [source]#
Solves the Riesz problem to determine the gradient.
This can be used for debugging, or code validation. The necessary solutions of the state and adjoint systems are carried out automatically.
- Returns:
A list consisting of the (components) of the gradient.
- Return type:
list[fenics.Function]
- compute_state_variables() None #
Solves the state system.
This can be used for debugging purposes and to validate the solver. Updates and overwrites the user input for the state variables.
- Return type:
None
- gradient_test(u: list[fenics.Function] | None = None, h: list[fenics.Function] | None = None, rng: RandomState | None = None) float [source]#
Performs a Taylor test to verify correctness of the computed gradient.
- Parameters:
u (list[fenics.Function] | None) – The point, at which the gradient shall be verified. If this is
None
, then the current controls of the optimization problem are used. Default isNone
.h (list[fenics.Function] | None) – The direction(s) for the directional (Gâteaux) derivative. If this is
None
, one random direction is chosen. Default isNone
.rng (RandomState | None) – A numpy random state for calculating a random direction.
- Returns:
The convergence order from the Taylor test. If this is (approximately) 2 or larger, everything works as expected.
- Return type:
float
- initialize_solve_parameters(algorithm: str | None = None, rtol: float | None = None, atol: float | None = None, max_iter: int | None = None) None #
Solves the optimization problem by the method specified in the config file.
- Parameters:
algorithm (str | None) – Selects the optimization algorithm. Valid choices are
'gradient_descent'
or'gd'
for a gradient descent method,'conjugate_gradient'
,'nonlinear_cg'
,'ncg'
or'cg'
for nonlinear conjugate gradient methods, and'lbfgs'
or'bfgs'
for limited memory BFGS methods. This overwrites the value specified in the config file. If this isNone
, then the value in the config file is used. Default isNone
. In addition, for optimal control problems, one can use'newton'
for a truncated Newton method.rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.
- Return type:
None
Notes
If either
rtol
oratol
are specified as arguments to the.solve
call, the termination criterion changes to:a purely relative one (if only
rtol
is specified), i.e.,
\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]a purely absolute one (if only
atol
is specified), i.e.,
\[|| \nabla J(u_K) || \leq \texttt{atol}.\]a combined one if both
rtol
andatol
are specified, i.e.,
\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
- inject_post_callback(function: Callable | None) None #
Changes the a-posteriori callback of the OptimizationProblem to function.
- Parameters:
function (Callable | None) – A custom function without arguments, which will be called after the computation of the gradient(s)
- Return type:
None
- inject_pre_callback(function: Callable | None) None #
Changes the a-priori callback of the OptimizationProblem to function.
- Parameters:
function (Callable | None) – A custom function without arguments, which will be called before each solve of the state system
- Return type:
None
- inject_pre_post_callback(pre_function: Callable | None, post_function: Callable | None) None #
Changes the a-priori (pre) and a-posteriori (post) callbacks.
- Parameters:
pre_function (Callable | None) – A function without arguments, which is to be called before each solve of the state system
post_function (Callable | None) – A function without arguments, which is to be called after each computation of the (shape) gradient
- Return type:
None
- shift_cost_functional(shift: float = 0.0) None #
Shifts the cost functional by a constant.
- Parameters:
shift (float) – The constant, by which the cost functional is shifted.
- Return type:
None
- solve(algorithm: str | None = None, rtol: float | None = None, atol: float | None = None, max_iter: int | None = None) None [source]#
Solves the optimization problem by the method specified in the config file.
Updates / overwrites states, controls, and adjoints according to the optimization method, i.e., the user-input
fenics.Function()
objects.- Parameters:
algorithm (str | None) – Selects the optimization algorithm. Valid choices are
'gradient_descent'
or'gd'
for a gradient descent method,'conjugate_gradient'
,'nonlinear_cg'
,'ncg'
or'cg'
for nonlinear conjugate gradient methods,'lbfgs'
or'bfgs'
for limited memory BFGS methods, and'newton'
for a truncated Newton method. This overwrites the value specified in the config file. If this isNone
, then the value in the config file is used. Default isNone
.rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is
None
, the value from the config file is taken. Default isNone
.
- Return type:
None
Notes
If either
rtol
oratol
are specified as arguments to the.solve
call, the termination criterion changes to: - a purely relative one (if onlyrtol
is specified), i.e.,\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]a purely absolute one (if only
atol
is specified), i.e.,
\[|| \nabla J(u_K) || \leq \texttt{atol}.\]a combined one if both
rtol
andatol
are specified, i.e.,
\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
- supply_adjoint_forms(adjoint_forms: ufl.Form | list[ufl.Form], adjoint_bcs_list: fenics.DirichletBC | list[fenics.DirichletBC] | list[list[fenics.DirichletBC]]) None #
Overwrites the computed weak forms of the adjoint system.
This allows the user to specify their own weak forms of the problems and to use cashocs merely as a solver for solving the optimization problems.
- Parameters:
adjoint_forms (ufl.Form | list[ufl.Form]) – The UFL forms of the adjoint system(s).
adjoint_bcs_list (fenics.DirichletBC | list[fenics.DirichletBC] | list[list[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).
- Return type:
None
- supply_custom_forms(derivatives: ufl.Form | list[ufl.Form], adjoint_forms: ufl.Form | list[ufl.Form], adjoint_bcs_list: fenics.DirichletBC | list[fenics.DirichletBC] | list[list[fenics.DirichletBC]]) None [source]#
Overrides both adjoint system and derivatives with user input.
This allows the user to specify both the derivatives of the reduced cost functional and the corresponding adjoint system, and allows them to use cashocs as a solver.
- Parameters:
derivatives (ufl.Form | list[ufl.Form]) – The derivatives of the reduced (!) cost functional w.r.t. the control variables.
adjoint_forms (ufl.Form | list[ufl.Form]) – The UFL forms of the adjoint system(s).
adjoint_bcs_list (fenics.DirichletBC | list[fenics.DirichletBC] | list[list[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).
- Return type:
None
- supply_derivatives(derivatives: ufl.Form | list[ufl.Form]) None [source]#
Overwrites the derivatives of the reduced cost functional w.r.t. controls.
This allows users to implement their own derivatives and use cashocs as a solver library only.
- Parameters:
derivatives (ufl.Form | list[ufl.Form]) – The derivatives of the reduced (!) cost functional w.r.t.
variables. (the control)
- Return type:
None
Attributes Documentation
- uses_custom_scalar_product: bool = False#
- gradient: list[fenics.Function]#
- reduced_cost_functional: cost_functional.ReducedCostFunctional#
- gradient_problem: _pde_problems.ControlGradientProblem#
- output_manager: io.OutputManager#
- form_handler: _forms.ControlFormHandler#
- optimization_variable_abstractions: ova.OptimizationVariableAbstractions#
- adjoint_problem: _pde_problems.AdjointProblem#
- state_problem: _pde_problems.StateProblem#
- line_search: ls.LineSearch#
- hessian_problem: _pde_problems.HessianProblem#
- solver: optimization_algorithms.OptimizationAlgorithm | topology_optimization_algorithm.TopologyOptimizationAlgorithm#
- initial_guess: list[fenics.Function] | None#
- cost_functional_list: list[_typing.CostFunctional]#