ROL
ROL_StdConstraint_Def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_STDEQUALITY_CONSTRAINT_DEF_H
45 #define ROL_STDEQUALITY_CONSTRAINT_DEF_H
46 
47 namespace ROL {
48 
49 template<typename Real>
50 void StdConstraint<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
51  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
52  update(*(xs.getVector()),flag,iter);
53 }
54 
55 template<typename Real>
56 void StdConstraint<Real>::update( const Vector<Real> &x, UpdateType type, int iter ) {
57  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
58  update(*(xs.getVector()),type,iter);
59 }
60 
61 template<typename Real>
62 void StdConstraint<Real>::value(Vector<Real> &c, const Vector<Real> &x, Real &tol) {
63  StdVector<Real> cs = dynamic_cast<StdVector<Real>&>(c);
64  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
65  value(*(cs.getVector()),*(xs.getVector()),tol);
66 }
67 
68 template<typename Real>
70  const Vector<Real> &v,
71  const Vector<Real> &x, Real &tol) {
72  StdVector<Real> jvs = dynamic_cast<StdVector<Real>&>(jv);
73  const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
74  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
75  try {
76  applyJacobian(*(jvs.getVector()),*(vs.getVector()),*(xs.getVector()),tol);
77  }
78  catch (std::exception &e ){
80  }
81 }
82 
83 template<typename Real>
84 void StdConstraint<Real>::applyJacobian( std::vector<Real> &jv,
85  const std::vector<Real> &v,
86  const std::vector<Real> &x, Real &tol ) {
87  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
88  ">>> ERROR (ROL::StdConstraint): applyJacobian not implemented!");
89 }
90 
91 template<typename Real>
93  const Vector<Real> &v,
94  const Vector<Real> &x, Real &tol) {
95  StdVector<Real> ajvs = dynamic_cast<StdVector<Real>&>(ajv);
96  const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
97  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
98  try {
99  applyAdjointJacobian(*(ajvs.getVector()),*(vs.getVector()),*(xs.getVector()),tol);
100  }
101  catch (std::exception &e ){
103  }
104 }
105 
106 template<typename Real>
107 void StdConstraint<Real>::applyAdjointJacobian( std::vector<Real> &ajv,
108  const std::vector<Real> &v,
109  const std::vector<Real> &x, Real &tol ) {
110  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
111  ">>> ERROR (ROL::StdConstraint): applyAdjointJacobian not implemented!");
112 }
113 
114 template<typename Real>
116  const Vector<Real> &u,
117  const Vector<Real> &v,
118  const Vector<Real> &x, Real &tol) {
119  StdVector<Real> ahuvs = dynamic_cast<StdVector<Real>&>(ahuv);
120  const StdVector<Real> us = dynamic_cast<const StdVector<Real>&>(u);
121  const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
122  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
123  try {
124  applyAdjointHessian( *(ahuvs.getVector()), *(us.getVector()), *(vs.getVector()),
125  *(xs.getVector()), tol );
126  }
127  catch (std::exception &e) {
129  }
130 }
131 
132 template<typename Real>
133 void StdConstraint<Real>::applyAdjointHessian( std::vector<Real> &ahuv,
134  const std::vector<Real> &u,
135  const std::vector<Real> &v,
136  const std::vector<Real> &x, Real &tol ) {
137  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
138  ">>> ERROR (ROL::StdConstraint) : applyAdjointHessian not implemented!");
139 }
140 
141 template<typename Real>
143  Vector<Real> &v2,
144  const Vector<Real> &b1,
145  const Vector<Real> &b2,
146  const Vector<Real> &x, Real &tol) {
147  StdVector<Real> v1s = dynamic_cast<StdVector<Real>&>(v1);
148  StdVector<Real> v2s = dynamic_cast<StdVector<Real>&>(v2);
149  const StdVector<Real> b1s = dynamic_cast<const StdVector<Real>&>(b1);
150  const StdVector<Real> b2s = dynamic_cast<const StdVector<Real>&>(b2);
151  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
152  try {
153  return solveAugmentedSystem( *(v1s.getVector()), *(v2s.getVector()), *(b1s.getVector()),
154  *(b2s.getVector()), *(xs.getVector()), tol );
155  }
156  catch (std::exception &e) {
157  return Constraint<Real>::solveAugmentedSystem(v1,v2,b1,b2,x,tol);
158  }
159 }
160 
161 template<typename Real>
162 std::vector<Real> StdConstraint<Real>::solveAugmentedSystem( std::vector<Real> &v1,
163  std::vector<Real> &v2,
164  const std::vector<Real> &b1,
165  const std::vector<Real> &b2,
166  const std::vector<Real> &x, Real tol ) {
167  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
168  ">>> ERROR (ROL::StdConstraint) : solveAugmentedSystem not implemented!");
169  return std::vector<Real>();
170 }
171 
172 template<typename Real>
174  const Vector<Real> &v,
175  const Vector<Real> &x,
176  const Vector<Real> &g, Real &tol) {
177  StdVector<Real> pvs = dynamic_cast<StdVector<Real>&>(pv);
178  const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
179  const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
180  const StdVector<Real> gs = dynamic_cast<const StdVector<Real>&>(g);
181  try {
182  applyPreconditioner( *(pvs.getVector()), *(vs.getVector()), *(xs.getVector()),
183  *(gs.getVector()), tol );
184  }
185  catch (std::exception &e) {
187  }
188 }
189 
190 template<typename Real>
191 void StdConstraint<Real>::applyPreconditioner( std::vector<Real> &pv,
192  const std::vector<Real> &v,
193  const std::vector<Real> &x,
194  const std::vector<Real> &g, Real &tol ) {
195  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
196  ">>> ERROR (ROL::StdConstraint) : applyPreconditioner not implemented!");
197 }
198 
199 } // namespace ROL
200 
201 #endif
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
virtual void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
Ptr< const std::vector< Element > > getVector() const
ROL::Objective_SimOpt value
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
void applyJacobian(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &sol)
std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol) override
Approximately solves the augmented system where , , , , is an identity or Riesz operator...
void update(const Vector< Real > &x, bool flag=true, int iter=-1) override
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
virtual void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ...
virtual void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
virtual std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
Approximately solves the augmented system where , , , , is an identity or Riesz operator...
void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol) override
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the constraint Jacobian at , , to vector .
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the adjoint of the the constraint Jacobian at , , to vector .