20 #include "Teuchos_GlobalMPISession.hpp"
21 #include "ROL_ParameterList.hpp"
26 int main(
int argc,
char *argv[]) {
28 typedef std::vector<RealT> vector;
34 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
37 int iprint = argc - 1;
38 ROL::Ptr<std::ostream> outStream;
41 outStream = ROL::makePtrFromRef(std::cout);
43 outStream = ROL::makePtrFromRef(bhs);
47 oldFormatState.copyfmt(std::cout);
61 ROL::Ptr<vector> x_ptr = ROL::makePtr<vector>(
dim,0.0);
62 ROL::Ptr<vector> g_ptr = ROL::makePtr<vector>(
dim,0.0);
63 ROL::Ptr<vector> v_ptr = ROL::makePtr<vector>(
dim,0.0);
64 ROL::Ptr<vector> hv_ptr = ROL::makePtr<vector>(
dim,0.0);
66 ROL::Ptr<vector> l_ptr = ROL::makePtr<vector>(
dim,1.0);
67 ROL::Ptr<vector> u_ptr = ROL::makePtr<vector>(
dim,2.0);
69 ROL::Ptr<vector> xlog0_ptr = ROL::makePtr<vector>(
dim,0.0);
70 ROL::Ptr<vector> xlog1_ptr = ROL::makePtr<vector>(
dim,0.0);
71 ROL::Ptr<vector> xlog2_ptr = ROL::makePtr<vector>(
dim,0.0);
73 ROL::Ptr<vector> xquad0_ptr = ROL::makePtr<vector>(
dim,0.0);
74 ROL::Ptr<vector> xquad1_ptr = ROL::makePtr<vector>(
dim,0.0);
75 ROL::Ptr<vector> xquad2_ptr = ROL::makePtr<vector>(
dim,0.0);
77 ROL::Ptr<vector> xdwell0_ptr = ROL::makePtr<vector>(
dim,0.0);
78 ROL::Ptr<vector> xdwell1_ptr = ROL::makePtr<vector>(
dim,0.0);
79 ROL::Ptr<vector> xdwell2_ptr = ROL::makePtr<vector>(
dim,0.0);
88 ROL::Ptr<SV> xlog0 = ROL::makePtr<SV>(xlog0_ptr);
89 ROL::Ptr<SV> xlog1 = ROL::makePtr<SV>(xlog1_ptr);
90 ROL::Ptr<SV> xlog2 = ROL::makePtr<SV>(xlog2_ptr);
92 ROL::Ptr<SV> xquad0 = ROL::makePtr<SV>(xquad0_ptr);
93 ROL::Ptr<SV> xquad1 = ROL::makePtr<SV>(xquad1_ptr);
94 ROL::Ptr<SV> xquad2 = ROL::makePtr<SV>(xquad2_ptr);
96 ROL::Ptr<SV> xdwell0 = ROL::makePtr<SV>(xdwell0_ptr);
97 ROL::Ptr<SV> xdwell1 = ROL::makePtr<SV>(xdwell1_ptr);
98 ROL::Ptr<SV> xdwell2 = ROL::makePtr<SV>(xdwell2_ptr);
100 ROL::Ptr<V> lo = ROL::makePtr<SV>(l_ptr);
101 ROL::Ptr<V> up = ROL::makePtr<SV>(u_ptr);
103 for(uint i=0; i<
dim; ++i) {
104 RealT t =
static_cast<RealT>(i)/static_cast<RealT>(dim-1);
105 (*x_ptr)[i] = xmin*(1-t) + xmax*t;
111 ROL::ParameterList logList;
112 ROL::ParameterList quadList;
113 ROL::ParameterList dwellList;
115 logList.sublist(
"Barrier Function").set(
"Type",
"Logarithmic");
116 quadList.sublist(
"Barrier Function").set(
"Type",
"Quadratic");
117 dwellList.sublist(
"Barrier Function").set(
"Type",
"Double Well");
136 quadObj.
value(x,tol);
146 dwellObj.
value(x,tol);
156 *outStream << std::setw(14) <<
"x"
157 << std::setw(14) <<
"log"
158 << std::setw(14) <<
"D(log)"
159 << std::setw(14) <<
"D2(log)"
160 << std::setw(14) <<
"quad"
161 << std::setw(14) <<
"D(quad)"
162 << std::setw(14) <<
"D2(quad)"
163 << std::setw(14) <<
"dwell"
164 << std::setw(14) <<
"D(dwell)"
165 << std::setw(14) <<
"D2(dwell)"
167 *outStream << std::string(140,
'-') << std::endl;
169 for(uint i=0; i<
dim; ++i) {
170 *outStream << std::setw(14) << (*x_ptr)[i]
171 << std::setw(14) << (*xlog0_ptr)[i]
172 << std::setw(14) << (*xlog1_ptr)[i]
173 << std::setw(14) << (*xlog2_ptr)[i]
174 << std::setw(14) << (*xquad0_ptr)[i]
175 << std::setw(14) << (*xquad1_ptr)[i]
176 << std::setw(14) << (*xquad2_ptr)[i]
177 << std::setw(14) << (*xdwell0_ptr)[i]
178 << std::setw(14) << (*xdwell1_ptr)[i]
179 << std::setw(14) << (*xdwell2_ptr)[i]
187 *outStream <<
"\n\n";
188 *outStream <<
"Test of logarithmic penalty objective" << std::endl;
189 logObj.
checkGradient(x,v,
true,*outStream); *outStream << std::endl;
190 logObj.
checkHessVec(x,v,
true,*outStream); *outStream << std::endl;
195 *outStream <<
"\n\n";
196 *outStream <<
"Test of piecewise quadratic penalty objective" << std::endl;
197 quadObj.
checkGradient(x,v,
true,*outStream); *outStream << std::endl;
198 quadObj.
checkHessVec(x,v,
true,*outStream); *outStream << std::endl;
201 *outStream <<
"\n\n";
202 *outStream <<
"Test of double well penalty objective" << std::endl;
203 dwellObj.
checkGradient(x,v,
true,*outStream); *outStream << std::endl;
204 dwellObj.
checkHessVec(x,v,
true,*outStream); *outStream << std::endl;
211 catch (std::logic_error& err) {
212 *outStream << err.what() <<
"\n";
217 std::cout <<
"End Result: TEST FAILED\n";
219 std::cout <<
"End Result: TEST PASSED\n";
typename PV< Real >::size_type size_type
ROL::Ptr< Vector< Real > > getBarrierVector(void)
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.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
basic_nullstream< char, std::char_traits< char >> nullstream
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Provides the elementwise interface to apply upper and lower bound constraints.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
int main(int argc, char *argv[])
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Create a penalty objective from upper and lower bound vectors.