55 #include "Teuchos_GlobalMPISession.hpp"
57 template<
typename Real>
60 void value(std::vector<Real> &c,
const std::vector<Real> &x, Real &tol) {
63 void applyJacobian(std::vector<Real> &jv,
const std::vector<Real> &v,
const std::vector<Real> &x, Real &tol) {
66 void applyAdjointJacobian(std::vector<Real> &ajv,
const std::vector<Real> &v,
const std::vector<Real> &x, Real &tol) {
74 int main(
int argc,
char *argv[]) {
76 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
80 int iprint = argc - 1;
81 ROL::Ptr<std::ostream> outStream;
84 outStream = ROL::makePtrFromRef(std::cout);
86 outStream = ROL::makePtrFromRef(bhs);
91 RealT tol = std::sqrt(ROL::ROL_EPSILON<RealT>());
93 ROL::Ptr<con2d<RealT>> con = ROL::makePtr<con2d<RealT>>();
95 ROL::Ptr<std::vector<RealT>> yptr = ROL::makePtr<std::vector<RealT>>(2);
96 (*yptr)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
97 (*yptr)[1] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
102 ROL::Ptr<std::vector<RealT>> xptr = ROL::makePtr<std::vector<RealT>>(2);
103 (*xptr)[0] = (*yptr)[0];
104 (*xptr)[1] = (*yptr)[1];
107 ROL::Ptr<std::vector<RealT>> Pxptr = ROL::makePtr<std::vector<RealT>>(2,0.0);
113 ROL::Ptr<std::vector<RealT>> x0ptr = ROL::makePtr<std::vector<RealT>>(2);
114 (*x0ptr)[0] = ((*yptr)[0]-(*yptr)[1])/static_cast<RealT>(2);
115 (*x0ptr)[1] = -(*x0ptr)[0];
120 *outStream << std::setprecision(6) << std::scientific << std::endl;
121 *outStream <<
" x[0] = " << (*xptr)[0] <<
" x[1] = " << (*xptr)[1] << std::endl;
122 *outStream <<
" Px[0] = " << (*Pxptr)[0] <<
" Px[1] = " << (*Pxptr)[1] << std::endl;
123 *outStream <<
" x*[0] = " << (*x0ptr)[0] <<
" x*[1] = " << (*x0ptr)[1] << std::endl;
125 e0.
set(x0); e0.
axpy(static_cast<RealT>(-1),Px);
127 *outStream <<
" Error in Euclidean Projection: " << err << std::endl;
129 e0.
set(x); e0.
axpy(static_cast<RealT>(-1),x0);
130 *outStream <<
" ||x*-x||^2 = " << e0.
norm() << std::endl;
132 e0.
set(x); e0.
axpy(static_cast<RealT>(-1),Px);
133 *outStream <<
" ||Px-x||^2 = " << e0.
norm() << std::endl << std::endl;
135 errorFlag += (err > tol);
137 ROL::Ptr<std::vector<RealT>> dptr = ROL::makePtr<std::vector<RealT>>(2);
138 (*dptr)[0] =
static_cast<RealT>(1)+static_cast<RealT>(2)*
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
139 (*dptr)[1] =
static_cast<RealT>(1)+static_cast<RealT>(5)*
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
141 ROL::Ptr<std::vector<RealT>> x1ptr = ROL::makePtr<std::vector<RealT>>(2);
142 (*x1ptr)[0] = ((*dptr)[0]*(*yptr)[0]-(*dptr)[1]*(*yptr)[1])/((*dptr)[0]+(*dptr)[1]);
143 (*x1ptr)[1] = -(*x1ptr)[0];
146 ROL::Ptr<std::vector<RealT>> zptr = ROL::makePtr<std::vector<RealT>>(2);
147 (*zptr)[0] = (*yptr)[0];
148 (*zptr)[1] = (*yptr)[1];
151 ROL::Ptr<std::vector<RealT>> Pzptr = ROL::makePtr<std::vector<RealT>>(2,0.0);
157 ROL::Ptr<std::vector<RealT>> e1ptr = ROL::makePtr<std::vector<RealT>>(2);
160 *outStream << std::endl;
161 *outStream <<
" x[0] = " << (*zptr)[0] <<
" x[1] = " << (*zptr)[1] << std::endl;
162 *outStream <<
" Px[0] = " << (*Pzptr)[0] <<
" Px[1] = " << (*Pzptr)[1] << std::endl;
163 *outStream <<
" x*[0] = " << (*x1ptr)[0] <<
" x*[1] = " << (*x1ptr)[1] << std::endl;
165 e1.set(x1); e1.axpy(static_cast<RealT>(-1),Pz);
167 *outStream <<
" Error in Euclidean Projection: " << err << std::endl;
169 e1.set(z); e1.axpy(static_cast<RealT>(-1),x1);
170 *outStream <<
" ||x*-x||^2 = " << e1.norm() << std::endl;
172 e1.set(z); e1.axpy(static_cast<RealT>(-1),Pz);
173 *outStream <<
" ||Px-x||^2 = " << e1.norm() << std::endl << std::endl;
175 errorFlag += (err > tol);
178 catch (std::logic_error& err) {
179 *outStream << err.what() <<
"\n";
184 std::cout <<
"End Result: TEST FAILED\n";
186 std::cout <<
"End Result: TEST PASSED\n";
void axpy(const Real alpha, const Vector< Real > &x)
Compute where .
Defines the equality constraint operator interface for StdVectors.
virtual void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Projects on to the null space of a linear constraint.
Real norm() const
Returns where .
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
void applyJacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
Provides the std::vector implementation of the ROL::Vector interface that handles scalings in the inn...
void set(const Vector< Real > &x)
Set where .
basic_nullstream< char, char_traits< char >> nullstream
int main(int argc, char *argv[])
void applyAdjointJacobian(std::vector< Real > &ajv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
void value(std::vector< Real > &c, const std::vector< Real > &x, Real &tol)