ROL
ROL_Constraint_Partitioned.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_H
42 #define ROL_CONSTRAINT_PARTITIONED_H
43 
45 #include "ROL_Constraint.hpp"
46 
47 namespace ROL {
48 
55 template<class Real>
56 class Constraint_Partitioned : public Constraint<Real> {
57 private:
58  std::vector<ROL::Ptr<Constraint<Real> > > cvec_;
59  std::vector<bool> isInequality_; // Label whether cvec_[i] is inequality
60  ROL::Ptr<Vector<Real> > scratch_; // Scratch vector for intermediate computation
61  int ncval_; // Number of constraint evaluations
62  bool initialized_; // Is scratch vector initialized?
63 
65  try {
66  return *dynamic_cast<PartitionedVector<Real>&>(xs).get(0);
67  }
68  catch (std::exception &e) {
69  return xs;
70  }
71  }
72 
73  const Vector<Real>& getOpt( const Vector<Real> &xs ) {
74  try {
75  return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(0);
76  }
77  catch (std::exception &e) {
78  return xs;
79  }
80  }
81 
82  Vector<Real>& getSlack( Vector<Real> &xs, const int ind ) {
83  return *dynamic_cast<PartitionedVector<Real>&>(xs).get(ind);
84  }
85 
86  const Vector<Real>& getSlack( const Vector<Real> &xs, const int ind ) {
87  return *dynamic_cast<const PartitionedVector<Real>&>(xs).get(ind);
88  }
89 
90 
91 public:
92  Constraint_Partitioned(const std::vector<ROL::Ptr<Constraint<Real> > > &cvec,
93  bool isInequality = false)
94  : cvec_(cvec),
95  scratch_(ROL::nullPtr), ncval_(0), initialized_(false) {
96  isInequality_.clear(); isInequality_.resize(cvec.size(),isInequality);
97  }
98 
99  Constraint_Partitioned(const std::vector<ROL::Ptr<Constraint<Real> > > &cvec,
100  const std::vector<bool> &isInequality)
101  : cvec_(cvec), isInequality_(isInequality),
102  scratch_(ROL::nullPtr), ncval_(0), initialized_(false) {}
103 
105  return ncval_;
106  }
107 
108  Ptr<Constraint<Real>> get(const int ind = 0) const {
109  if (ind < 0 || ind > static_cast<int>(cvec_.size())) {
110  throw Exception::NotImplemented(">>> Constraint_Partitioned::get : Index out of bounds!");
111  }
112  return cvec_[ind];
113  }
114 
115  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
116  const int ncon = static_cast<int>(cvec_.size());
117  for (int i = 0; i < ncon; ++i) {
118  cvec_[i]->update(getOpt(x),flag,iter);
119  }
120  }
121 
122  void value( Vector<Real> &c, const Vector<Real> &x, Real &tol ) {
124  = dynamic_cast<PartitionedVector<Real>&>(c);
125 
126  const int ncon = static_cast<int>(cvec_.size());
127  int cnt = 1;
128  for (int i = 0; i < ncon; ++i) {
129  cvec_[i]->value(*cpv.get(i), getOpt(x), tol);
130  if (isInequality_[i]) {
131  cpv.get(i)->axpy(static_cast<Real>(-1),getSlack(x,cnt));
132  cnt++;
133  }
134  }
135  ++ncval_;
136  }
137 
139  const Vector<Real> &v,
140  const Vector<Real> &x,
141  Real &tol ) {
143  = dynamic_cast<PartitionedVector<Real>&>(jv);
144 
145  const int ncon = static_cast<int>(cvec_.size());
146  int cnt = 1;
147  for (int i = 0; i < ncon; ++i) {
148  cvec_[i]->applyJacobian(*jvpv.get(i), getOpt(v), getOpt(x), tol);
149  if (isInequality_[i]) {
150  jvpv.get(i)->axpy(static_cast<Real>(-1),getSlack(v,cnt));
151  cnt++;
152  }
153  }
154  }
155 
158  const Vector<Real> &v,
159  const Vector<Real> &x,
160  Real &tol ) {
161  if (!initialized_) {
162  scratch_ = getOpt(ajv).clone();
163  initialized_ = true;
164  }
165 
166  const PartitionedVector<Real> &vpv
167  = dynamic_cast<const PartitionedVector<Real>&>(v);
168 
169  const int ncon = static_cast<int>(cvec_.size());
170  int cnt = 1;
171  getOpt(ajv).zero();
172  for (int i = 0; i < ncon; ++i) {
173  scratch_->zero();
174  cvec_[i]->applyAdjointJacobian(*scratch_, *vpv.get(i), getOpt(x), tol);
175  getOpt(ajv).plus(*scratch_);
176  if (isInequality_[i]) {
177  getSlack(ajv,cnt).set(*vpv.get(i));
178  getSlack(ajv,cnt).scale(static_cast<Real>(-1));
179  cnt++;
180  }
181  }
182  }
183 
185  const Vector<Real> &u,
186  const Vector<Real> &v,
187  const Vector<Real> &x,
188  Real &tol ) {
189  if (!initialized_) {
190  scratch_ = getOpt(ahuv).clone();
191  initialized_ = true;
192  }
193 
194  const PartitionedVector<Real> &upv
195  = dynamic_cast<const PartitionedVector<Real>&>(u);
196 
197  const int ncon = static_cast<int>(cvec_.size());
198  int cnt = 1;
199  getOpt(ahuv).zero();
200  for (int i = 0; i < ncon; ++i) {
201  ROL::Ptr<const Vector<Real> > ui = upv.get(i);
202  scratch_->zero();
203  cvec_[i]->applyAdjointHessian(*scratch_, *ui, getOpt(v), getOpt(x), tol);
204  getOpt(ahuv).plus(*scratch_);
205  if (isInequality_[i]) {
206  getSlack(ahuv,cnt).zero();
207  cnt++;
208  }
209  }
210  }
211 
213  const Vector<Real> &v,
214  const Vector<Real> &x,
215  const Vector<Real> &g,
216  Real &tol) {
218  = dynamic_cast<PartitionedVector<Real>&>(pv);
219  const PartitionedVector<Real> &vpv
220  = dynamic_cast<const PartitionedVector<Real>&>(v);
221 
222  const int ncon = static_cast<int>(cvec_.size());
223  for (int i = 0; i < ncon; ++i) {
224  cvec_[i]->applyPreconditioner(*pvpv.get(i), *vpv.get(i), getOpt(x), getOpt(g), tol);
225  }
226  }
227 
228 // Definitions for parametrized (stochastic) equality constraints
229 public:
230  void setParameter(const std::vector<Real> &param) {
232  const int ncon = static_cast<int>(cvec_.size());
233  for (int i = 0; i < ncon; ++i) {
234  cvec_[i]->setParameter(param);
235  }
236  }
237 }; // class Constraint_Partitioned
238 
239 } // namespace ROL
240 
241 #endif
ROL::Ptr< const Vector< Real > > get(size_type i) const
Vector< Real > & getSlack(Vector< Real > &xs, const int ind)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
Defines the linear algebra of vector space on a generic partitioned vector.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Constraint_Partitioned(const std::vector< ROL::Ptr< Constraint< Real > > > &cvec, const std::vector< bool > &isInequality)
const Vector< Real > & getSlack(const Vector< Real > &xs, const int ind)
Has both inequality and equality constraints. Treat inequality constraint as equality with slack vari...
virtual void setParameter(const std::vector< Real > &param)
ROL::Ptr< Vector< Real > > scratch_
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...
Vector< Real > & getOpt(Vector< Real > &xs)
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 .
const Vector< Real > & getOpt(const Vector< Real > &xs)
Constraint_Partitioned(const std::vector< ROL::Ptr< Constraint< Real > > > &cvec, bool isInequality=false)
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
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 ...
std::vector< ROL::Ptr< Constraint< Real > > > cvec_
void setParameter(const std::vector< Real > &param)
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
Defines the general constraint operator interface.