ROL
function/constraint/test_02.cpp
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 
16 #include "ROL_StdVector.hpp"
17 #include "ROL_RandomVector.hpp"
18 #include "ROL_Stream.hpp"
19 #include "Teuchos_GlobalMPISession.hpp"
20 
21 #include <iostream>
22 
23 typedef double RealT;
24 
25 int main(int argc, char *argv[]) {
26 
27  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
28 
29  using V = ROL::Vector<RealT>;
30 
31 
32 
33  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
34  int iprint = argc - 1;
35  ROL::Ptr<std::ostream> outStream;
36  ROL::nullstream bhs; // outputs nothing
37  if (iprint > 0)
38  outStream = ROL::makePtrFromRef(std::cout);
39  else
40  outStream = ROL::makePtrFromRef(bhs);
41 
42  // Save the format state of the original std::cout.
43  ROL::nullstream oldFormatState;
44  oldFormatState.copyfmt(std::cout);
45 
46  int errorFlag = 0;
47 
48  // *** Test body.
49 
50  try {
51 
52  int xdim = 5;
53 
54  // Optimization vector
55  ROL::Ptr<V> a = ROL::makePtr<ROL::StdVector<RealT>>( ROL::makePtr<std::vector<RealT>>(xdim) );
56  ROL::Ptr<V> c = ROL::makePtr<ROL::SingletonVector<RealT>>( 0.0 );
57 
58  ROL::Ptr<V> x = a->clone();
59  ROL::Ptr<V> d = x->clone();
60  ROL::Ptr<V> v = c->clone();
61 
62  RealT b = 0.5;
63 
68 
69  std::cout << "a = "; a->print(*outStream); std::cout << std::endl;
70  std::cout << "x = "; x->print(*outStream); std::cout << std::endl;
71  std::cout << "d = "; d->print(*outStream); std::cout << std::endl;
72  std::cout << "v = "; v->print(*outStream); std::cout << std::endl;
73 
75 
76  con.checkApplyJacobian(*x,*d,*c,true,*outStream);
77  con.checkAdjointConsistencyJacobian(*v,*d,*x,true,*outStream);
78  con.checkApplyAdjointHessian(*x,*v,*d,*x,true,*outStream);
79 
80 
81 
82  }
83  catch (std::logic_error& err) {
84  *outStream << err.what() << "\n";
85  errorFlag = -1000;
86  }; // end try
87 
88  if (errorFlag != 0)
89  std::cout << "End Result: TEST FAILED\n";
90  else
91  std::cout << "End Result: TEST PASSED\n";
92 
93  // reset format state of std::cout
94  std::cout.copyfmt(oldFormatState);
95 
96  return 0;
97 
98 }
99 
virtual Real checkAdjointConsistencyJacobian(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout)
void RandomizeVector(Vector< Real > &x, const Real &lower=0.0, const Real &upper=1.0)
Fill a ROL::Vector with uniformly-distributed random numbers in the interval [lower,upper].
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
This equality constraint defines an affine hyperplane.
Vector< Real > V
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian. ...
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the constraint Jacobian application.
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:38
int main(int argc, char *argv[])