tlm_adjoint.optimization
Module Contents
- tlm_adjoint.optimization.minimize_scipy(forward, M0, *, manager=None, **kwargs)
Minimize a functional using
scipy.optimize.minimize()
.- Parameters:
- forwardcallable
Accepts one or more variable arguments, and returns a scalar variable defining the functional.
- M0variable or Sequence[variable, …]
The initial guess.
- manager
EquationManager
Used to create an internal manager via
EquationManager.new()
. manager() is used if not supplied.- kwargs
Passed to
scipy.optimize.minimize()
.
- Returns:
- Mvariable or Sequence[variable, …]
The result of the minimization
- minimize_return_value
The return value from
scipy.optimize.minimize()
.
Warning
Data is gathered onto the root process (process zero) so that the serial
scipy.optimize.minimize()
can be used.
- class tlm_adjoint.optimization.TAOSolver(forward, space, *, solver_parameters=None, H_0_action=None, M_inv_action=None, manager=None)
Functional minimization using TAO.
- Parameters:
- forwardcallable
Accepts one or more variable arguments, and returns a scalar variable defining the functional.
- spacespace or Sequence[space, …]
The control space.
- solver_parametersMapping
TAO solver parameters.
- H_0_actioncallable
Defines the initial Hessian inverse approximation. Accepts one or more variables as arguments, defining a direction, and returns a variable or a
Sequence
of variables defining the (conjugate of) the action of an initial Hessian inverse approximation on this direction. Arguments should not be modified.- M_inv_actioncallable
Defines a dual space norm. Accepts one or more variables as arguments, defining a direction, and returns a variable or a
Sequence
of variables defining the (conjugate of) the action of a Hermitian and positive definite matrix on this direction. Arguments should not be modified. H_0_action is used if not supplied.- manager
EquationManager
Used to create an internal manager via
EquationManager.new()
. manager() is used if not supplied.
- property tao
The
petsc4py.PETSc.TAO
used to solve the optimization problem.
- solve(M)
Solve the optimization problem.
- Parameters:
- Mvariable or Sequence[variable, …]
Defines the initial guess, and stores the result of the minimization.
- tlm_adjoint.optimization.minimize_tao(forward, M0, *args, **kwargs)
Minimize a functional using TAO.
- Parameters:
- forwardcallable
Accepts one or more variable arguments, and returns a scalar variable defining the functional.
- M0variable or Sequence[variable, …]
The initial guess.
- args, kwargs
Passed to the
TAOSolver
constructor.
- Returns:
- variable or Sequence[variable, …]
The result of the minimization.