Class ae108::cpppetsc::TAOSolver
template <class Policy>
ClassList > ae108 > cpppetsc > TAOSolver
Classes
| Type | Name | 
|---|---|
| struct | EqualityConstraints | 
Public Types
| Type | Name | 
|---|---|
| typedef std::function< void(const distributed< vector_type > &, distributed< vector_type > *)> | GradientFunctor | 
| typedef std::function< void(const distributed< vector_type > &, matrix_type *)> | HessianFunctor | 
| typedef std::function< void(const distributed< vector_type > &, real_type *)> | ObjectiveFunctor | 
| enum std::uint8_t | Type | 
| typedef Matrix< Policy > | matrix_type | 
| typedef typename Vector< Policy >::real_type | real_type | 
| typedef typename Vector< Policy >::size_type | size_type | 
| typedef typename Vector< Policy >::value_type | value_type | 
| typedef Vector< Policy > | vector_type | 
Public Static Attributes
| Type | Name | 
|---|---|
| constexpr value_type | no_lower_bound   = PETSC\_NINFINITY | 
| constexpr value_type | no_upper_bound   = PETSC\_INFINITY | 
Public Functions
| Type | Name | 
|---|---|
| TAOSolver (matrix_type bufferMatrix) Initializes the nonlinear solver with the problem dimension. | |
| TAOSolver (matrix_type bufferMatrix, const Type defaultType) Initializes the nonlinear solver with the problem dimension. | |
| void | setBounds (const distributed< vector_type > & lowerBounds, const distributed< vector_type > & upperBounds) Define upper and lower bounds for the variables. | 
| void | setConstraints (EqualityConstraints constraints) Set the equality constraints of the problem. | 
| void | setType (const Type type) Set the type of the optimizer. | 
| distributed< vector_type > | solve (ObjectiveFunctor objective, GradientFunctor gradient, HessianFunctor hessian, distributed< vector_type > guess) const Minimizes the objective function using optimization. | 
Public Types Documentation
typedef GradientFunctor
using ae108::cpppetsc::TAOSolver< Policy >::GradientFunctor =  std::function<void(const distributed<vector_type> &, distributed<vector_type> *)>;
Remark:
The first parameter is the input, the second parameter is the output.
typedef HessianFunctor
using ae108::cpppetsc::TAOSolver< Policy >::HessianFunctor =  std::function<void(const distributed<vector_type> &, matrix_type *)>;
Remark:
The first parameter is the input, the second parameter is the output.
typedef ObjectiveFunctor
using ae108::cpppetsc::TAOSolver< Policy >::ObjectiveFunctor =  std::function<void(const distributed<vector_type> &, real_type *)>;
Remark:
The first parameter is the input, the second parameter is the output.
enum Type
enum ae108::cpppetsc::TAOSolver::Type {
    lmvm,
    nls,
    ntr,
    ntl,
    cg,
    tron,
    owlqn,
    bmrm,
    blmvm,
    bqnls,
    bncg,
    bnls,
    bntr,
    bqnkls,
    bqnktr,
    bqnktl,
    bqpip,
    gpcg,
    nm,
    pounders,
    brgn,
    lcl,
    ssils,
    ssfls,
    asils,
    asfls,
    ipm,
    pdipm,
    shell,
    admm,
    almm
};
typedef matrix_type
using ae108::cpppetsc::TAOSolver< Policy >::matrix_type =  Matrix<Policy>;
typedef real_type
using ae108::cpppetsc::TAOSolver< Policy >::real_type =  typename Vector<Policy>::real_type;
typedef size_type
using ae108::cpppetsc::TAOSolver< Policy >::size_type =  typename Vector<Policy>::size_type;
typedef value_type
using ae108::cpppetsc::TAOSolver< Policy >::value_type =  typename Vector<Policy>::value_type;
typedef vector_type
using ae108::cpppetsc::TAOSolver< Policy >::vector_type =  Vector<Policy>;
Public Static Attributes Documentation
variable no_lower_bound
constexpr value_type ae108::cpppetsc::TAOSolver< Policy >::no_lower_bound;
variable no_upper_bound
constexpr value_type ae108::cpppetsc::TAOSolver< Policy >::no_upper_bound;
Public Functions Documentation
function TAOSolver [1/2]
Initializes the nonlinear solver with the problem dimension.
explicit ae108::cpppetsc::TAOSolver::TAOSolver (
    matrix_type bufferMatrix
) 
Parameters:
- Thematrix to use as an internal buffer.
function TAOSolver [2/2]
Initializes the nonlinear solver with the problem dimension.
explicit ae108::cpppetsc::TAOSolver::TAOSolver (
    matrix_type bufferMatrix,
    const Type defaultType
) 
Parameters:
- Thematrix to use as an internal buffer.
- defaultTypeThe default type of solver used. This type will not be used if e.g. a command line option specifying the type is present (see PETSc documentation).
function setBounds
Define upper and lower bounds for the variables.
void ae108::cpppetsc::TAOSolver::setBounds (
    const distributed < vector_type > & lowerBounds,
    const distributed < vector_type > & upperBounds
) 
Parameters:
- lowerBoundsA vector of constraints for every variable. Use no_lower_bound to set no constraints for a variable.
- upperBoundsA vector of constraints for every variable. Use no_upper_bound to set no constraints for a variable.
function setConstraints
Set the equality constraints of the problem.
void ae108::cpppetsc::TAOSolver::setConstraints (
    EqualityConstraints constraints
) 
Note:
Not all solvers consider these constraints.
function setType
Set the type of the optimizer.
void ae108::cpppetsc::TAOSolver::setType (
    const Type type
) 
Parameters:
- typeThe type of the optimizer. See the definition of Type and the PETSc documentation for more information about the possible choices.
function solve
Minimizes the objective function using optimization.
distributed < vector_type > ae108::cpppetsc::TAOSolver::solve (
    ObjectiveFunctor objective,
    GradientFunctor gradient,
    HessianFunctor hessian,
    distributed < vector_type > guess
) const
Parameters:
- objectiveThe function to minimize.
- gradientThe gradient of the objective.
- hessianThe hessian of the objective.
- guessThe starting point for the optimization.
Exception:
- TAOSolverDivergedException if solve did not converge.
Returns:
The minimizer.
The documentation for this class was generated from the following file cpppetsc/src/include/ae108/cpppetsc/TAOSolver.h