56 #include "Teuchos_oblackholestream.hpp" 57 #include "Teuchos_GlobalMPISession.hpp" 58 #include "Teuchos_XMLParameterListHelpers.hpp" 64 int main(
int argc,
char *argv[]) {
66 using namespace Teuchos;
68 typedef std::vector<RealT> vector;
72 GlobalMPISession mpiSession(&argc, &argv);
75 int iprint = argc - 1;
76 Teuchos::RCP<std::ostream> outStream;
77 Teuchos::oblackholestream bhs;
79 outStream = Teuchos::rcp(&std::cout,
false);
81 outStream = Teuchos::rcp(&bhs,
false);
91 RCP<ParameterList> parlist = rcp(
new ParameterList());
92 std::string paramfile =
"parameters.xml";
93 updateParametersFromXmlFile(paramfile,parlist.ptr());
99 RCP<vector> x_rcp = rcp(
new vector(dim, 0.0) );
102 RCP<vector> k_rcp = rcp(
new vector(dim, 0.0) );
105 RCP<vector> xtest_rcp = rcp(
new vector(dim, 0.0) );
106 RCP<vector> d_rcp = rcp(
new vector(dim, 0.0) );
107 RCP<vector> v_rcp = rcp(
new vector(dim, 0.0) );
108 RCP<vector> hv_rcp = rcp(
new vector(dim, 0.0) );
109 RCP<vector> ihhv_rcp = rcp(
new vector(dim, 0.0) );
111 RealT left = -1e0, right = 1e0;
112 for (
int i=0; i<dim; i++) {
116 (*xtest_rcp)[i] = ( (
RealT)rand() / (
RealT)RAND_MAX ) * (right - left) + left;
117 (*d_rcp)[i] = ( (
RealT)rand() / (
RealT)RAND_MAX ) * (right - left) + left;
118 (*v_rcp)[i] = ( (
RealT)rand() / (
RealT)RAND_MAX ) * (right - left) + left;
121 RCP<V> k = rcp(
new SV(k_rcp) );
133 obj.
checkGradient(xtest, d,
true, *outStream); *outStream <<
"\n";
134 obj.
checkHessVec(xtest, v,
true, *outStream); *outStream <<
"\n";
135 obj.
checkHessSym(xtest, d, v,
true, *outStream); *outStream <<
"\n";
142 *outStream <<
"Checking inverse Hessian" << std::endl;
143 *outStream <<
"||H^{-1}Hv-v|| = " << ihhv.norm() << std::endl;
147 algo.run(x, obj,
true, *outStream);
150 RCP<vector> xtrue_rcp = rcp(
new vector(dim, 0.0) );
156 RealT abserr = x.norm();
157 *outStream << std::scientific <<
"\n Absolute Error: " << abserr << std::endl;
162 catch (std::logic_error err) {
163 *outStream << err.what() <<
"\n";
168 std::cout <<
"End Result: TEST FAILED\n";
170 std::cout <<
"End Result: TEST PASSED\n";
void invHessVec(Vector< Real > &ihv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Defines the linear algebra or vector space interface.
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.
Provides the std::vector implementation of the ROL::Vector interface.
Provides an interface to run optimization algorithms.
Contains definitions for the Zakharov function as evaluated using only the ROL::Vector interface...
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.
virtual std::vector< Real > checkHessSym(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
Hessian symmetry check.
static const double ROL_EPSILON
Platform-dependent machine epsilon.