ROL
ROL_Constraint_Partitioned_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_CONSTRAINT_PARTITIONED_DEF_H
11#define ROL_CONSTRAINT_PARTITIONED_DEF_H
12
13namespace ROL {
14
15template<typename Real>
17 bool isInequality,
18 int offset)
19 : cvec_(cvec), offset_(offset),
20 scratch_(nullPtr), ncval_(0), initialized_(false) {
21 isInequality_.clear(); isInequality_.resize(cvec.size(),isInequality);
22}
23
24template<typename Real>
26 std::vector<bool> isInequality,
27 int offset)
28 : cvec_(cvec), isInequality_(isInequality), offset_(offset),
29 scratch_(nullPtr), ncval_(0), initialized_(false) {}
30
31template<typename Real>
35
36template<typename Real>
37Ptr<Constraint<Real>> Constraint_Partitioned<Real>::get(int ind) const {
38 if (ind < 0 || ind >= static_cast<int>(cvec_.size())) {
39 throw Exception::NotImplemented(">>> Constraint_Partitioned::get : Index out of bounds!");
40 }
41 return cvec_[ind];
42}
43
44template<typename Real>
46 const int ncon = static_cast<int>(cvec_.size());
47 for (int i = 0; i < ncon; ++i) {
48 cvec_[i]->update(getOpt(x),type,iter);
49 }
50}
51
52template<typename Real>
53void Constraint_Partitioned<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
54 const int ncon = static_cast<int>(cvec_.size());
55 for (int i = 0; i < ncon; ++i) {
56 cvec_[i]->update(getOpt(x),flag,iter);
57 }
58}
59
60template<typename Real>
63 = dynamic_cast<PartitionedVector<Real>&>(c);
64
65 const int ncon = static_cast<int>(cvec_.size());
66 int cnt = offset_+1;
67 for (int i = 0; i < ncon; ++i) {
68 cvec_[i]->value(*cpv.get(i), getOpt(x), tol);
69 if (isInequality_[i]) {
70 cpv.get(i)->axpy(static_cast<Real>(-1),getSlack(x,cnt));
71 cnt++;
72 }
73 }
74 ++ncval_;
75}
76
77template<typename Real>
79 const Vector<Real> &v,
80 const Vector<Real> &x,
81 Real &tol ) {
83 = dynamic_cast<PartitionedVector<Real>&>(jv);
84
85 const int ncon = static_cast<int>(cvec_.size());
86 int cnt = offset_+1;
87 for (int i = 0; i < ncon; ++i) {
88 cvec_[i]->applyJacobian(*jvpv.get(i), getOpt(v), getOpt(x), tol);
89 if (isInequality_[i]) {
90 jvpv.get(i)->axpy(static_cast<Real>(-1),getSlack(v,cnt));
91 cnt++;
92 }
93 }
94}
95
96template<typename Real>
98 const Vector<Real> &v,
99 const Vector<Real> &x,
100 Real &tol ) {
101 if (!initialized_) {
102 scratch_ = getOpt(ajv).clone();
103 initialized_ = true;
104 }
105
106 const PartitionedVector<Real> &vpv
107 = dynamic_cast<const PartitionedVector<Real>&>(v);
108
109 const int ncon = static_cast<int>(cvec_.size());
110 int cnt = offset_+1;
111 getOpt(ajv).zero();
112 for (int i = 0; i < ncon; ++i) {
113 scratch_->zero();
114 cvec_[i]->applyAdjointJacobian(*scratch_, *vpv.get(i), getOpt(x), tol);
115 getOpt(ajv).plus(*scratch_);
116 if (isInequality_[i]) {
117 getSlack(ajv,cnt).set(*vpv.get(i));
118 getSlack(ajv,cnt).scale(static_cast<Real>(-1));
119 cnt++;
120 }
121 }
122}
123
124template<typename Real>
126 const Vector<Real> &u,
127 const Vector<Real> &v,
128 const Vector<Real> &x,
129 Real &tol ) {
130 if (!initialized_) {
131 scratch_ = getOpt(ahuv).clone();
132 initialized_ = true;
133 }
134
135 const PartitionedVector<Real> &upv
136 = dynamic_cast<const PartitionedVector<Real>&>(u);
137
138 const int ncon = static_cast<int>(cvec_.size());
139 int cnt = offset_+1;
140 getOpt(ahuv).zero();
141 for (int i = 0; i < ncon; ++i) {
142 Ptr<const Vector<Real>> ui = upv.get(i);
143 scratch_->zero();
144 cvec_[i]->applyAdjointHessian(*scratch_, *ui, getOpt(v), getOpt(x), tol);
145 getOpt(ahuv).plus(*scratch_);
146 if (isInequality_[i]) {
147 getSlack(ahuv,cnt).zero();
148 cnt++;
149 }
150 }
151}
152
153template<typename Real>
155 const Vector<Real> &v,
156 const Vector<Real> &x,
157 const Vector<Real> &g,
158 Real &tol) {
160 = dynamic_cast<PartitionedVector<Real>&>(pv);
161 const PartitionedVector<Real> &vpv
162 = dynamic_cast<const PartitionedVector<Real>&>(v);
163
164 const int ncon = static_cast<int>(cvec_.size());
165 for (int i = 0; i < ncon; ++i) {
166 cvec_[i]->applyPreconditioner(*pvpv.get(i), *vpv.get(i), getOpt(x), getOpt(g), tol);
167 }
168}
169
170template<typename Real>
171void Constraint_Partitioned<Real>::setParameter(const std::vector<Real> &param) {
173 const int ncon = static_cast<int>(cvec_.size());
174 for (int i = 0; i < ncon; ++i) {
175 cvec_[i]->setParameter(param);
176 }
177}
178
179template<typename Real>
181 try {
182 return *dynamic_cast<PartitionedVector<Real>&>(xs).get(0);
183 }
184 catch (std::exception &e) {
185 return xs;
186 }
187}
188
189template<typename Real>
191 try {
192 return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(0);
193 }
194 catch (std::exception &e) {
195 return xs;
196 }
197}
198
199template<typename Real>
201 return *dynamic_cast<PartitionedVector<Real>&>(xs).get(ind);
202}
203
204template<typename Real>
206 return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(ind);
207}
208
209} // namespace ROL
210
211#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.