61 #include "Teuchos_GlobalMPISession.hpp"
90 Teuchos::SerialDenseMatrix<int, Real>
H_;
108 hu_ = 1.0/((Real)
nu_+1.0);
112 void apply_mass(std::vector<Real> &Mz,
const std::vector<Real> &z ) {
116 Mz[i] =
hu_/6.0*(2.0*z[i] + z[i+1]);
118 else if ( i == nu_-1 ) {
119 Mz[i] =
hu_/6.0*(z[i-1] + 2.0*z[i]);
122 Mz[i] =
hu_/6.0*(z[i-1] + 4.0*z[i] + z[i+1]);
128 return (x <= 0.5) ? 1.0 : 0.0;
132 const std::vector<Real> &d,
const std::vector<Real> &u,
133 bool addBC =
true ) {
136 for (
uint i = 0; i <
nu_; i++) {
138 Bd[i] = 1.0/
hu_*( u[i]*d[i] + (u[i]-u[i+1])*d[i+1] );
140 else if ( i == nu_-1 ) {
141 Bd[i] = 1.0/
hu_*( (u[i]-u[i-1])*d[i] + u[i]*d[i+1] );
144 Bd[i] = 1.0/
hu_*( (u[i]-u[i-1])*d[i] + (u[i]-u[i+1])*d[i+1] );
154 const std::vector<Real> &d,
const std::vector<Real> &u,
155 bool addBC =
true ) {
158 for (
uint i = 0; i <
nz_; i++) {
160 Bd[i] = 1.0/
hu_*u[i]*d[i];
162 else if ( i == nz_-1 ) {
163 Bd[i] = 1.0/
hu_*u[i-1]*d[i-1];
166 Bd[i] = 1.0/
hu_*( (u[i]-u[i-1])*(d[i]-d[i-1]) );
178 std::vector<Real> d(
nu_,1.0);
179 std::vector<Real> o(
nu_-1,1.0);
180 for (
uint i = 0; i <
nu_; i++ ) {
181 d[i] = (z[i] + z[i+1])/
hu_;
192 Teuchos::LAPACK<int,Real> lp;
196 lp.PTTRF(nu_,&d[0],&o[0],&info);
197 lp.PTTRS(nu_,nhrs,&d[0],&o[0],&u[0],ldb,&info);
204 for (
uint i = 0; i <
nu_; i++ ) {
205 d[i] = (z[i] + z[i+1])/
hu_;
212 for (
uint i = 0; i <
nu_; i++) {
219 Teuchos::LAPACK<int,Real> lp;
223 lp.PTTRF(nu_,&d[0],&o[0],&info);
224 lp.PTTRS(nu_,nhrs,&d[0],&o[0],&p[0],ldb,&info);
228 const std::vector<Real> &u,
const std::vector<Real> &z) {
232 for (
uint i = 0; i <
nu_; i++ ) {
233 d[i] = (z[i] + z[i+1])/
hu_;
243 Teuchos::LAPACK<int,Real> lp;
247 lp.PTTRF(nu_,&d[0],&o[0],&info);
248 lp.PTTRS(nu_,nhrs,&d[0],&o[0],&w[0],ldb,&info);
252 const std::vector<Real> &v,
const std::vector<Real> &p,
253 const std::vector<Real> &u,
const std::vector<Real> &z) {
257 for (
uint i = 0; i <
nu_; i++ ) {
258 d[i] = (z[i] + z[i+1])/
hu_;
267 std::vector<Real> res(nu_,0.0);
269 for (
uint i = 0; i <
nu_; i++) {
273 Teuchos::LAPACK<int,Real> lp;
277 lp.PTTRF(nu_,&d[0],&o[0],&info);
278 lp.PTTRS(nu_,nhrs,&d[0],&o[0],&q[0],ldb,&info);
286 Real tol = std::sqrt(ROL::ROL_EPSILON<Real>());
288 ROL::Ptr<V> e = z.
clone();
289 ROL::Ptr<V> h = z.
clone();
290 for (
uint i = 0; i <
nz_; i++ ) {
293 for (
uint j = 0; j <
nz_; j++ ) {
295 (
H_)(j,i) = e->dot(*h);
298 std::vector<vector> eigenvals = ROL::computeEigenvalues<Real>(
H_);
299 std::sort((eigenvals[0]).begin(), (eigenvals[0]).end());
300 Real inertia = (eigenvals[0])[0];
301 Real correction = 0.0;
302 if ( inertia <= 0.0 ) {
303 correction = (1.0+std::sqrt(ROL::ROL_EPSILON<Real>()))*std::abs(inertia);
304 if ( inertia == 0.0 ) {
306 while ( eigenvals[0][cnt] == 0.0 ) {
309 correction = std::sqrt(ROL::ROL_EPSILON<Real>())*eigenvals[0][cnt];
310 if ( cnt == nz_-1 ) {
314 for (
uint i = 0; i <
nz_; i++ ) {
315 (
H_)(i,i) += correction;
325 ROL::Ptr<const vector> zp =
getVector(z);
344 res =
hu_/6.0*(2.0*res1 + res2)*res1;
346 else if ( i == nu_-1 ) {
349 res = hu_/6.0*(res1 + 2.0*res2)*res2;
355 res = hu_/6.0*(res1 + 4.0*res2 + res3)*res2;
366 ROL::Ptr<const vector> zp =
getVector(z);
380 for (
uint i = 0; i <
nz_; i++ ) {
406 ROL::Ptr<const vector> vp =
getVector(v);
407 ROL::Ptr<const vector> zp =
getVector(z);
429 std::vector<Real> tmp(
nz_,0.0);
443 using ROL::constPtrCast;
448 ROL::Ptr<vector> vp = ROL::constPtrCast<vector>(
getVector(v));
450 Teuchos::SerialDenseVector<int, Real> hv_teuchos(Teuchos::View, &((*hvp)[0]), static_cast<int>(
nz_));
451 Teuchos::SerialDenseVector<int, Real> v_teuchos(Teuchos::View, &(( *vp)[0]), static_cast<int>(
nz_));
452 hv_teuchos.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0,
H_, v_teuchos, 0.0);
461 int main(
int argc,
char *argv[]) {
463 typedef std::vector<RealT> vector;
471 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
474 int iprint = argc - 1;
475 ROL::Ptr<std::ostream> outStream;
478 outStream = ROL::makePtrFromRef(std::cout);
480 outStream = ROL::makePtrFromRef(bhs);
493 ROL::Ptr<vector> x_ptr = ROL::makePtr<vector>(
dim, 0.0);
494 ROL::Ptr<vector> y_ptr = ROL::makePtr<vector>(
dim, 0.0);
497 for (uint i=0; i<
dim; i++) {
498 (*x_ptr)[i] = (
RealT)rand()/(
RealT)RAND_MAX + 1.e2;
499 (*y_ptr)[i] = (
RealT)rand()/(
RealT)RAND_MAX + 1.e2;
508 ROL::Ptr<vector> l_ptr = ROL::makePtr<vector>(
dim,1.0);
509 ROL::Ptr<vector> u_ptr = ROL::makePtr<vector>(
dim,10.0);
511 ROL::Ptr<V> lo = ROL::makePtr<SV>(l_ptr);
512 ROL::Ptr<V> up = ROL::makePtr<SV>(u_ptr);
516 ROL::ParameterList parlist;
519 parlist.sublist(
"General").sublist(
"Krylov").set(
"Absolute Tolerance",1.e-8);
520 parlist.sublist(
"General").sublist(
"Krylov").set(
"Relative Tolerance",1.e-4);
521 parlist.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit",static_cast<int>(dim));
523 parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").set(
"Relative Step Tolerance",1.e-8);
524 parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").set(
"Relative Gradient Tolerance",1.e-6);
525 parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").set(
"Iteration Limit", 10);
526 parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").set(
"Dual Scaling",(alpha>0.0)?alpha:1.e-4);
527 parlist.sublist(
"General").sublist(
"Secant").set(
"Use as Hessian",
true);
529 parlist.sublist(
"Status Test").set(
"Gradient Tolerance",1.e-12);
530 parlist.sublist(
"Status Test").set(
"Step Tolerance",1.e-14);
531 parlist.sublist(
"Status Test").set(
"Iteration Limit",100);
533 ROL::Ptr<ROL::Step<RealT>> step;
534 ROL::Ptr<ROL::StatusTest<RealT>> status;
537 step = ROL::makePtr<ROL::PrimalDualActiveSetStep<RealT>>(parlist);
538 status = ROL::makePtr<ROL::StatusTest<RealT>>(parlist);
543 algo.run(x,obj,icon,
true,*outStream);
547 file.open(
"control_PDAS.txt");
548 for ( uint i = 0; i <
dim; i++ ) {
549 file << (*x_ptr)[i] <<
"\n";
555 parlist.sublist(
"General").sublist(
"Secant").set(
"Use as Hessian",
false);
556 parlist.sublist(
"Step").sublist(
"Trust Region").set(
"Subproblem Solver",
"Truncated CG");
557 parlist.sublist(
"Step").sublist(
"Trust Region").set(
"Initial Radius", 1e3);
558 parlist.sublist(
"Step").sublist(
"Trust Region").set(
"Maximum Radius", 1e8);
559 step = ROL::makePtr<ROL::TrustRegionStep<RealT>>(parlist);
560 status = ROL::makePtr<ROL::StatusTest<RealT>>(parlist);
565 algo_tr.run(y,obj,icon,
true,*outStream);
567 std::ofstream file_tr;
568 file_tr.open(
"control_TR.txt");
569 for ( uint i = 0; i <
dim; i++ ) {
570 file_tr << (*y_ptr)[i] <<
"\n";
574 std::vector<RealT> u;
576 std::ofstream file_u;
577 file_u.open(
"state.txt");
578 for ( uint i = 0; i < (dim-1); i++ ) {
579 file_u << u[i] <<
"\n";
583 ROL::Ptr<V> diff = x.clone();
586 RealT error = diff->norm()/std::sqrt((
RealT)dim-1.0);
587 std::cout <<
"\nError between PDAS solution and TR solution is " << error <<
"\n";
588 errorFlag = ((error > 1.e2*std::sqrt(ROL::ROL_EPSILON<RealT>())) ? 1 : 0);
591 catch (std::logic_error& err) {
592 *outStream << err.what() <<
"\n";
597 std::cout <<
"End Result: TEST FAILED\n";
599 std::cout <<
"End Result: TEST PASSED\n";
Provides the interface to evaluate objective functions.
typename PV< Real >::size_type size_type
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
std::vector< Real > vector
virtual ROL::Ptr< Vector > basis(const int i) const
Return i-th basis vector.
void solve_state_equation(std::vector< Real > &u, const std::vector< Real > &z)
void apply_transposed_linearized_control_operator(std::vector< Real > &Bd, const std::vector< Real > &z, const std::vector< Real > &d, const std::vector< Real > &u, bool addBC=true)
Contains definitions of custom data types in ROL.
Contains definitions for helper functions in ROL.
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...
void solve_adjoint_sensitivity_equation(std::vector< Real > &q, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &p, const std::vector< Real > &u, const std::vector< Real > &z)
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.
void apply_mass(std::vector< Real > &Mz, const std::vector< Real > &z)
Real value(const ROL::Vector< Real > &z, Real &tol)
Compute value.
void update(const ROL::Vector< Real > &z, bool flag, int iter)
Update objective function.
void deactivateInertia(void)
Provides an interface to run optimization algorithms.
Provides the elementwise interface to apply upper and lower bound constraints.
void hessVec_true(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
ROL::StdVector< Real > SV
void solve_adjoint_equation(std::vector< Real > &p, const std::vector< Real > &u, const std::vector< Real > &z)
void solve_state_sensitivity_equation(std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
basic_nullstream< char, char_traits< char >> nullstream
void hessVec(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
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 apply_linearized_control_operator(std::vector< Real > &Bd, const std::vector< Real > &z, const std::vector< Real > &d, const std::vector< Real > &u, bool addBC=true)
void activateInertia(void)
void hessVec_inertia(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
Objective_PoissonInversion(int nz=32, Real alpha=1.e-4)
void gradient(ROL::Vector< Real > &g, const ROL::Vector< Real > &z, Real &tol)
Compute gradient.
Real evaluate_target(Real x)
ROL::Ptr< const vector > getVector(const V &x)
Teuchos::SerialDenseMatrix< int, Real > H_
ROL::Ptr< vector > getVector(V &x)