ROL
|
Defines the reduced time-dependent objective function interface for simulation-based optimization. More...
#include <ROL_ReducedDynamicObjective.hpp>
Public Member Functions | |
ReducedDynamicObjective (const Ptr< DynamicObjective< Real >> &obj, const Ptr< DynamicConstraint< Real >> &con, const Ptr< Vector< Real >> &u0, const Ptr< Vector< Real >> &zvec, const Ptr< Vector< Real >> &cvec, const std::vector< TimeStamp< Real >> &timeStamp, ROL::ParameterList &pl, const Ptr< std::ostream > &stream=nullPtr) | |
Ptr< Vector< Real > > | makeDynamicVector (const Vector< Real > &x) const |
void | update (const Vector< Real > &x, bool flag=true, int iter=-1) |
Update objective function. More... | |
Real | value (const Vector< Real > &x, Real &tol) |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol) |
Compute gradient. More... | |
void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply Hessian approximation to vector. More... | |
Public Member Functions inherited from ROL::Objective< Real > | |
virtual | ~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 | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply preconditioner to vector. More... | |
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... | |
virtual void | setParameter (const std::vector< Real > ¶m) |
Private Types | |
using | size_type = typename std::vector< Real >::size_type |
Private Member Functions | |
PartitionedVector< Real > & | partition (Vector< Real > &x) const |
const PartitionedVector< Real > & | partition (const Vector< Real > &x) const |
void | throwError (const int err, const std::string &sfunc, const std::string &func, const int line) const |
Real | solveState (const Vector< Real > &x) |
Real | updateSketch (const Vector< Real > &x, const Real tol) |
void | setTerminalCondition (Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | computeAdjointRHS (Vector< Real > &rhs, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | advanceAdjoint (Vector< Real > &l, Vector< Real > &rhs, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | solveAdjoint (const Vector< Real > &x) |
void | updateGradient (Vector< Real > &g, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | advanceStateSens (Vector< Real > &wnew, const Vector< Real > &v, const Vector< Real > &wold, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | setTerminalConditionHess (Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | computeOldStateHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false) |
void | computeOldMixedHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false) |
void | computeNewStateJacobian (Vector< Real > &Hv, const Vector< Real > &p, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false) |
void | computeNewStateHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false) |
void | computeNewMixedHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false) |
void | advanceAdjointSens (Vector< Real > &p, Vector< Real > &rhs, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | computeControlHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | addMixedHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
void | addAdjointSens (Vector< Real > &Hv, const Vector< Real > &p, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts) |
Private Attributes | |
const Ptr< DynamicObjective < Real > > | obj_ |
const Ptr< DynamicConstraint < Real > > | con_ |
const Ptr< Vector< Real > > | u0_ |
const std::vector< TimeStamp < Real > > | timeStamp_ |
const size_type | Nt_ |
const bool | useSketch_ |
size_type | rankState_ |
Ptr< Sketch< Real > > | stateSketch_ |
size_type | rankAdjoint_ |
Ptr< Sketch< Real > > | adjointSketch_ |
size_type | rankStateSens_ |
Ptr< Sketch< Real > > | stateSensSketch_ |
const bool | useInexact_ |
const size_type | updateFactor_ |
const size_type | maxRank_ |
const bool | syncHessRank_ |
std::vector< Ptr< Vector< Real > > > | uhist_ |
std::vector< Ptr< Vector< Real > > > | lhist_ |
std::vector< Ptr< Vector< Real > > > | whist_ |
std::vector< Ptr< Vector< Real > > > | phist_ |
Ptr< Vector< Real > > | cprimal_ |
Ptr< Vector< Real > > | crhs_ |
Ptr< Vector< Real > > | udual_ |
Ptr< Vector< Real > > | rhs_ |
Ptr< Vector< Real > > | zdual_ |
Real | val_ |
bool | isValueComputed_ |
bool | isStateComputed_ |
bool | isAdjointComputed_ |
bool | useHessian_ |
bool | useSymHess_ |
Ptr< std::ostream > | stream_ |
const bool | print_ |
const int | freq_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Objective< Real > | |
const std::vector< Real > | getParameter (void) const |
Defines the reduced time-dependent objective function interface for simulation-based optimization.
This objective function implements the implicitly defined objective function given by
\[ F(z) := \sum_{n=1}^{N_t} f_n(u_{n-1}(z),u_n(z),z_n) \]
where \(f_n:\mathcal{U}\times\mathcal{U}\times\mathcal{Z}\to\mathbb{R}\), and \(u_n\in\mathcal{U}\) solves the system of equations
\[ c_n(u_{n-1},u_n,z_n) = 0,\quad n=1,\ldots,N_t \]
with \(u_0\) provided.
Disclaimer: This is currently only set up for single step time integrators and piecewise-constant-in-time controls.
Definition at line 80 of file ROL_ReducedDynamicObjective.hpp.
|
private |
Definition at line 81 of file ROL_ReducedDynamicObjective.hpp.
|
inline |
Definition at line 160 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::cprimal_, ROL::ReducedDynamicObjective< Real >::crhs_, ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::phist_, ROL::ReducedDynamicObjective< Real >::rankAdjoint_, ROL::ReducedDynamicObjective< Real >::rankState_, ROL::ReducedDynamicObjective< Real >::rankStateSens_, ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::stateSensSketch_, ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::syncHessRank_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::udual_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::useHessian_, ROL::ReducedDynamicObjective< Real >::useSketch_, ROL::ReducedDynamicObjective< Real >::whist_, and ROL::ReducedDynamicObjective< Real >::zdual_.
|
inlineprivate |
Definition at line 127 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
inlineprivate |
Definition at line 131 of file ROL_ReducedDynamicObjective.hpp.
|
inlineprivate |
Definition at line 135 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
inline |
Definition at line 239 of file ROL_ReducedDynamicObjective.hpp.
References ROL::PartitionedVector< Real >::create(), and ROL::ReducedDynamicObjective< Real >::Nt_.
|
inlinevirtual |
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 243 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::freq_, ROL::ReducedDynamicObjective< Real >::isAdjointComputed_, ROL::ReducedDynamicObjective< Real >::isStateComputed_, ROL::ReducedDynamicObjective< Real >::isValueComputed_, ROL::Vector< Real >::print(), ROL::ReducedDynamicObjective< Real >::print_, ROL::ReducedDynamicObjective< Real >::stateSensSketch_, ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::useHessian_, ROL::ReducedDynamicObjective< Real >::useSketch_, and ROL::ReducedDynamicObjective< Real >::val_.
|
inlinevirtual |
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 272 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::cprimal_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isStateComputed_, ROL::ReducedDynamicObjective< Real >::isValueComputed_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::useSketch_, and ROL::ReducedDynamicObjective< Real >::val_.
|
inlinevirtual |
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 315 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::con_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isAdjointComputed_, ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::rankState_, ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::stream_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::updateGradient(), ROL::ReducedDynamicObjective< Real >::updateSketch(), ROL::ReducedDynamicObjective< Real >::useInexact_, and ROL::ReducedDynamicObjective< Real >::useSketch_.
|
inlinevirtual |
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 410 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::addAdjointSens(), ROL::ReducedDynamicObjective< Real >::addMixedHessLag(), ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::advanceAdjointSens(), ROL::ReducedDynamicObjective< Real >::advanceStateSens(), ROL::ReducedDynamicObjective< Real >::computeControlHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateJacobian(), ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag(), ROL::ReducedDynamicObjective< Real >::con_, ROL::PartitionedVector< Real >::get(), ROL::Objective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::phist_, ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSensSketch_, ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::useHessian_, ROL::ReducedDynamicObjective< Real >::useSketch_, ROL::ReducedDynamicObjective< Real >::useSymHess_, and ROL::ReducedDynamicObjective< Real >::whist_.
|
inlineprivate |
Definition at line 581 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::cprimal_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isStateComputed_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, and ROL::ReducedDynamicObjective< Real >::useSketch_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
inlineprivate |
Definition at line 617 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::cprimal_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isAdjointComputed_, ROL::ReducedDynamicObjective< Real >::isStateComputed_, ROL::ReducedDynamicObjective< Real >::maxRank_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::rankAdjoint_, ROL::ReducedDynamicObjective< Real >::rankState_, ROL::ReducedDynamicObjective< Real >::rankStateSens_, ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSensSketch_, ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::stream_, ROL::ReducedDynamicObjective< Real >::syncHessRank_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, and ROL::ReducedDynamicObjective< Real >::updateFactor_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient().
|
inlineprivate |
Definition at line 680 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), and ROL::ReducedDynamicObjective< Real >::solveAdjoint().
|
inlineprivate |
Definition at line 687 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), and ROL::ReducedDynamicObjective< Real >::solveAdjoint().
|
inlineprivate |
Definition at line 696 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), and ROL::ReducedDynamicObjective< Real >::solveAdjoint().
|
inlineprivate |
Definition at line 704 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::con_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isAdjointComputed_, ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, and ROL::ReducedDynamicObjective< Real >::useSketch_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 768 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, and ROL::ReducedDynamicObjective< Real >::zdual_.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient().
|
inlineprivate |
Definition at line 780 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::cprimal_, and ROL::ReducedDynamicObjective< Real >::crhs_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 794 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::rhs_, and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 817 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 840 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 858 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::Vector< Real >::scale(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 873 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 896 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::udual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 914 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 924 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, and ROL::ReducedDynamicObjective< Real >::zdual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 935 of file ROL_ReducedDynamicObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::ReducedDynamicObjective< Real >::con_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::zdual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 952 of file ROL_ReducedDynamicObjective.hpp.
References ROL::ReducedDynamicObjective< Real >::con_, ROL::Vector< Real >::plus(), and ROL::ReducedDynamicObjective< Real >::zdual_.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
private |
Definition at line 84 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::addMixedHessLag(), ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::computeControlHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag(), ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::updateGradient(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 85 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::addAdjointSens(), ROL::ReducedDynamicObjective< Real >::addMixedHessLag(), ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::advanceAdjointSens(), ROL::ReducedDynamicObjective< Real >::advanceStateSens(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::computeControlHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateJacobian(), ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag(), ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateGradient(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 86 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 87 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 88 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::makeDynamicVector(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 90 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::update(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 92 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 93 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::update(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 95 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 96 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::update(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 98 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
|
private |
Definition at line 101 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient().
|
private |
Definition at line 102 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 103 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 104 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 106 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::update(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
|
private |
Definition at line 108 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec(), and ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective().
|
private |
Definition at line 109 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec(), and ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective().
|
private |
Definition at line 110 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::advanceStateSens(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::updateSketch(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 111 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::advanceStateSens(), and ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective().
|
private |
Definition at line 112 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateJacobian(), ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), and ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess().
|
private |
Definition at line 113 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess(), and ROL::ReducedDynamicObjective< Real >::solveAdjoint().
|
private |
Definition at line 114 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::addAdjointSens(), ROL::ReducedDynamicObjective< Real >::addMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeControlHessLag(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::updateGradient().
|
private |
Definition at line 115 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::update(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
Definition at line 117 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::update(), and ROL::ReducedDynamicObjective< Real >::value().
|
private |
|
private |
|
private |
Definition at line 120 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective(), and ROL::ReducedDynamicObjective< Real >::update().
|
private |
Definition at line 121 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::hessVec().
|
private |
Definition at line 123 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::gradient(), and ROL::ReducedDynamicObjective< Real >::updateSketch().
|
private |
Definition at line 124 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::update().
|
private |
Definition at line 125 of file ROL_ReducedDynamicObjective.hpp.
Referenced by ROL::ReducedDynamicObjective< Real >::update().