ROL
|
Provides the interface to evaluate nonlinear least squares objective functions. More...
#include <ROL_NonlinearLeastSquaresObjective.hpp>
Public Member Functions | |
NonlinearLeastSquaresObjective (const Ptr< Constraint< Real > > &con, const Vector< Real > &optvec, const Vector< Real > &convec, const bool GNH=false) | |
Constructor. More... | |
void | update (const Vector< Real > &x, UpdateType type, int iter=-1) override |
Update objective function. More... | |
void | update (const Vector< Real > &x, bool flag=true, int iter=-1) override |
Update objective function. More... | |
Real | value (const Vector< Real > &x, Real &tol) override |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol) override |
Compute gradient. More... | |
void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply Hessian approximation to vector. More... | |
void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply preconditioner to vector. More... | |
void | setParameter (const std::vector< Real > ¶m) override |
Public Member Functions inherited from ROL::Objective< Real > | |
virtual | ~Objective () |
Objective () | |
virtual Real | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
Compute directional derivative. More... | |
virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply inverse Hessian approximation to vector. More... | |
virtual void | prox (Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol) |
virtual std::vector < std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference gradient check. More... | |
virtual std::vector < std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference gradient check. More... | |
virtual std::vector < std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference gradient check with specified step sizes. More... | |
virtual std::vector < std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference gradient check with specified step sizes. More... | |
virtual std::vector < std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference Hessian-applied-to-vector check. More... | |
virtual std::vector < std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference Hessian-applied-to-vector check. More... | |
virtual std::vector < std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference Hessian-applied-to-vector check with specified step sizes. More... | |
virtual std::vector < std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference Hessian-applied-to-vector check with specified step sizes. More... | |
virtual std::vector< Real > | checkHessSym (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout) |
Hessian symmetry check. More... | |
virtual std::vector< Real > | checkHessSym (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout) |
Hessian symmetry check. More... | |
Private Attributes | |
const Ptr< Constraint< Real > > | con_ |
const bool | GaussNewtonHessian_ |
Ptr< Vector< Real > > | c1_ |
Ptr< Vector< Real > > | c2_ |
Ptr< Vector< Real > > | c1dual_ |
Ptr< Vector< Real > > | x_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Objective< Real > | |
const std::vector< Real > | getParameter (void) const |
Provides the interface to evaluate nonlinear least squares objective functions.
ROL's nonlinear least squares objective function interface constructs the the nonlinear least squares objective function associated with the equality constraint \(c(x)=0\). That is,
\[ J(x) = \langle \mathfrak{R} c(x),c(x) \rangle_{\mathcal{C}^*,\mathcal{C}} \]
where \(c:\mathcal{X}\to\mathcal{C}\) and \(\mathfrak{R}\in\mathcal{L}( \mathcal{C},\mathcal{C}^*)\) denotes the Riesz map from \(\mathcal{C}\) into \(\mathcal{C}^*\).
Definition at line 39 of file ROL_NonlinearLeastSquaresObjective.hpp.
ROL::NonlinearLeastSquaresObjective< Real >::NonlinearLeastSquaresObjective | ( | const Ptr< Constraint< Real > > & | con, |
const Vector< Real > & | optvec, | ||
const Vector< Real > & | convec, | ||
const bool | GNH = false |
||
) |
Constructor.
This function constructs a nonlinear least squares objective function.
[in] | con | is the nonlinear equation to be solved. |
[in] | vec | is a constraint space vector used for cloning. |
[in] | GHN | is a flag dictating whether or not to use the Gauss-Newton Hessian. |
Definition at line 16 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
References ROL::NonlinearLeastSquaresObjective< Real >::c1_, ROL::NonlinearLeastSquaresObjective< Real >::c1dual_, ROL::NonlinearLeastSquaresObjective< Real >::c2_, ROL::Vector< Real >::clone(), ROL::Vector< Real >::dual(), and ROL::NonlinearLeastSquaresObjective< Real >::x_.
|
overridevirtual |
Update objective function.
This function updates the objective function at new iterations.
[in] | x | is the new iterate. |
[in] | type | is the type of update requested. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Objective< Real >.
Definition at line 27 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
|
overridevirtual |
Update objective function.
This function updates the objective function at new iterations.
[in] | x | is the new iterate. |
[in] | flag | is true if the iterate has changed. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Objective< Real >.
Definition at line 35 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
|
overridevirtual |
Compute value.
This function returns the objective function value.
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Implements ROL::Objective< Real >.
Definition at line 43 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
|
overridevirtual |
Compute gradient.
This function returns the objective function gradient.
[out] | g | is the gradient. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
The default implementation is a finite-difference approximation based on the function value. This requires the definition of a basis \(\{\phi_i\}\) for the optimization vectors x and the definition of a basis \(\{\psi_j\}\) for the dual optimization vectors (gradient vectors g). The bases must be related through the Riesz map, i.e., \( R \{\phi_i\} = \{\psi_j\}\), and this must be reflected in the implementation of the ROL::Vector::dual() method.
Reimplemented from ROL::Objective< Real >.
Definition at line 49 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
|
overridevirtual |
Apply Hessian approximation to vector.
This function applies the Hessian of the objective function to the vector \(v\).
[out] | hv | is the the action of the Hessian on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::Objective< Real >.
Definition at line 54 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
References ROL::Vector< Real >::dual(), and ROL::Vector< Real >::plus().
|
overridevirtual |
Apply preconditioner to vector.
This function applies a preconditioner for the Hessian of the objective function to the vector \(v\).
[out] | Pv | is the action of the Hessian preconditioner on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::Objective< Real >.
Definition at line 64 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
References ROL::Vector< Real >::dual().
|
overridevirtual |
Reimplemented from ROL::Objective< Real >.
Definition at line 69 of file ROL_NonlinearLeastSquaresObjective_Def.hpp.
References ROL::Objective< Real >::setParameter().
|
private |
Definition at line 41 of file ROL_NonlinearLeastSquaresObjective.hpp.
|
private |
Definition at line 42 of file ROL_NonlinearLeastSquaresObjective.hpp.
|
private |
Definition at line 44 of file ROL_NonlinearLeastSquaresObjective.hpp.
Referenced by ROL::NonlinearLeastSquaresObjective< Real >::NonlinearLeastSquaresObjective().
|
private |
Definition at line 44 of file ROL_NonlinearLeastSquaresObjective.hpp.
Referenced by ROL::NonlinearLeastSquaresObjective< Real >::NonlinearLeastSquaresObjective().
|
private |
Definition at line 44 of file ROL_NonlinearLeastSquaresObjective.hpp.
Referenced by ROL::NonlinearLeastSquaresObjective< Real >::NonlinearLeastSquaresObjective().
|
private |
Definition at line 44 of file ROL_NonlinearLeastSquaresObjective.hpp.
Referenced by ROL::NonlinearLeastSquaresObjective< Real >::NonlinearLeastSquaresObjective().