ROL
Loading...
Searching...
No Matches
ROL_SeparableConstraint_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_SEPARABLE_CONSTRAINT_DEF_H
11#define ROL_SEPARABLE_CONSTRAINT_DEF_H
12
13namespace ROL {
14
15template<typename Real>
17 : con_(cvec), size_(cvec.size()) {}
18
19template<typename Real>
20Ptr<Constraint<Real>> SeparableConstraint<Real>::get(unsigned ind) const {
21 ROL_TEST_FOR_EXCEPTION( (size_ <= ind ) , std::invalid_argument,
22 ">>> ERROR (ROL_SeparableConstraint, get): "
23 "Dimension mismatch.");
24 return con_[ind];
25}
26
27template<typename Real>
29 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
30 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
31 ">>> ERROR (ROL_SeparableConstraint, update): "
32 "Dimension mismatch.");
33 for(unsigned i = 0; i < size_; ++i) con_[i]->update(*xp.get(i),type,iter);
34}
35
36template<typename Real>
37void SeparableConstraint<Real>::update( const Vector<Real> &x, bool flag, int iter ) {
38 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
39 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
40 ">>> ERROR (ROL_SeparableConstraint, update): "
41 "Dimension mismatch.");
42 for(unsigned i = 0; i < size_; ++i) con_[i]->update(*xp.get(i),flag,iter);
43}
44
45template<typename Real>
48 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
49 ROL_TEST_FOR_EXCEPTION( (size_ != cp.numVectors() ) , std::invalid_argument,
50 ">>> ERROR (ROL_SeparableConstraint, value): "
51 "Dimension mismatch.");
52 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
53 ">>> ERROR (ROL_SeparableConstraint, value): "
54 "Dimension mismatch.");
55 for(unsigned i = 0; i < size_; ++i) con_[i]->value(*cp.get(i),*xp.get(i),tol);
56}
57
58template<typename Real>
60 const Vector<Real> &v,
61 const Vector<Real> &x,
62 Real &tol ) {
63 PartitionedVector<Real> &jvp = dynamic_cast<PartitionedVector<Real>&>(jv);
64 const PartitionedVector<Real> &vp = dynamic_cast<const PartitionedVector<Real>&>(v);
65 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
66 ROL_TEST_FOR_EXCEPTION( (size_ != jvp.numVectors() ) , std::invalid_argument,
67 ">>> ERROR (ROL_SeparableConstraint, applyJacobian): "
68 "Dimension mismatch.");
69 ROL_TEST_FOR_EXCEPTION( (size_ != vp.numVectors() ) , std::invalid_argument,
70 ">>> ERROR (ROL_SeparableConstraint, applyJacobian): "
71 "Dimension mismatch.");
72 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
73 ">>> ERROR (ROL_SeparableConstraint, applyJacobian): "
74 "Dimension mismatch.");
75 for(unsigned i = 0; i < size_; ++i) con_[i]->applyJacobian(*jvp.get(i),*vp.get(i),*xp.get(i),tol);
76}
77
78template<typename Real>
80 const Vector<Real> &v,
81 const Vector<Real> &x,
82 Real &tol ) {
83 PartitionedVector<Real> &ajvp = dynamic_cast<PartitionedVector<Real>&>(ajv);
84 const PartitionedVector<Real> &vp = dynamic_cast<const PartitionedVector<Real>&>(v);
85 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
86 ROL_TEST_FOR_EXCEPTION( (size_ != ajvp.numVectors() ) , std::invalid_argument,
87 ">>> ERROR (ROL_SeparableConstraint, applyAdjointJacobian): "
88 "Dimension mismatch.");
89 ROL_TEST_FOR_EXCEPTION( (size_ != vp.numVectors() ) , std::invalid_argument,
90 ">>> ERROR (ROL_SeparableConstraint, applyAdjointJacobian): "
91 "Dimension mismatch.");
92 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
93 ">>> ERROR (ROL_SeparableConstraint, applyAdjointJacobian): "
94 "Dimension mismatch.");
95 for(unsigned i = 0; i < size_; ++i) con_[i]->applyAdjointJacobian(*ajvp.get(i),*vp.get(i),*xp.get(i),tol);
96}
97
98template<typename Real>
100 const Vector<Real> &u,
101 const Vector<Real> &v,
102 const Vector<Real> &x,
103 Real &tol ) {
104 PartitionedVector<Real> &ahuvp = dynamic_cast<PartitionedVector<Real>&>(ahuv);
105 const PartitionedVector<Real> &up = dynamic_cast<const PartitionedVector<Real>&>(u);
106 const PartitionedVector<Real> &vp = dynamic_cast<const PartitionedVector<Real>&>(v);
107 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
108 ROL_TEST_FOR_EXCEPTION( (size_ != ahuvp.numVectors() ) , std::invalid_argument,
109 ">>> ERROR (ROL_SeparableConstraint, applyAdjointHessian): "
110 "Dimension mismatch.");
111 ROL_TEST_FOR_EXCEPTION( (size_ != up.numVectors() ) , std::invalid_argument,
112 ">>> ERROR (ROL_SeparableConstraint, applyAdjointHessian): "
113 "Dimension mismatch.");
114 ROL_TEST_FOR_EXCEPTION( (size_ != vp.numVectors() ) , std::invalid_argument,
115 ">>> ERROR (ROL_SeparableConstraint, applyAdjointHessian): "
116 "Dimension mismatch.");
117 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
118 ">>> ERROR (ROL_SeparableConstraint, applyAdjointHessian): "
119 "Dimension mismatch.");
120 for(unsigned i = 0; i < size_; ++i) con_[i]->applyAdjointHessian(*ahuvp.get(i),*up.get(i),*vp.get(i),*xp.get(i),tol);
121}
122
123template<typename Real>
125 const Vector<Real> &v,
126 const Vector<Real> &x,
127 const Vector<Real> &g,
128 Real &tol) {
129 PartitionedVector<Real> &pvp = dynamic_cast<PartitionedVector<Real>&>(pv);
130 const PartitionedVector<Real> &vp = dynamic_cast<const PartitionedVector<Real>&>(v);
131 const PartitionedVector<Real> &xp = dynamic_cast<const PartitionedVector<Real>&>(x);
132 const PartitionedVector<Real> &gp = dynamic_cast<const PartitionedVector<Real>&>(g);
133 ROL_TEST_FOR_EXCEPTION( (size_ != pvp.numVectors() ) , std::invalid_argument,
134 ">>> ERROR (ROL_SeparableConstraint, applyPreconditioner): "
135 "Dimension mismatch.");
136 ROL_TEST_FOR_EXCEPTION( (size_ != vp.numVectors() ) , std::invalid_argument,
137 ">>> ERROR (ROL_SeparableConstraint, applyPreconditioner): "
138 "Dimension mismatch.");
139 ROL_TEST_FOR_EXCEPTION( (size_ != xp.numVectors() ) , std::invalid_argument,
140 ">>> ERROR (ROL_SeparableConstraint, applyPreconditioner): "
141 "Dimension mismatch.");
142 ROL_TEST_FOR_EXCEPTION( (size_ != gp.numVectors() ) , std::invalid_argument,
143 ">>> ERROR (ROL_SeparableConstraint, applyPreconditioner): "
144 "Dimension mismatch.");
145 for(unsigned i = 0; i < size_; ++i) con_[i]->applyPreconditioner(*pvp.get(i),*vp.get(i),*xp.get(i),*gp.get(i),tol);
146}
147
148template<typename Real>
149void SeparableConstraint<Real>::setParameter(const std::vector<Real> &param) {
151 for(unsigned i = 0; i < size_; ++i) con_[i]->setParameter(param);
152}
153
154} // namespace ROL
155
156#endif
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
Ptr< Constraint< Real > > get(unsigned ind=0) const
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 setParameter(const std::vector< Real > &param) override
void update(const Vector< Real > &x, UpdateType type, int iter=-1) override
Update constraint function.
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 ,...
SeparableConstraint(const std::vector< Ptr< Constraint< Real > > > &cvec)
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 applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the constraint Jacobian at , , to vector .
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
Defines the linear algebra or vector space interface.
void setParameter(ROL::ParameterList &parlist, const std::vector< std::string > &location, const std::vector< std::string >::iterator iter, ParameterType value)
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
void applyJacobian(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &sol)