25#include "ROL_ParameterList.hpp"
28#include "ROL_GlobalMPISession.hpp"
52 return std::sqrt(
dot(r,r));
55 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
57 Real c = ((x.size()==
nx_) ? 4.0 : 2.0);
58 for (
unsigned i = 0; i < x.size(); i++) {
60 ip +=
dx_/6.0*(c*x[i] + x[i+1])*y[i];
62 else if ( i == x.size()-1 ) {
63 ip +=
dx_/6.0*(x[i-1] + c*x[i])*y[i];
66 ip +=
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
74 void update(std::vector<Real> &u,
const std::vector<Real> &s,
const Real alpha=1.0) {
75 for (
unsigned i = 0; i < u.size(); i++) {
80 void scale(std::vector<Real> &u,
const Real alpha=0.0) {
81 for (
unsigned i = 0; i < u.size(); i++) {
86 void compute_residual(std::vector<Real> &r,
const std::vector<Real> &uold,
const std::vector<Real> &zold,
87 const std::vector<Real> &unew,
const std::vector<Real> &znew) {
90 for (
unsigned n = 0; n <
nx_; n++) {
104 r[n] -= 0.5*
dt_*unew[n-1]*(unew[n-1]+unew[n])/6.0;
105 r[n] -= 0.5*
dt_*uold[n-1]*(uold[n-1]+uold[n])/6.0;
108 r[n] += 0.5*
dt_*unew[n+1]*(unew[n]+unew[n+1])/6.0;
109 r[n] += 0.5*
dt_*uold[n+1]*(uold[n]+uold[n+1])/6.0;
112 r[n] -= 0.5*
dt_*
dx_/6.0*(znew[n]+4.0*znew[n+1]+znew[n+2]);
113 r[n] -= 0.5*
dt_*
dx_/6.0*(zold[n]+4.0*zold[n+1]+zold[n+2]);
123 const std::vector<Real> &u) {
132 for (
unsigned n = 0; n <
nx_; n++) {
134 dl[n] += 0.5*
dt_*(-2.0*u[n]-u[n+1])/6.0;
135 d[n] += 0.5*
dt_*u[n+1]/6.0;
138 d[n] -= 0.5*
dt_*u[n-1]/6.0;
139 du[n-1] += 0.5*
dt_*(u[n-1]+2.0*u[n])/6.0;
148 bool adjoint =
false) {
152 for (
unsigned n = 0; n <
nx_; n++) {
157 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]);
160 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]);
166 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]);
169 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]);
173 jv[0] -= 0.5*
dt_*
u0_/6.0*v[0];
178 bool adjoint =
false) {
182 for (
unsigned n = 0; n <
nx_; n++) {
187 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]);
190 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]);
196 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]);
199 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]);
203 jv[0] -= 0.5*
dt_*
u0_/6.0*v[0];
207 void apply_pde_jacobian(std::vector<Real> &jv,
const std::vector<Real> &vold,
const std::vector<Real> &uold,
208 const std::vector<Real> &vnew,
const std::vector<Real> unew,
bool adjoint =
false) {
212 for (
unsigned n = 0; n <
nx_; n++) {
219 jv[n] -= 0.5*
dt_*(unew[n-1]/6.0*vnew[n]-(unew[n-1]+2.0*unew[n])/6.0*vnew[n-1]);
220 jv[n] -= 0.5*
dt_*(uold[n-1]/6.0*vold[n]-(uold[n-1]+2.0*uold[n])/6.0*vold[n-1]);
223 jv[n] -= 0.5*
dt_*(unew[n-1]/6.0*vnew[n]+(unew[n]+2.0*unew[n-1])/6.0*vnew[n-1]);
224 jv[n] -= 0.5*
dt_*(uold[n-1]/6.0*vold[n]+(uold[n]+2.0*uold[n-1])/6.0*vold[n-1]);
231 jv[n] += 0.5*
dt_*(unew[n+1]/6.0*vnew[n]-(unew[n+1]+2.0*unew[n])/6.0*vnew[n+1]);
232 jv[n] += 0.5*
dt_*(uold[n+1]/6.0*vold[n]-(uold[n+1]+2.0*uold[n])/6.0*vold[n+1]);
235 jv[n] += 0.5*
dt_*(unew[n+1]/6.0*vnew[n]+(unew[n]+2.0*unew[n+1])/6.0*vnew[n+1]);
236 jv[n] += 0.5*
dt_*(uold[n+1]/6.0*vold[n]+(uold[n]+2.0*uold[n+1])/6.0*vold[n+1]);
240 jv[0] -= 0.5*
dt_*
u0_/6.0*vnew[0];
241 jv[0] -= 0.5*
dt_*
u0_/6.0*vold[0];
246 void apply_pde_hessian(std::vector<Real> &hv,
const std::vector<Real> &wold,
const std::vector<Real> &vold,
247 const std::vector<Real> &wnew,
const std::vector<Real> &vnew) {
250 for (
unsigned n = 0; n <
nx_; n++) {
252 hv[n] += 0.5*
dt_*((wnew[n-1]*(vnew[n-1]+2.0*vnew[n]) - wnew[n]*vnew[n-1])/6.0);
253 hv[n] += 0.5*
dt_*((wold[n-1]*(vold[n-1]+2.0*vold[n]) - wold[n]*vold[n-1])/6.0);
256 hv[n] += 0.5*
dt_*((wnew[n]*vnew[n+1] - wnew[n+1]*(2.0*vnew[n]+vnew[n+1]))/6.0);
257 hv[n] += 0.5*
dt_*((wold[n]*vold[n+1] - wold[n+1]*(2.0*vold[n]+vold[n+1]))/6.0);
264 unsigned dim = ((adjoint ==
true) ?
nx_+2 :
nx_);
266 for (
unsigned n = 0; n <
dim; n++) {
269 jv[n] = -0.5*
dt_*
dx_/6.0*v[n];
272 jv[n] = -0.5*
dt_*
dx_/6.0*(4.0*v[n-1]+v[n]);
274 else if ( n ==
nx_ ) {
275 jv[n] = -0.5*
dt_*
dx_/6.0*(4.0*v[n-1]+v[n-2]);
277 else if ( n ==
nx_+1 ) {
278 jv[n] = -0.5*
dt_*
dx_/6.0*v[n-2];
281 jv[n] = -0.5*
dt_*
dx_/6.0*(v[n-2]+4.0*v[n-1]+v[n]);
285 jv[n] -= 0.5*
dt_*
dx_/6.0*(v[n]+4.0*v[n+1]+v[n+2]);
290 void run_newton(std::vector<Real> &u,
const std::vector<Real> &znew,
291 const std::vector<Real> &uold,
const std::vector<Real> &zold) {
295 std::vector<Real> r(
nx_,0.0);
299 Real rtol = 1.e2*ROL::ROL_EPSILON<Real>();
300 unsigned maxit = 500;
302 std::vector<Real> d(
nx_,0.0);
303 std::vector<Real> dl(
nx_-1,0.0);
304 std::vector<Real> du(
nx_-1,0.0);
306 Real alpha = 1.0, tmp = 0.0;
307 std::vector<Real> s(
nx_,0.0);
308 std::vector<Real> utmp(
nx_,0.0);
309 for (
unsigned i = 0; i < maxit; i++) {
318 utmp.assign(u.begin(),u.end());
322 while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON<Real>()) ) {
324 utmp.assign(u.begin(),u.end());
330 u.assign(utmp.begin(),utmp.end());
331 if ( rnorm < rtol ) {
338 const std::vector<Real> &dl,
const std::vector<Real> &d,
const std::vector<Real> &du,
339 const std::vector<Real> &r,
const bool transpose =
false) {
340 bool useLAPACK =
true;
342 u.assign(r.begin(),r.end());
344 std::vector<Real> Dl(dl);
345 std::vector<Real> Du(du);
346 std::vector<Real> D(d);
348 ROL::LAPACK<int,Real> lp;
349 std::vector<Real> Du2(
nx_-2,0.0);
350 std::vector<int> ipiv(
nx_,0);
354 lp.GTTRF(
nx_,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&info);
355 char trans = ((transpose ==
true) ?
'T' :
'N');
356 lp.GTTRS(trans,
nx_,nhrs,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&u[0],ldb,&info);
361 unsigned maxit = 100;
362 Real rtol = std::min(1.e-12,1.e-4*std::sqrt(
dot(r,r)));
365 Real rnorm = 10.0*rtol;
366 for (
unsigned i = 0; i < maxit; i++) {
367 for (
unsigned n = 0; n <
nx_; n++) {
370 u[n] -= ((transpose ==
false) ? du[n] : dl[n])*u[n+1]/d[n];
373 u[n] -= ((transpose ==
false) ? dl[n-1] : du[n-1])*u[n-1]/d[n];
378 for (
unsigned n = 0; n <
nx_; n++) {
379 resid = r[n] - d[n]*u[n];
381 resid -= ((transpose ==
false) ? du[n] : dl[n])*u[n+1];
384 resid -= ((transpose ==
false) ? dl[n-1] : du[n-1])*u[n-1];
386 rnorm += resid*resid;
388 rnorm = std::sqrt(rnorm);
389 if ( rnorm < rtol ) {
400 Real nu = 1.e-2, Real u0 = 0.0, Real u1 = 0.0, Real f = 0.0)
402 dx_ = 1.0/((Real)nx+1.0);
407 for (
unsigned n = 0; n <
nx_; n++) {
409 u_init_[n] = ((x <= 0.5) ? 1.0 : 0.0);
414 ROL::Ptr<std::vector<Real> > cp =
416 ROL::Ptr<const std::vector<Real> > up =
418 ROL::Ptr<const std::vector<Real> > zp =
421 std::vector<Real> C(
nx_,0.0);
422 std::vector<Real> uold(
u_init_);
423 std::vector<Real> unew(
u_init_);
424 std::vector<Real> zold(
nx_+2,0.0);
425 std::vector<Real> znew(
nx_+2,0.0);
427 for (
unsigned n = 0; n <
nx_+2; n++) {
431 for (
unsigned t = 0; t <
nt_; t++) {
433 for (
unsigned n = 0; n <
nx_; n++) {
434 unew[n] = (*up)[t*
nx_+n];
437 for (
unsigned n = 0; n <
nx_+2; n++) {
438 znew[n] = (*zp)[(t+1)*(
nx_+2)+n];
443 for (
unsigned n = 0; n <
nx_; n++) {
444 (*cp)[t*
nx_+n] = C[n];
447 uold.assign(unew.begin(),unew.end());
448 zold.assign(znew.begin(),znew.end());
453 ROL::Ptr<std::vector<Real> > up =
455 up->assign(up->size(),z.
norm()/up->size());
456 ROL::Ptr<const std::vector<Real> > zp =
459 std::vector<Real> uold(
u_init_);
460 std::vector<Real> unew(
u_init_);
461 std::vector<Real> zold(
nx_+2,0.0);
462 std::vector<Real> znew(
nx_+2,0.0);
464 for (
unsigned n = 0; n <
nx_+2; n++) {
468 for (
unsigned t = 0; t <
nt_; t++) {
470 for (
unsigned n = 0; n <
nx_+2; n++) {
471 znew[n] = (*zp)[(t+1)*(
nx_+2)+n];
476 for (
unsigned n = 0; n <
nx_; n++) {
477 (*up)[t*
nx_+n] = unew[n];
480 uold.assign(unew.begin(),unew.end());
481 zold.assign(znew.begin(),znew.end());
488 ROL::Ptr<std::vector<Real> > jvp =
490 ROL::Ptr<const std::vector<Real> > vp =
492 ROL::Ptr<const std::vector<Real> > up =
494 std::vector<Real> J(
nx_,0.0);
495 std::vector<Real> vold(
nx_,0.0);
496 std::vector<Real> vnew(
nx_,0.0);
497 std::vector<Real> uold(
nx_,0.0);
498 std::vector<Real> unew(
nx_,0.0);
499 for (
unsigned t = 0; t <
nt_; t++) {
500 for (
unsigned n = 0; n <
nx_; n++) {
501 unew[n] = (*up)[t*
nx_+n];
502 vnew[n] = (*vp)[t*
nx_+n];
505 for (
unsigned n = 0; n <
nx_; n++) {
506 (*jvp)[t*
nx_+n] = J[n];
508 vold.assign(vnew.begin(),vnew.end());
509 uold.assign(unew.begin(),unew.end());
515 ROL::Ptr<std::vector<Real> > jvp =
517 ROL::Ptr<const std::vector<Real> > vp =
519 ROL::Ptr<const std::vector<Real> > zp =
521 std::vector<Real> vnew(
nx_+2,0.0);
522 std::vector<Real> vold(
nx_+2,0.0);
523 std::vector<Real> jnew(
nx_,0.0);
524 std::vector<Real> jold(
nx_,0.0);
525 for (
unsigned n = 0; n <
nx_+2; n++) {
529 for (
unsigned t = 0; t <
nt_; t++) {
530 for (
unsigned n = 0; n <
nx_+2; n++) {
531 vnew[n] = (*vp)[(t+1)*(
nx_+2)+n];
534 for (
unsigned n = 0; n <
nx_; n++) {
536 (*jvp)[t*
nx_+n] = jnew[n] + jold[n];
538 jold.assign(jnew.begin(),jnew.end());
544 ROL::Ptr<std::vector<Real> > ijvp =
546 ROL::Ptr<const std::vector<Real> > vp =
548 ROL::Ptr<const std::vector<Real> > up =
550 std::vector<Real> J(
nx_,0.0);
551 std::vector<Real> r(
nx_,0.0);
552 std::vector<Real> s(
nx_,0.0);
553 std::vector<Real> uold(
nx_,0.0);
554 std::vector<Real> unew(
nx_,0.0);
555 std::vector<Real> d(
nx_,0.0);
556 std::vector<Real> dl(
nx_-1,0.0);
557 std::vector<Real> du(
nx_-1,0.0);
558 for (
unsigned t = 0; t <
nt_; t++) {
560 for (
unsigned n = 0; n <
nx_; n++) {
561 r[n] = (*vp)[t*
nx_+n];
562 unew[n] = (*up)[t*
nx_+n];
571 for (
unsigned n = 0; n <
nx_; n++) {
572 (*ijvp)[t*
nx_+n] = s[n];
574 uold.assign(unew.begin(),unew.end());
580 ROL::Ptr<std::vector<Real> > jvp =
582 ROL::Ptr<const std::vector<Real> > vp =
584 ROL::Ptr<const std::vector<Real> > up =
586 std::vector<Real> J(
nx_,0.0);
587 std::vector<Real> vold(
nx_,0.0);
588 std::vector<Real> vnew(
nx_,0.0);
589 std::vector<Real> unew(
nx_,0.0);
590 for (
unsigned t =
nt_; t > 0; t--) {
591 for (
unsigned n = 0; n <
nx_; n++) {
592 unew[n] = (*up)[(t-1)*
nx_+n];
593 vnew[n] = (*vp)[(t-1)*
nx_+n];
596 for (
unsigned n = 0; n <
nx_; n++) {
597 (*jvp)[(t-1)*
nx_+n] = J[n];
599 vold.assign(vnew.begin(),vnew.end());
605 ROL::Ptr<std::vector<Real> > jvp =
607 ROL::Ptr<const std::vector<Real> > vp =
609 ROL::Ptr<const std::vector<Real> > zp =
611 std::vector<Real> vnew(
nx_,0.0);
612 std::vector<Real> vold(
nx_,0.0);
613 std::vector<Real> jnew(
nx_+2,0.0);
614 std::vector<Real> jold(
nx_+2,0.0);
615 for (
unsigned t =
nt_+1; t > 0; t--) {
616 for (
unsigned n = 0; n <
nx_; n++) {
618 vnew[n] = (*vp)[(t-2)*
nx_+n];
625 for (
unsigned n = 0; n <
nx_+2; n++) {
627 (*jvp)[(t-1)*(
nx_+2)+n] = jnew[n] + jold[n];
629 jold.assign(jnew.begin(),jnew.end());
635 ROL::Ptr<std::vector<Real> > ijvp =
637 ROL::Ptr<const std::vector<Real> > vp =
639 ROL::Ptr<const std::vector<Real> > up =
641 std::vector<Real> J(
nx_,0.0);
642 std::vector<Real> r(
nx_,0.0);
643 std::vector<Real> s(
nx_,0.0);
644 std::vector<Real> unew(
nx_,0.0);
645 std::vector<Real> d(
nx_,0.0);
646 std::vector<Real> dl(
nx_-1,0.0);
647 std::vector<Real> du(
nx_-1,0.0);
648 for (
unsigned t =
nt_; t > 0; t--) {
650 for (
unsigned n = 0; n <
nx_; n++) {
651 r[n] = (*vp)[(t-1)*
nx_+n];
652 unew[n] = (*up)[(t-1)*
nx_+n];
661 for (
unsigned n = 0; n <
nx_; n++) {
662 (*ijvp)[(t-1)*
nx_+n] = s[n];
669 ROL::Ptr<std::vector<Real> > hwvp =
671 ROL::Ptr<const std::vector<Real> > wp =
673 ROL::Ptr<const std::vector<Real> > vp =
675 std::vector<Real> snew(
nx_,0.0);
676 std::vector<Real> wnew(
nx_,0.0);
677 std::vector<Real> wold(
nx_,0.0);
678 std::vector<Real> vnew(
nx_,0.0);
679 for (
unsigned t =
nt_; t > 0; t--) {
680 for (
unsigned n = 0; n <
nx_; n++) {
681 vnew[n] = (*vp)[(t-1)*
nx_+n];
682 wnew[n] = (*wp)[(t-1)*
nx_+n];
685 for (
unsigned n = 0; n <
nx_; n++) {
686 (*hwvp)[(t-1)*
nx_+n] = snew[n];
688 wold.assign(wnew.begin(),wnew.end());
724 case 1: val = ((x<=0.5) ? 1.0 : 0.0);
break;
725 case 2: val = 1.0;
break;
726 case 3: val = std::abs(std::sin(8.0*M_PI*x));
break;
727 case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5));
break;
732 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
734 Real c = ((x.size()==
nx_) ? 4.0 : 2.0);
735 for (
unsigned i=0; i<x.size(); i++) {
737 ip +=
dx_/6.0*(c*x[i] + x[i+1])*y[i];
739 else if ( i == x.size()-1 ) {
740 ip +=
dx_/6.0*(x[i-1] + c*x[i])*y[i];
743 ip +=
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
749 void apply_mass(std::vector<Real> &Mu,
const std::vector<Real> &u ) {
750 Mu.resize(u.size(),0.0);
751 Real c = ((u.size()==
nx_) ? 4.0 : 2.0);
752 for (
unsigned i=0; i<u.size(); i++) {
754 Mu[i] =
dx_/6.0*(c*u[i] + u[i+1]);
756 else if ( i == u.size()-1 ) {
757 Mu[i] =
dx_/6.0*(u[i-1] + c*u[i]);
760 Mu[i] =
dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
772 dx_ = 1.0/((Real)nx+1.0);
777 ROL::Ptr<const std::vector<Real> > up =
779 ROL::Ptr<const std::vector<Real> > zp =
782 std::vector<Real> U(
nx_,0.0);
783 std::vector<Real> G(
nx_,0.0);
784 std::vector<Real> Z(
nx_+2,0.0);
785 for (
unsigned n = 0; n <
nx_+2; n++) {
790 for (
unsigned t = 0; t <
nt_; t++) {
792 for (
unsigned n = 0; n <
nx_; n++) {
796 val += 0.5*ss*
dot(U,U);
797 val -= 0.5*ss*
dot(G,G);
798 for (
unsigned n = 0; n <
nx_+2; n++) {
799 Z[n] = (*zp)[(t+1)*(
nx_+2)+n];
807 ROL::Ptr<std::vector<Real> > gp =
809 ROL::Ptr<const std::vector<Real> > up =
812 std::vector<Real> U(
nx_,0.0);
813 std::vector<Real> M(
nx_,0.0);
815 for (
unsigned t = 0; t <
nt_; t++) {
817 for (
unsigned n = 0; n <
nx_; n++) {
821 for (
unsigned n = 0; n <
nx_; n++) {
822 (*gp)[t*
nx_+n] = ss*M[n];
828 ROL::Ptr<std::vector<Real> > gp =
830 ROL::Ptr<const std::vector<Real> > zp =
833 std::vector<Real> Z(
nx_+2,0.0);
834 std::vector<Real> M(
nx_+2,0.0);
836 for (
unsigned t = 0; t <
nt_+1; t++) {
837 ss = ((t < nt_ && t > 0) ?
dt_ : 0.5*
dt_);
838 for (
unsigned n = 0; n <
nx_+2; n++) {
839 Z[n] = (*zp)[t*(
nx_+2)+n];
842 for (
unsigned n = 0; n <
nx_+2; n++) {
850 ROL::Ptr<std::vector<Real> > hvp =
852 ROL::Ptr<const std::vector<Real> > vp =
855 std::vector<Real>
V(
nx_,0.0);
856 std::vector<Real> M(
nx_,0.0);
858 for (
unsigned t = 0; t <
nt_; t++) {
860 for (
unsigned n = 0; n <
nx_; n++) {
861 V[n] = (*vp)[t*
nx_+n];
864 for (
unsigned n = 0; n <
nx_; n++) {
865 (*hvp)[t*
nx_+n] = ss*M[n];
882 ROL::Ptr<std::vector<Real> > hvp = ROL::constPtrCast<std::vector<Real> >(
884 ROL::Ptr<const std::vector<Real> > vp =
887 std::vector<Real>
V(
nx_+2,0.0);
888 std::vector<Real> M(
nx_+2,0.0);
890 for (
unsigned t = 0; t <
nt_+1; t++) {
891 ss = ((t < nt_ && t > 0) ?
dt_ : 0.5*
dt_);
892 for (
unsigned n = 0; n <
nx_+2; n++) {
893 V[n] = (*vp)[t*(
nx_+2)+n];
896 for (
unsigned n = 0; n <
nx_+2; n++) {
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Contains definitions of custom data types in ROL.
void scale(std::vector< Real > &u, const Real alpha=0.0)
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 ...
void compute_residual(std::vector< Real > &r, const std::vector< Real > &uold, const std::vector< Real > &zold, const std::vector< Real > &unew, const std::vector< Real > &znew)
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
void apply_pde_jacobian_new(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
void applyAdjointHessian_11(ROL::Vector< Real > &hwv, 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 apply_pde_jacobian(std::vector< Real > &jv, const std::vector< Real > &vold, const std::vector< Real > &uold, const std::vector< Real > &vnew, const std::vector< Real > unew, bool adjoint=false)
void applyInverseAdjointJacobian_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 of the adjoint of the partial constraint Jacobian at , , to the vector .
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...
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...
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
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 .
std::vector< Real > u_init_
Real compute_norm(const std::vector< Real > &r)
void apply_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, bool adjoint=false)
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 update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
void apply_pde_hessian(std::vector< Real > &hv, const std::vector< Real > &wold, const std::vector< Real > &vold, const std::vector< Real > &wnew, const std::vector< Real > &vnew)
void run_newton(std::vector< Real > &u, const std::vector< Real > &znew, const std::vector< Real > &uold, const std::vector< Real > &zold)
void linear_solve(std::vector< Real > &u, const std::vector< Real > &dl, const std::vector< Real > &d, const std::vector< Real > &du, const std::vector< Real > &r, const bool transpose=false)
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 .
Constraint_BurgersControl(int nx=128, int nt=100, Real T=1, Real nu=1.e-2, Real u0=0.0, Real u1=0.0, Real f=0.0)
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
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...
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)
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 apply_pde_jacobian_old(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
void solve(ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
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 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.
Real evaluate_target(Real x)
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 hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
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.
void hessVec_12(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.
Objective_BurgersControl(Real alpha=1.e-4, int nx=128, int nt=100, Real T=1.0)
Defines the constraint operator interface for simulation-based optimization.
Provides the interface to evaluate simulation-based objective functions.
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void zero()
Set to zero vector.