ROL
ROL_Constraint_Partitioned_Def.hpp
Go to the documentation of this file.
1// Rapid Optimization Library (ROL) Package
2// Copyright (2014) Sandia Corporation
3//
4// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
5// license for use of this work by or on behalf of the U.S. Government.
6//
7// Redistribution and use in source and binary forms, with or without
8// modification, are permitted provided that the following conditions are
9// met:
10//
11// 1. Redistributions of source code must retain the above copyright
12// notice, this list of conditions and the following disclaimer.
13//
14// 2. Redistributions in binary form must reproduce the above copyright
15// notice, this list of conditions and the following disclaimer in the
16// documentation and/or other materials provided with the distribution.
17//
18// 3. Neither the name of the Corporation nor the names of the
19// contributors may be used to endorse or promote products derived from
20// this software without specific prior written permission.
21//
22// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
23// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
24// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
25// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
26// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
27// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
28// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
29// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
30// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
31// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
32// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33//
34// Questions? Contact lead developers:
35// Drew Kouri (dpkouri@sandia.gov) and
36// Denis Ridzal (dridzal@sandia.gov)
37//
38// ************************************************************************
39// @HEADER
40
41#ifndef ROL_CONSTRAINT_PARTITIONED_DEF_H
42#define ROL_CONSTRAINT_PARTITIONED_DEF_H
43
44namespace ROL {
45
46template<typename Real>
48 bool isInequality,
49 int offset)
50 : cvec_(cvec), offset_(offset),
51 scratch_(nullPtr), ncval_(0), initialized_(false) {
52 isInequality_.clear(); isInequality_.resize(cvec.size(),isInequality);
53}
54
55template<typename Real>
57 std::vector<bool> isInequality,
58 int offset)
59 : cvec_(cvec), isInequality_(isInequality), offset_(offset),
60 scratch_(nullPtr), ncval_(0), initialized_(false) {}
61
62template<typename Real>
64 return ncval_;
65}
66
67template<typename Real>
68Ptr<Constraint<Real>> Constraint_Partitioned<Real>::get(int ind) const {
69 if (ind < 0 || ind > static_cast<int>(cvec_.size())) {
70 throw Exception::NotImplemented(">>> Constraint_Partitioned::get : Index out of bounds!");
71 }
72 return cvec_[ind];
73}
74
75template<typename Real>
77 const int ncon = static_cast<int>(cvec_.size());
78 for (int i = 0; i < ncon; ++i) {
79 cvec_[i]->update(getOpt(x),type,iter);
80 }
81}
82
83template<typename Real>
84void Constraint_Partitioned<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
85 const int ncon = static_cast<int>(cvec_.size());
86 for (int i = 0; i < ncon; ++i) {
87 cvec_[i]->update(getOpt(x),flag,iter);
88 }
89}
90
91template<typename Real>
94 = dynamic_cast<PartitionedVector<Real>&>(c);
95
96 const int ncon = static_cast<int>(cvec_.size());
97 int cnt = offset_+1;
98 for (int i = 0; i < ncon; ++i) {
99 cvec_[i]->value(*cpv.get(i), getOpt(x), tol);
100 if (isInequality_[i]) {
101 cpv.get(i)->axpy(static_cast<Real>(-1),getSlack(x,cnt));
102 cnt++;
103 }
104 }
105 ++ncval_;
106}
107
108template<typename Real>
110 const Vector<Real> &v,
111 const Vector<Real> &x,
112 Real &tol ) {
114 = dynamic_cast<PartitionedVector<Real>&>(jv);
115
116 const int ncon = static_cast<int>(cvec_.size());
117 int cnt = offset_+1;
118 for (int i = 0; i < ncon; ++i) {
119 cvec_[i]->applyJacobian(*jvpv.get(i), getOpt(v), getOpt(x), tol);
120 if (isInequality_[i]) {
121 jvpv.get(i)->axpy(static_cast<Real>(-1),getSlack(v,cnt));
122 cnt++;
123 }
124 }
125}
126
127template<typename Real>
129 const Vector<Real> &v,
130 const Vector<Real> &x,
131 Real &tol ) {
132 if (!initialized_) {
133 scratch_ = getOpt(ajv).clone();
134 initialized_ = true;
135 }
136
137 const PartitionedVector<Real> &vpv
138 = dynamic_cast<const PartitionedVector<Real>&>(v);
139
140 const int ncon = static_cast<int>(cvec_.size());
141 int cnt = offset_+1;
142 getOpt(ajv).zero();
143 for (int i = 0; i < ncon; ++i) {
144 scratch_->zero();
145 cvec_[i]->applyAdjointJacobian(*scratch_, *vpv.get(i), getOpt(x), tol);
146 getOpt(ajv).plus(*scratch_);
147 if (isInequality_[i]) {
148 getSlack(ajv,cnt).set(*vpv.get(i));
149 getSlack(ajv,cnt).scale(static_cast<Real>(-1));
150 cnt++;
151 }
152 }
153}
154
155template<typename Real>
157 const Vector<Real> &u,
158 const Vector<Real> &v,
159 const Vector<Real> &x,
160 Real &tol ) {
161 if (!initialized_) {
162 scratch_ = getOpt(ahuv).clone();
163 initialized_ = true;
164 }
165
166 const PartitionedVector<Real> &upv
167 = dynamic_cast<const PartitionedVector<Real>&>(u);
168
169 const int ncon = static_cast<int>(cvec_.size());
170 int cnt = offset_+1;
171 getOpt(ahuv).zero();
172 for (int i = 0; i < ncon; ++i) {
173 Ptr<const Vector<Real>> ui = upv.get(i);
174 scratch_->zero();
175 cvec_[i]->applyAdjointHessian(*scratch_, *ui, getOpt(v), getOpt(x), tol);
176 getOpt(ahuv).plus(*scratch_);
177 if (isInequality_[i]) {
178 getSlack(ahuv,cnt).zero();
179 cnt++;
180 }
181 }
182}
183
184template<typename Real>
186 const Vector<Real> &v,
187 const Vector<Real> &x,
188 const Vector<Real> &g,
189 Real &tol) {
191 = dynamic_cast<PartitionedVector<Real>&>(pv);
192 const PartitionedVector<Real> &vpv
193 = dynamic_cast<const PartitionedVector<Real>&>(v);
194
195 const int ncon = static_cast<int>(cvec_.size());
196 for (int i = 0; i < ncon; ++i) {
197 cvec_[i]->applyPreconditioner(*pvpv.get(i), *vpv.get(i), getOpt(x), getOpt(g), tol);
198 }
199}
200
201template<typename Real>
202void Constraint_Partitioned<Real>::setParameter(const std::vector<Real> &param) {
204 const int ncon = static_cast<int>(cvec_.size());
205 for (int i = 0; i < ncon; ++i) {
206 cvec_[i]->setParameter(param);
207 }
208}
209
210template<typename Real>
212 try {
213 return *dynamic_cast<PartitionedVector<Real>&>(xs).get(0);
214 }
215 catch (std::exception &e) {
216 return xs;
217 }
218}
219
220template<typename Real>
222 try {
223 return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(0);
224 }
225 catch (std::exception &e) {
226 return xs;
227 }
228}
229
230template<typename Real>
232 return *dynamic_cast<PartitionedVector<Real>&>(xs).get(ind);
233}
234
235template<typename Real>
237 return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(ind);
238}
239
240} // namespace ROL
241
242#endif
void setParameter(const std::vector< Real > &param) override
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 .
Vector< Real > & getSlack(Vector< Real > &xs, int ind) const
Vector< Real > & getOpt(Vector< Real > &xs) const
Constraint_Partitioned(const std::vector< Ptr< Constraint< Real > > > &cvec, bool isInequality=false, int offset=0)
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
Ptr< Constraint< Real > > get(int ind=0) const
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 .
virtual 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 update(const Vector< Real > &x, UpdateType type, int iter=-1) override
Update constraint function.
Defines the general constraint operator interface.
virtual void setParameter(const std::vector< Real > &param)
Defines the linear algebra of vector space on a generic partitioned vector.
ROL::Ptr< const Vector< Real > > get(size_type i) const
void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:84