53 #include "ROL_ParameterList.hpp" 
   56 #include "Teuchos_GlobalMPISession.hpp" 
   71 int main(
int argc, 
char *argv[]) {
 
   73   Teuchos::GlobalMPISession mpiSession(&argc, &argv);
 
   75   int iprint     = argc - 1;
 
   76   ROL::Ptr<std::ostream> outStream;
 
   79     outStream = ROL::makePtrFromRef(std::cout);
 
   81     outStream = ROL::makePtrFromRef(bhs);
 
   99     ROL::Ptr<BurgersFEM<RealT> > fem
 
  100       = ROL::makePtr<BurgersFEM<RealT>>(nx,nu,nl,u0,u1,f,cH1,cL2);
 
  101     fem->test_inverse_mass(*outStream);
 
  102     fem->test_inverse_H1(*outStream);
 
  106     ROL::Ptr<std::vector<RealT> > ud_ptr
 
  107       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  108     ROL::Ptr<ROL::Vector<RealT> > ud
 
  109       = ROL::makePtr<L2VectorPrimal<RealT>>(ud_ptr,fem);
 
  114     bool useEChessian = 
true;
 
  120     std::vector<RealT> Ulo(nx, 0.), Uhi(nx, 1.);
 
  122     ROL::Ptr<ROL::BoundConstraint<RealT> > Ubnd
 
  123        = ROL::makePtr<H1BoundConstraint<RealT>>(Ulo,Uhi,fem);
 
  127     std::vector<RealT> Zlo(nx+2,0.), Zhi(nx+2,2.);
 
  128     ROL::Ptr<ROL::BoundConstraint<RealT> > Zbnd
 
  129       = ROL::makePtr<L2BoundConstraint<RealT>>(Zlo,Zhi,fem);
 
  138     ROL::Ptr<std::vector<RealT> > z_ptr
 
  139       = ROL::makePtr<std::vector<RealT>>(nx+2, 0.);
 
  140     ROL::Ptr<std::vector<RealT> > zrand_ptr
 
  141       = ROL::makePtr<std::vector<RealT>>(nx+2, 1.);
 
  142     ROL::Ptr<std::vector<RealT> > gz_ptr
 
  143       = ROL::makePtr<std::vector<RealT>>(nx+2, 1.);
 
  144     ROL::Ptr<std::vector<RealT> > yz_ptr
 
  145       = ROL::makePtr<std::vector<RealT>>(nx+2, 1.);
 
  146     for (
int i=0; i<nx+2; i++) {
 
  147       (*zrand_ptr)[i] = 10.*(
RealT)rand()/(
RealT)RAND_MAX-5.;
 
  148       (*yz_ptr)[i] = 10.*(
RealT)rand()/(
RealT)RAND_MAX-5.;
 
  150     ROL::Ptr<ROL::Vector<RealT> > zp
 
  151       = ROL::makePtr<PrimalControlVector>(z_ptr,fem);
 
  152     ROL::Ptr<ROL::Vector<RealT> > zrandp
 
  153       = ROL::makePtr<PrimalControlVector>(zrand_ptr,fem);
 
  154     ROL::Ptr<ROL::Vector<RealT> > gzp
 
  155       = ROL::makePtr<DualControlVector>(gz_ptr,fem);
 
  156     ROL::Ptr<ROL::Vector<RealT> > yzp
 
  157       = ROL::makePtr<PrimalControlVector>(yz_ptr,fem);
 
  159     ROL::Ptr<std::vector<RealT> > u_ptr
 
  160       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  161     ROL::Ptr<std::vector<RealT> > gu_ptr
 
  162       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  163     ROL::Ptr<std::vector<RealT> > yu_ptr
 
  164       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  165     for (
int i=0; i<nx; i++) {
 
  166       (*yu_ptr)[i] = 10.*(
RealT)rand()/(
RealT)RAND_MAX-5.;
 
  168     ROL::Ptr<ROL::Vector<RealT> > up
 
  169       = ROL::makePtr<PrimalStateVector>(u_ptr,fem);
 
  170     ROL::Ptr<ROL::Vector<RealT> > gup
 
  171       = ROL::makePtr<DualStateVector>(gu_ptr,fem);
 
  172     ROL::Ptr<ROL::Vector<RealT> > yup
 
  173       = ROL::makePtr<PrimalStateVector>(yu_ptr,fem);
 
  175     ROL::Ptr<std::vector<RealT> > c_ptr
 
  176       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  177     ROL::Ptr<std::vector<RealT> > l_ptr
 
  178       = ROL::makePtr<std::vector<RealT>>(nx, 1.);
 
  179     for (
int i=0; i<nx; i++) {
 
  189     std::string filename = 
"input.xml";
 
  190     auto parlist = ROL::getParametersFromXmlFile( filename );
 
  197     obj.checkGradient(x,g,y,
true,*outStream);
 
  198     obj.checkHessVec(x,g,y,
true,*outStream);
 
  210     ROL::Ptr<ROL::Objective<RealT> > obj_ptr = ROL::makePtrFromRef(obj);
 
  211     ROL::Ptr<ROL::Constraint<RealT> > con_ptr = ROL::makePtrFromRef(con);
 
  212     ROL::Ptr<ROL::BoundConstraint<RealT> > bnd_ptr = ROL::makePtrFromRef(bnd);
 
  223     ROL::Ptr<ROL::Step<RealT>>
 
  224       stepMY = ROL::makePtr<ROL::MoreauYosidaPenaltyStep<RealT>>(*parlist);
 
  225     ROL::Ptr<ROL::StatusTest<RealT>>
 
  226       statusMY = ROL::makePtr<ROL::ConstraintStatusTest<RealT>>(*parlist);
 
  229     RealT zerotol = std::sqrt(ROL::ROL_EPSILON<RealT>());
 
  230     con.
solve(c,*up,*zp,zerotol);
 
  231     obj.gradient_1(*gup,*up,*zp,zerotol);
 
  234     gup->zero(); c.
zero();
 
  235     algoMY.
run(x, g, l, c, myPen, con, bnd, 
true, *outStream);
 
  236     ROL::Ptr<ROL::Vector<RealT> > xMY = x.
clone();
 
  239     ROL::Ptr<ROL::Step<RealT>>
 
  240       stepAL = ROL::makePtr<ROL::AugmentedLagrangianStep<RealT>>(*parlist);
 
  241     ROL::Ptr<ROL::StatusTest<RealT>>
 
  242       statusAL = ROL::makePtr<ROL::ConstraintStatusTest<RealT>>(*parlist);
 
  245     con.
solve(c,*up,*zp,zerotol);
 
  246     obj.gradient_1(*gup,*up,*zp,zerotol);
 
  249     gup->zero(); c.
zero();
 
  250     algoAL.
run(x, g, l, c, myAugLag, con, bnd, 
true, *outStream);
 
  252     ROL::Ptr<ROL::Vector<RealT> > err = x.
clone();
 
  253     err->set(x); err->axpy(-1.,*xMY);
 
  254     errorFlag += ((err->norm() > 1.e-7*x.
norm()) ? 1 : 0);
 
  256   catch (std::logic_error& err) {
 
  257     *outStream << err.what() << 
"\n";
 
  262     std::cout << 
"End Result: TEST FAILED\n";
 
  264     std::cout << 
"End Result: TEST PASSED\n";
 
Provides the interface to evaluate the augmented Lagrangian. 
virtual Real checkSolve(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
Defines the linear algebra or vector space interface for simulation-based optimization. 
void solve(ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve  for . 
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout, bool printVectors=false, std::ostream &vectorStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface. 
virtual Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface. 
L2VectorPrimal< RealT > PrimalControlVector
void applyInverseAdjointJacobian_1(ROL::Vector< Real > &iajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
Real norm() const 
Returns  where . 
virtual void zero()
Set to zero vector. 
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. 
ROL::Ptr< Vector< Real > > clone() const 
Clone to make a new (uninitialized) vector. 
virtual Real checkInverseJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
H1VectorDual< RealT > DualStateVector
L2VectorDual< RealT > DualControlVector
Provides an interface to run optimization algorithms. 
Provides the interface to evaluate the Moreau-Yosida penalty function. 
H1VectorDual< RealT > PrimalConstraintVector
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
int main(int argc, char *argv[])
Provides definitions of equality constraint and objective for example_04. 
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 Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface. 
H1VectorPrimal< RealT > PrimalStateVector
virtual Real checkInverseAdjointJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
H1VectorPrimal< RealT > DualConstraintVector