ROL
ROL_StdConstraint_Def.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_STDEQUALITY_CONSTRAINT_DEF_H
11#define ROL_STDEQUALITY_CONSTRAINT_DEF_H
12
13namespace ROL {
14
15template<typename Real>
16void StdConstraint<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
17 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
18 update(*(xs.getVector()),flag,iter);
19}
20
21template<typename Real>
22void StdConstraint<Real>::update( const Vector<Real> &x, UpdateType type, int iter ) {
23 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
24 update(*(xs.getVector()),type,iter);
25}
26
27template<typename Real>
29 StdVector<Real> cs = dynamic_cast<StdVector<Real>&>(c);
30 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
31 value(*(cs.getVector()),*(xs.getVector()),tol);
33
34template<typename Real>
36 const Vector<Real> &v,
37 const Vector<Real> &x, Real &tol) {
38 StdVector<Real> jvs = dynamic_cast<StdVector<Real>&>(jv);
39 const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
40 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
41 try {
42 applyJacobian(*(jvs.getVector()),*(vs.getVector()),*(xs.getVector()),tol);
43 }
44 catch (std::exception &e ){
46 }
47}
49template<typename Real>
50void StdConstraint<Real>::applyJacobian( std::vector<Real> &jv,
51 const std::vector<Real> &v,
52 const std::vector<Real> &x, Real &tol ) {
53 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
54 ">>> ERROR (ROL::StdConstraint): applyJacobian not implemented!");
55}
56
57template<typename Real>
59 const Vector<Real> &v,
60 const Vector<Real> &x, Real &tol) {
61 StdVector<Real> ajvs = dynamic_cast<StdVector<Real>&>(ajv);
62 const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
63 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
64 try {
65 applyAdjointJacobian(*(ajvs.getVector()),*(vs.getVector()),*(xs.getVector()),tol);
66 }
67 catch (std::exception &e ){
69 }
70}
71
72template<typename Real>
73void StdConstraint<Real>::applyAdjointJacobian( std::vector<Real> &ajv,
74 const std::vector<Real> &v,
75 const std::vector<Real> &x, Real &tol ) {
76 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
77 ">>> ERROR (ROL::StdConstraint): applyAdjointJacobian not implemented!");
78}
79
80template<typename Real>
82 const Vector<Real> &u,
83 const Vector<Real> &v,
84 const Vector<Real> &x, Real &tol) {
85 StdVector<Real> ahuvs = dynamic_cast<StdVector<Real>&>(ahuv);
86 const StdVector<Real> us = dynamic_cast<const StdVector<Real>&>(u);
87 const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
88 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
89 try {
90 applyAdjointHessian( *(ahuvs.getVector()), *(us.getVector()), *(vs.getVector()),
91 *(xs.getVector()), tol );
92 }
93 catch (std::exception &e) {
95 }
96}
97
98template<typename Real>
99void StdConstraint<Real>::applyAdjointHessian( std::vector<Real> &ahuv,
100 const std::vector<Real> &u,
101 const std::vector<Real> &v,
102 const std::vector<Real> &x, Real &tol ) {
103 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
104 ">>> ERROR (ROL::StdConstraint) : applyAdjointHessian not implemented!");
105}
106
107template<typename Real>
109 Vector<Real> &v2,
110 const Vector<Real> &b1,
111 const Vector<Real> &b2,
112 const Vector<Real> &x, Real &tol) {
113 StdVector<Real> v1s = dynamic_cast<StdVector<Real>&>(v1);
114 StdVector<Real> v2s = dynamic_cast<StdVector<Real>&>(v2);
115 const StdVector<Real> b1s = dynamic_cast<const StdVector<Real>&>(b1);
116 const StdVector<Real> b2s = dynamic_cast<const StdVector<Real>&>(b2);
117 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
118 try {
119 return solveAugmentedSystem( *(v1s.getVector()), *(v2s.getVector()), *(b1s.getVector()),
120 *(b2s.getVector()), *(xs.getVector()), tol );
121 }
122 catch (std::exception &e) {
123 return Constraint<Real>::solveAugmentedSystem(v1,v2,b1,b2,x,tol);
124 }
125}
126
127template<typename Real>
128std::vector<Real> StdConstraint<Real>::solveAugmentedSystem( std::vector<Real> &v1,
129 std::vector<Real> &v2,
130 const std::vector<Real> &b1,
131 const std::vector<Real> &b2,
132 const std::vector<Real> &x, Real tol ) {
133 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
134 ">>> ERROR (ROL::StdConstraint) : solveAugmentedSystem not implemented!");
135 return std::vector<Real>();
136}
137
138template<typename Real>
140 const Vector<Real> &v,
141 const Vector<Real> &x,
142 const Vector<Real> &g, Real &tol) {
143 StdVector<Real> pvs = dynamic_cast<StdVector<Real>&>(pv);
144 const StdVector<Real> vs = dynamic_cast<const StdVector<Real>&>(v);
145 const StdVector<Real> xs = dynamic_cast<const StdVector<Real>&>(x);
146 const StdVector<Real> gs = dynamic_cast<const StdVector<Real>&>(g);
147 try {
148 applyPreconditioner( *(pvs.getVector()), *(vs.getVector()), *(xs.getVector()),
149 *(gs.getVector()), tol );
150 }
151 catch (std::exception &e) {
153 }
154}
155
156template<typename Real>
158 const std::vector<Real> &v,
159 const std::vector<Real> &x,
160 const std::vector<Real> &g, Real &tol ) {
161 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
162 ">>> ERROR (ROL::StdConstraint) : applyPreconditioner not implemented!");
163}
164
165} // namespace ROL
166
167#endif
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 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 applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply 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
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
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 i...
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 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 .
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...
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
virtual void applyJacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the constraint Jacobian at , , to vector .
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Ptr< const std::vector< Element > > getVector() const
Defines the linear algebra or vector space interface.
ROL::Objective_SerialSimOpt Objective_SimOpt value(const V &u, const V &z, Real &tol) override
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override