ROL
ROL_NonlinearLeastSquaresObjective_Dynamic.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_NONLINEARLEASTSQUARESOBJECTIVE_DYNAMIC_H
11#define ROL_NONLINEARLEASTSQUARESOBJECTIVE_DYNAMIC_H
12
13#include "ROL_Objective.hpp"
15#include "ROL_Types.hpp"
16
36namespace ROL {
37
38template <class Real>
39class DynamicConstraint;
40
41template <class Real>
43private:
44 const Ptr<DynamicConstraint<Real>> con_;
45 const Ptr<const Vector<Real>> uo_;
46 const Ptr<const Vector<Real>> z_;
47 const Ptr<const TimeStamp<Real>> ts_;
49
50 Ptr<Vector<Real> > c1_, c2_, cdual_, udual_;
51
52public:
61 const Vector<Real> &c,
62 const Ptr<const Vector<Real>> &uo,
63 const Ptr<const Vector<Real>> &z,
64 const Ptr<const TimeStamp<Real>> &ts,
65 const bool GNH = false)
66 : con_(con), uo_(uo), z_(z), ts_(ts), GaussNewtonHessian_(GNH) {
67 c1_ = c.clone();
68 c2_ = c.clone();
69 cdual_ = c.dual().clone();
70 udual_ = uo->dual().clone();
71 }
72
73 void update( const Vector<Real> &u, bool flag = true, int iter = -1 ) {
74 //con_->update_un(u,*ts_);
75 con_->update(*uo_,u,*z_,*ts_);
76 con_->value(*c1_,*uo_,u,*z_,*ts_);
77 cdual_->set(c1_->dual());
78 }
79
80 Real value( const Vector<Real> &x, Real &tol ) {
81 Real half(0.5);
82 return half*(c1_->dot(*cdual_));
83 }
84
85 void gradient( Vector<Real> &g, const Vector<Real> &u, Real &tol ) {
86 con_->applyAdjointJacobian_un(g,*cdual_,*uo_,u,*z_,*ts_);
87 }
88
89 void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &u, Real &tol ) {
90 con_->applyJacobian_un(*c2_,v,*uo_,u,*z_,*ts_);
91 con_->applyAdjointJacobian_un(hv,c2_->dual(),*uo_,u,*z_,*ts_);
92 if ( !GaussNewtonHessian_ ) {
93 con_->applyAdjointHessian_un_un(*udual_,*cdual_,v,*uo_,u,*z_,*ts_);
94 hv.plus(*udual_);
95 }
96 }
97
98 void precond( Vector<Real> &pv, const Vector<Real> &v, const Vector<Real> &u, Real &tol ) {
99 con_->applyInverseAdjointJacobian_un(*cdual_,v,*uo_,u,*z_,*ts_);
100 con_->applyInverseJacobian_un(pv,cdual_->dual(),*uo_,u,*z_,*ts_);
101 }
102
103// Definitions for parametrized (stochastic) equality constraints
104//public:
105// void setParameter(const std::vector<Real> &param) {
106// Objective<Real>::setParameter(param);
107// con_->setParameter(param);
108// }
109};
110
111} // namespace ROL
112
113#endif
Contains definitions of custom data types in ROL.
Defines the time-dependent constraint operator interface for simulation-based optimization.
Provides the interface to evaluate nonlinear least squares objective functions.
void update(const Vector< Real > &u, bool flag=true, int iter=-1)
Update objective function.
void precond(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &u, Real &tol)
Apply preconditioner to vector.
NonlinearLeastSquaresObjective_Dynamic(const Ptr< DynamicConstraint< Real > > &con, const Vector< Real > &c, const Ptr< const Vector< Real > > &uo, const Ptr< const Vector< Real > > &z, const Ptr< const TimeStamp< Real > > &ts, const bool GNH=false)
Constructor.
void gradient(Vector< Real > &g, const Vector< Real > &u, Real &tol)
Compute gradient.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, Real &tol)
Apply Hessian approximation to vector.
Provides the interface to evaluate objective functions.
Defines the linear algebra or vector space interface.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Contains local time step information.