51 #include "ROL_ParameterList.hpp"
53 #include "Teuchos_GlobalMPISession.hpp"
54 #include "Teuchos_LAPACK.hpp"
77 return std::sqrt(this->
dot(r,r));
80 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
82 Real c = (((int)x.size()==this->
nx_) ? 4.0 : 2.0);
83 for (
unsigned i=0; i<x.size(); i++) {
85 ip += this->
dx_/6.0*(c*x[i] + x[i+1])*y[i];
87 else if ( i == x.size()-1 ) {
88 ip += this->
dx_/6.0*(x[i-1] + c*x[i])*y[i];
91 ip += this->
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
99 void update(std::vector<Real> &u,
const std::vector<Real> &s,
const Real alpha=1.0) {
100 for (
unsigned i=0; i<u.size(); i++) {
105 void scale(std::vector<Real> &u,
const Real alpha=0.0) {
106 for (
unsigned i=0; i<u.size(); i++) {
112 const std::vector<Real> &z) {
114 r.resize(this->
nx_,0.0);
115 for (
int i=0; i<this->
nx_; i++) {
118 r[i] = this->
nu_/this->
dx_*(2.0*u[i]-u[i+1]);
120 else if (i==this->nx_-1) {
121 r[i] = this->
nu_/this->
dx_*(2.0*u[i]-u[i-1]);
124 r[i] = this->
nu_/this->
dx_*(2.0*u[i]-u[i-1]-u[i+1]);
128 r[i] += u[i+1]*(u[i]+u[i+1])/6.0;
131 r[i] -= u[i-1]*(u[i-1]+u[i])/6.0;
134 r[i] -= this->
dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
136 r[i] -= this->
dx_*this->
f_;
139 r[0] -= this->
u0_*u[ 0]/6.0 + this->
u0_*this->
u0_/6.0 + this->
nu_*this->
u0_/this->
dx_;
140 r[this->nx_-1] += this->
u1_*u[this->nx_-1]/6.0 + this->
u1_*this->
u1_/6.0 - this->
nu_*this->
u1_/this->
dx_;
144 const std::vector<Real> &u) {
147 d.resize(this->
nx_,this->
nu_*2.0/this->
dx_);
149 dl.resize(this->
nx_-1,-this->
nu_/this->
dx_);
151 du.resize(this->
nx_-1,-this->
nu_/this->
dx_);
153 for (
int i=0; i<this->
nx_; i++) {
155 dl[i] += (-2.0*u[i]-u[i+1])/6.0;
160 du[i-1] += (u[i-1]+2.0*u[i])/6.0;
164 d[0] -= this->
u0_/6.0;
165 d[this->nx_-1] += this->
u1_/6.0;
168 void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
169 const std::vector<Real> &r,
const bool transpose =
false) {
170 u.assign(r.begin(),r.end());
172 Teuchos::LAPACK<int,Real> lp;
173 std::vector<Real> du2(this->
nx_-2,0.0);
174 std::vector<int> ipiv(this->
nx_,0);
178 lp.GTTRF(this->
nx_,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
183 lp.GTTRS(trans,this->
nx_,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
190 dx_ = 1.0/((Real)nx+1.0);
195 ROL::Ptr<std::vector<Real> > cp =
197 ROL::Ptr<const std::vector<Real> > up =
199 ROL::Ptr<const std::vector<Real> > zp =
206 ROL::Ptr<std::vector<Real> > jvp =
208 ROL::Ptr<const std::vector<Real> > vp =
210 ROL::Ptr<const std::vector<Real> > up =
212 ROL::Ptr<const std::vector<Real> > zp =
215 for (
int i = 0; i < this->
nx_; i++) {
216 (*jvp)[i] = this->
nu_/this->
dx_*2.0*(*vp)[i];
218 (*jvp)[i] += -this->
nu_/this->
dx_*(*vp)[i-1]
219 -(*up)[i-1]/6.0*(*vp)[i]
220 -((*up)[i]+2.0*(*up)[i-1])/6.0*(*vp)[i-1];
222 if ( i < this->nx_-1 ) {
223 (*jvp)[i] += -this->
nu_/this->
dx_*(*vp)[i+1]
224 +(*up)[i+1]/6.0*(*vp)[i]
225 +((*up)[i]+2.0*(*up)[i+1])/6.0*(*vp)[i+1];
228 (*jvp)[0] -= this->
u0_/6.0*(*vp)[0];
229 (*jvp)[this->nx_-1] += this->
u1_/6.0*(*vp)[this->nx_-1];
234 ROL::Ptr<std::vector<Real> > jvp =
236 ROL::Ptr<const std::vector<Real> > vp =
238 ROL::Ptr<const std::vector<Real> > up =
240 ROL::Ptr<const std::vector<Real> > zp =
242 for (
int i=0; i<this->
nx_; i++) {
244 (*jvp)[i] = -this->
dx_/6.0*((*vp)[i]+4.0*(*vp)[i+1]+(*vp)[i+2]);
250 ROL::Ptr<std::vector<Real> > ijvp =
252 ROL::Ptr<const std::vector<Real> > vp =
254 ROL::Ptr<const std::vector<Real> > up =
256 ROL::Ptr<const std::vector<Real> > zp =
259 std::vector<Real> d(this->
nx_,0.0);
260 std::vector<Real> dl(this->
nx_-1,0.0);
261 std::vector<Real> du(this->
nx_-1,0.0);
269 ROL::Ptr<std::vector<Real> > jvp =
271 ROL::Ptr<const std::vector<Real> > vp =
273 ROL::Ptr<const std::vector<Real> > up =
275 ROL::Ptr<const std::vector<Real> > zp =
278 for (
int i = 0; i < this->
nx_; i++) {
279 (*jvp)[i] = this->
nu_/this->
dx_*2.0*(*vp)[i];
281 (*jvp)[i] += -this->
nu_/this->
dx_*(*vp)[i-1]
282 -(*up)[i-1]/6.0*(*vp)[i]
283 +((*up)[i-1]+2.0*(*up)[i])/6.0*(*vp)[i-1];
285 if ( i < this->nx_-1 ) {
286 (*jvp)[i] += -this->
nu_/this->
dx_*(*vp)[i+1]
287 +(*up)[i+1]/6.0*(*vp)[i]
288 -((*up)[i+1]+2.0*(*up)[i])/6.0*(*vp)[i+1];
291 (*jvp)[0] -= this->
u0_/6.0*(*vp)[0];
292 (*jvp)[this->nx_-1] += this->
u1_/6.0*(*vp)[this->nx_-1];
297 ROL::Ptr<std::vector<Real> > jvp =
299 ROL::Ptr<const std::vector<Real> > vp =
301 ROL::Ptr<const std::vector<Real> > up =
303 ROL::Ptr<const std::vector<Real> > zp =
305 for (
int i=0; i<this->
nx_+2; i++) {
307 (*jvp)[i] = -this->
dx_/6.0*(*vp)[i];
310 (*jvp)[i] = -this->
dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i]);
312 else if ( i == this->nx_ ) {
313 (*jvp)[i] = -this->
dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i-2]);
315 else if ( i == this->nx_+1 ) {
316 (*jvp)[i] = -this->
dx_/6.0*(*vp)[i-2];
319 (*jvp)[i] = -this->
dx_/6.0*((*vp)[i-2]+4.0*(*vp)[i-1]+(*vp)[i]);
326 ROL::Ptr<std::vector<Real> > iajvp =
328 ROL::Ptr<const std::vector<Real> > vp =
330 ROL::Ptr<const std::vector<Real> > up =
333 std::vector<Real> d(this->
nx_,0.0);
334 std::vector<Real> du(this->
nx_-1,0.0);
335 std::vector<Real> dl(this->
nx_-1,0.0);
343 ROL::Ptr<std::vector<Real> > ahwvp =
345 ROL::Ptr<const std::vector<Real> > wp =
347 ROL::Ptr<const std::vector<Real> > vp =
349 ROL::Ptr<const std::vector<Real> > up =
351 ROL::Ptr<const std::vector<Real> > zp =
353 for (
int i=0; i<this->
nx_; i++) {
357 (*ahwvp)[i] += ((*wp)[i]*(*vp)[i+1] - (*wp)[i+1]*(2.0*(*vp)[i]+(*vp)[i+1]))/6.0;
360 (*ahwvp)[i] += ((*wp)[i-1]*((*vp)[i-1]+2.0*(*vp)[i]) - (*wp)[i]*(*vp)[i-1])/6.0;
397 case 1: val = ((x<0.5) ? 1.0 : 0.0);
break;
398 case 2: val = 1.0;
break;
399 case 3: val = std::abs(std::sin(8.0*M_PI*x));
break;
400 case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5));
break;
405 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
407 Real c = (((int)x.size()==this->
nx_) ? 4.0 : 2.0);
408 for (
unsigned i=0; i<x.size(); i++) {
410 ip += this->
dx_/6.0*(c*x[i] + x[i+1])*y[i];
412 else if ( i == x.size()-1 ) {
413 ip += this->
dx_/6.0*(x[i-1] + c*x[i])*y[i];
416 ip += this->
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
422 void apply_mass(std::vector<Real> &Mu,
const std::vector<Real> &u ) {
423 Mu.resize(u.size(),0.0);
424 Real c = (((int)u.size()==this->
nx_) ? 4.0 : 2.0);
425 for (
unsigned i=0; i<u.size(); i++) {
427 Mu[i] = this->
dx_/6.0*(c*u[i] + u[i+1]);
429 else if ( i == u.size()-1 ) {
430 Mu[i] = this->
dx_/6.0*(u[i-1] + c*u[i]);
433 Mu[i] = this->
dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
444 dx_ = 1.0/((Real)nx+1.0);
448 ROL::Ptr<const std::vector<Real> > up =
450 ROL::Ptr<const std::vector<Real> > zp =
453 Real res1 = 0.0, res2 = 0.0, res3 = 0.0;
454 Real valu = 0.0, valz = this->
dot(*zp,*zp);
455 for (
int i=0; i<this->
nx_; i++) {
459 valu += this->
dx_/6.0*(4.0*res1 + res2)*res1;
461 else if ( i == this->nx_-1 ) {
464 valu += this->dx_/6.0*(res1 + 4.0*res2)*res2;
470 valu += this->dx_/6.0*(res1 + 4.0*res2 + res3)*res2;
473 return 0.5*(valu + this->
alpha_*valz);
478 ROL::Ptr<std::vector<Real> > gup = ROL::constPtrCast<std::vector<Real> >(
481 ROL::Ptr<const std::vector<Real> > up =
483 ROL::Ptr<const std::vector<Real> > zp =
486 std::vector<Real> diff(this->
nx_,0.0);
487 for (
int i=0; i<this->
nx_; i++) {
495 ROL::Ptr<std::vector<Real> > gzp = ROL::constPtrCast<std::vector<Real> >(
498 ROL::Ptr<const std::vector<Real> > up =
500 ROL::Ptr<const std::vector<Real> > zp =
503 for (
int i=0; i<this->
nx_+2; i++) {
505 (*gzp)[i] = this->
alpha_*this->
dx_/6.0*(2.0*(*zp)[i]+(*zp)[i+1]);
507 else if (i==this->nx_+1) {
508 (*gzp)[i] = this->
alpha_*this->
dx_/6.0*(2.0*(*zp)[i]+(*zp)[i-1]);
511 (*gzp)[i] = this->
alpha_*this->
dx_/6.0*((*zp)[i-1]+4.0*(*zp)[i]+(*zp)[i+1]);
518 ROL::Ptr<std::vector<Real> > hvup =
521 ROL::Ptr<const std::vector<Real> > vup =
539 ROL::Ptr<std::vector<Real> > hvzp =
542 ROL::Ptr<const std::vector<Real> > vzp =
545 for (
int i=0; i<this->
nx_+2; i++) {
547 (*hvzp)[i] = this->
alpha_*this->
dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i+1]);
549 else if (i==this->nx_+1) {
550 (*hvzp)[i] = this->
alpha_*this->
dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i-1]);
553 (*hvzp)[i] = this->
alpha_*this->
dx_/6.0*((*vzp)[i-1]+4.0*(*vzp)[i]+(*vzp)[i+1]);
Provides the interface to evaluate simulation-based objective functions.
void applyAdjointJacobian_1(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
Real evaluate_target(Real x)
void applyAdjointHessian_12(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
Contains definitions of custom data types in ROL.
void applyAdjointHessian_11(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
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 ...
virtual void zero()
Set to zero vector.
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 hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.
void gradient_1(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
void applyAdjointHessian_21(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Objective_BurgersControl(Real alpha=1.e-4, int nx=128)
void applyJacobian_1(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Constraint_BurgersControl(int nx=128, Real nu=1.e-2, Real u0=1.0, Real u1=0.0, Real f=0.0)
void gradient_2(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
void applyInverseJacobian_1(ROL::Vector< Real > &ijv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse partial constraint Jacobian at , , to the vector .
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
void scale(std::vector< Real > &u, const Real alpha=0.0)
void hessVec_11(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
void applyJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
void applyAdjointHessian_22(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
Real compute_norm(const std::vector< Real > &r)
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
void applyAdjointJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
Defines the constraint operator interface for simulation-based optimization.
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
void linear_solve(std::vector< Real > &u, std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &r, const bool transpose=false)