10#ifndef ROL_SEPARABLE_CONSTRAINT_DEF_H
11#define ROL_SEPARABLE_CONSTRAINT_DEF_H
15template<
typename Real>
17 : con_(cvec), size_(cvec.size()) {}
19template<
typename Real>
21 ROL_TEST_FOR_EXCEPTION( (size_ <= ind ) , std::invalid_argument,
22 ">>> ERROR (ROL_SeparableConstraint, get): "
23 "Dimension mismatch.");
27template<
typename Real>
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);
36template<
typename Real>
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);
45template<
typename Real>
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);
58template<
typename Real>
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.");
78template<
typename Real>
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);
98template<
typename Real>
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);
123template<
typename Real>
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);
148template<
typename Real>
151 for(
unsigned i = 0; i < size_; ++i) con_[i]->
setParameter(param);
Defines the general constraint operator interface.
virtual void setParameter(const std::vector< Real > ¶m)
Defines the linear algebra of vector space on a generic partitioned vector.
ROL::Ptr< const Vector< Real > > get(size_type i) const
size_type numVectors() 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 > ¶m) 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)