61 #include "Teuchos_oblackholestream.hpp"
62 #include "Teuchos_GlobalMPISession.hpp"
63 #include "Teuchos_XMLParameterListHelpers.hpp"
64 #include "Teuchos_LAPACK.hpp"
88 return std::sqrt(
dot(r,r));
91 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
93 Real c = ((x.size()==
nx_) ? 4.0 : 2.0);
94 for (
unsigned i = 0; i < x.size(); i++) {
96 ip +=
dx_/6.0*(c*x[i] + x[i+1])*y[i];
98 else if ( i == x.size()-1 ) {
99 ip +=
dx_/6.0*(x[i-1] + c*x[i])*y[i];
102 ip +=
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
108 void update(std::vector<Real> &u,
const std::vector<Real> &s,
const Real alpha=1.0) {
109 for (
unsigned i = 0; i < u.size(); i++) {
114 void scale(std::vector<Real> &u,
const Real alpha=0.0) {
115 for (
unsigned i = 0; i < u.size(); i++) {
120 void compute_residual(std::vector<Real> &r,
const std::vector<Real> &uold,
const std::vector<Real> &zold,
121 const std::vector<Real> &unew,
const std::vector<Real> &znew) {
124 for (
unsigned n = 0; n <
nx_; n++) {
138 r[n] -= 0.5*
dt_*unew[n-1]*(unew[n-1]+unew[n])/6.0;
139 r[n] -= 0.5*
dt_*uold[n-1]*(uold[n-1]+uold[n])/6.0;
142 r[n] += 0.5*
dt_*unew[n+1]*(unew[n]+unew[n+1])/6.0;
143 r[n] += 0.5*
dt_*uold[n+1]*(uold[n]+uold[n+1])/6.0;
146 r[n] -= 0.5*
dt_*
dx_/6.0*(znew[n]+4.0*znew[n+1]+znew[n+2]);
147 r[n] -= 0.5*
dt_*
dx_/6.0*(zold[n]+4.0*zold[n+1]+zold[n+2]);
157 const std::vector<Real> &u) {
166 for (
unsigned n = 0; n <
nx_; n++) {
168 dl[n] += 0.5*
dt_*(-2.0*u[n]-u[n+1])/6.0;
169 d[n] += 0.5*
dt_*u[n+1]/6.0;
172 d[n] -= 0.5*
dt_*u[n-1]/6.0;
173 du[n-1] += 0.5*
dt_*(u[n-1]+2.0*u[n])/6.0;
178 d[nx_-1] += 0.5*
dt_*
u1_/6.0;
182 bool adjoint =
false) {
186 for (
unsigned n = 0; n <
nx_; n++) {
191 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]);
194 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]);
200 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]);
203 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]);
207 jv[0] -= 0.5*
dt_*
u0_/6.0*v[0];
208 jv[nx_-1] += 0.5*
dt_*
u1_/6.0*v[nx_-1];
212 bool adjoint =
false) {
216 for (
unsigned n = 0; n <
nx_; n++) {
221 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]);
224 jv[n] -= 0.5*
dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]);
230 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]);
233 jv[n] += 0.5*
dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]);
237 jv[0] -= 0.5*
dt_*
u0_/6.0*v[0];
238 jv[nx_-1] += 0.5*
dt_*
u1_/6.0*v[nx_-1];
241 void apply_pde_jacobian(std::vector<Real> &jv,
const std::vector<Real> &vold,
const std::vector<Real> &uold,
242 const std::vector<Real> &vnew,
const std::vector<Real> unew,
bool adjoint =
false) {
246 for (
unsigned n = 0; n <
nx_; n++) {
253 jv[n] -= 0.5*
dt_*(unew[n-1]/6.0*vnew[n]-(unew[n-1]+2.0*unew[n])/6.0*vnew[n-1]);
254 jv[n] -= 0.5*
dt_*(uold[n-1]/6.0*vold[n]-(uold[n-1]+2.0*uold[n])/6.0*vold[n-1]);
257 jv[n] -= 0.5*
dt_*(unew[n-1]/6.0*vnew[n]+(unew[n]+2.0*unew[n-1])/6.0*vnew[n-1]);
258 jv[n] -= 0.5*
dt_*(uold[n-1]/6.0*vold[n]+(uold[n]+2.0*uold[n-1])/6.0*vold[n-1]);
265 jv[n] += 0.5*
dt_*(unew[n+1]/6.0*vnew[n]-(unew[n+1]+2.0*unew[n])/6.0*vnew[n+1]);
266 jv[n] += 0.5*
dt_*(uold[n+1]/6.0*vold[n]-(uold[n+1]+2.0*uold[n])/6.0*vold[n+1]);
269 jv[n] += 0.5*
dt_*(unew[n+1]/6.0*vnew[n]+(unew[n]+2.0*unew[n+1])/6.0*vnew[n+1]);
270 jv[n] += 0.5*
dt_*(uold[n+1]/6.0*vold[n]+(uold[n]+2.0*uold[n+1])/6.0*vold[n+1]);
274 jv[0] -= 0.5*
dt_*
u0_/6.0*vnew[0];
275 jv[0] -= 0.5*
dt_*
u0_/6.0*vold[0];
276 jv[nx_-1] += 0.5*
dt_*
u1_/6.0*vnew[nx_-1];
277 jv[nx_-1] += 0.5*
dt_*
u1_/6.0*vold[nx_-1];
280 void apply_pde_hessian(std::vector<Real> &hv,
const std::vector<Real> &wold,
const std::vector<Real> &vold,
281 const std::vector<Real> &wnew,
const std::vector<Real> &vnew) {
284 for (
unsigned n = 0; n <
nx_; n++) {
286 hv[n] += 0.5*
dt_*((wnew[n-1]*(vnew[n-1]+2.0*vnew[n]) - wnew[n]*vnew[n-1])/6.0);
287 hv[n] += 0.5*
dt_*((wold[n-1]*(vold[n-1]+2.0*vold[n]) - wold[n]*vold[n-1])/6.0);
290 hv[n] += 0.5*
dt_*((wnew[n]*vnew[n+1] - wnew[n+1]*(2.0*vnew[n]+vnew[n+1]))/6.0);
291 hv[n] += 0.5*
dt_*((wold[n]*vold[n+1] - wold[n+1]*(2.0*vold[n]+vold[n+1]))/6.0);
298 unsigned dim = ((adjoint ==
true) ?
nx_+2 :
nx_);
300 for (
unsigned n = 0; n < dim; n++) {
303 jv[n] = -0.5*
dt_*
dx_/6.0*v[n];
306 jv[n] = -0.5*
dt_*
dx_/6.0*(4.0*v[n-1]+v[n]);
308 else if ( n ==
nx_ ) {
309 jv[n] = -0.5*
dt_*
dx_/6.0*(4.0*v[n-1]+v[n-2]);
311 else if ( n ==
nx_+1 ) {
312 jv[n] = -0.5*
dt_*
dx_/6.0*v[n-2];
315 jv[n] = -0.5*
dt_*
dx_/6.0*(v[n-2]+4.0*v[n-1]+v[n]);
319 jv[n] -= 0.5*
dt_*
dx_/6.0*(v[n]+4.0*v[n+1]+v[n+2]);
324 void run_newton(std::vector<Real> &u,
const std::vector<Real> &znew,
325 const std::vector<Real> &uold,
const std::vector<Real> &zold) {
329 std::vector<Real> r(
nx_,0.0);
334 unsigned maxit = 500;
336 std::vector<Real> d(
nx_,0.0);
337 std::vector<Real> dl(
nx_-1,0.0);
338 std::vector<Real> du(
nx_-1,0.0);
340 Real alpha = 1.0, tmp = 0.0;
341 std::vector<Real> s(
nx_,0.0);
342 std::vector<Real> utmp(
nx_,0.0);
343 for (
unsigned i = 0; i < maxit; i++) {
352 utmp.assign(u.begin(),u.end());
356 while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(
ROL::ROL_EPSILON) ) {
358 utmp.assign(u.begin(),u.end());
364 u.assign(utmp.begin(),utmp.end());
365 if ( rnorm < rtol ) {
372 const std::vector<Real> &dl,
const std::vector<Real> &d,
const std::vector<Real> &du,
373 const std::vector<Real> &r,
const bool transpose =
false) {
374 bool useLAPACK =
false;
376 u.assign(r.begin(),r.end());
378 std::vector<Real> Dl(dl);
379 std::vector<Real> Du(du);
380 std::vector<Real> D(d);
382 Teuchos::LAPACK<int,Real> lp;
383 std::vector<Real> Du2(
nx_-2,0.0);
384 std::vector<int> ipiv(
nx_,0);
388 lp.GTTRF(
nx_,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&info);
389 char trans = ((transpose ==
true) ?
'T' :
'N');
390 lp.GTTRS(trans,
nx_,nhrs,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&u[0],ldb,&info);
395 unsigned maxit = 100;
396 Real rtol = std::min(1.e-12,1.e-4*std::sqrt(
dot(r,r)));
398 Real rnorm = 10.0*rtol;
399 for (
unsigned i = 0; i < maxit; i++) {
400 for (
unsigned n = 0; n <
nx_; n++) {
403 u[n] -= ((transpose ==
false) ? du[n] : dl[n])*u[n+1]/d[n];
406 u[n] -= ((transpose ==
false) ? dl[n-1] : du[n-1])*u[n-1]/d[n];
411 for (
unsigned n = 0; n <
nx_; n++) {
412 resid = r[n] - d[n]*u[n];
414 resid -= ((transpose ==
false) ? du[n] : dl[n])*u[n+1];
417 resid -= ((transpose ==
false) ? dl[n-1] : du[n-1])*u[n-1];
419 rnorm += resid*resid;
421 rnorm = std::sqrt(rnorm);
422 if ( rnorm < rtol ) {
433 Real nu = 1.e-2, Real u0 = 0.0, Real u1 = 0.0, Real f = 0.0)
435 dx_ = 1.0/((Real)nx+1.0);
440 for (
unsigned n = 0; n <
nx_; n++) {
442 u_init_[n] = ((x <= 0.5) ? 1.0 : 0.0);
447 Teuchos::RCP<std::vector<Real> > cp =
448 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(c)).getVector());
449 Teuchos::RCP<const std::vector<Real> > up =
451 Teuchos::RCP<const std::vector<Real> > zp =
454 std::vector<Real> C(
nx_,0.0);
455 std::vector<Real> uold(
u_init_);
456 std::vector<Real> unew(
u_init_);
457 std::vector<Real> zold(
nx_+2,0.0);
458 std::vector<Real> znew(
nx_+2,0.0);
460 for (
unsigned n = 0; n <
nx_+2; n++) {
464 for (
unsigned t = 0; t <
nt_; t++) {
466 for (
unsigned n = 0; n <
nx_; n++) {
467 unew[n] = (*up)[t*nx_+n];
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 (*cp)[t*nx_+n] = C[n];
480 uold.assign(unew.begin(),unew.end());
481 zold.assign(znew.begin(),znew.end());
486 Teuchos::RCP<std::vector<Real> > up =
487 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(u)).getVector());
488 up->assign(up->size(),z.
norm()/up->size());
489 Teuchos::RCP<const std::vector<Real> > zp =
492 std::vector<Real> uold(
u_init_);
493 std::vector<Real> unew(
u_init_);
494 std::vector<Real> zold(
nx_+2,0.0);
495 std::vector<Real> znew(
nx_+2,0.0);
497 for (
unsigned n = 0; n <
nx_+2; n++) {
501 for (
unsigned t = 0; t <
nt_; t++) {
503 for (
unsigned n = 0; n < nx_+2; n++) {
504 znew[n] = (*zp)[(t+1)*(nx_+2)+n];
509 for (
unsigned n = 0; n <
nx_; n++) {
510 (*up)[t*nx_+n] = unew[n];
513 uold.assign(unew.begin(),unew.end());
514 zold.assign(znew.begin(),znew.end());
520 Teuchos::RCP<std::vector<Real> > jvp =
521 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(jv)).getVector());
522 Teuchos::RCP<const std::vector<Real> > vp =
524 Teuchos::RCP<const std::vector<Real> > up =
526 std::vector<Real> J(
nx_,0.0);
527 std::vector<Real> vold(
nx_,0.0);
528 std::vector<Real> vnew(
nx_,0.0);
529 std::vector<Real> uold(
nx_,0.0);
530 std::vector<Real> unew(
nx_,0.0);
531 for (
unsigned t = 0; t <
nt_; t++) {
532 for (
unsigned n = 0; n <
nx_; n++) {
533 unew[n] = (*up)[t*nx_+n];
534 vnew[n] = (*vp)[t*nx_+n];
537 for (
unsigned n = 0; n <
nx_; n++) {
538 (*jvp)[t*nx_+n] = J[n];
540 vold.assign(vnew.begin(),vnew.end());
541 uold.assign(unew.begin(),unew.end());
547 Teuchos::RCP<std::vector<Real> > jvp =
548 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(jv)).getVector());
549 Teuchos::RCP<const std::vector<Real> > vp =
551 Teuchos::RCP<const std::vector<Real> > zp =
553 std::vector<Real> vnew(
nx_+2,0.0);
554 std::vector<Real> vold(
nx_+2,0.0);
555 std::vector<Real> jnew(
nx_,0.0);
556 std::vector<Real> jold(
nx_,0.0);
557 for (
unsigned n = 0; n <
nx_+2; n++) {
561 for (
unsigned t = 0; t <
nt_; t++) {
562 for (
unsigned n = 0; n < nx_+2; n++) {
563 vnew[n] = (*vp)[(t+1)*(nx_+2)+n];
566 for (
unsigned n = 0; n <
nx_; n++) {
568 (*jvp)[t*nx_+n] = jnew[n] + jold[n];
570 jold.assign(jnew.begin(),jnew.end());
576 Teuchos::RCP<std::vector<Real> > ijvp =
577 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(ijv)).getVector());
578 Teuchos::RCP<const std::vector<Real> > vp =
580 Teuchos::RCP<const std::vector<Real> > up =
582 std::vector<Real> J(
nx_,0.0);
583 std::vector<Real> r(
nx_,0.0);
584 std::vector<Real> s(
nx_,0.0);
585 std::vector<Real> uold(
nx_,0.0);
586 std::vector<Real> unew(
nx_,0.0);
587 std::vector<Real> d(
nx_,0.0);
588 std::vector<Real> dl(
nx_-1,0.0);
589 std::vector<Real> du(
nx_-1,0.0);
590 for (
unsigned t = 0; t <
nt_; t++) {
592 for (
unsigned n = 0; n <
nx_; n++) {
593 r[n] = (*vp)[t*nx_+n];
594 unew[n] = (*up)[t*nx_+n];
603 for (
unsigned n = 0; n <
nx_; n++) {
604 (*ijvp)[t*nx_+n] = s[n];
606 uold.assign(unew.begin(),unew.end());
612 Teuchos::RCP<std::vector<Real> > jvp =
613 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(ajv)).getVector());
614 Teuchos::RCP<const std::vector<Real> > vp =
616 Teuchos::RCP<const std::vector<Real> > up =
618 std::vector<Real> J(
nx_,0.0);
619 std::vector<Real> vold(
nx_,0.0);
620 std::vector<Real> vnew(
nx_,0.0);
621 std::vector<Real> unew(
nx_,0.0);
622 for (
unsigned t =
nt_; t > 0; t--) {
623 for (
unsigned n = 0; n <
nx_; n++) {
624 unew[n] = (*up)[(t-1)*nx_+n];
625 vnew[n] = (*vp)[(t-1)*nx_+n];
628 for (
unsigned n = 0; n <
nx_; n++) {
629 (*jvp)[(t-1)*nx_+n] = J[n];
631 vold.assign(vnew.begin(),vnew.end());
637 Teuchos::RCP<std::vector<Real> > jvp =
638 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(jv)).getVector());
639 Teuchos::RCP<const std::vector<Real> > vp =
641 Teuchos::RCP<const std::vector<Real> > zp =
643 std::vector<Real> vnew(
nx_,0.0);
644 std::vector<Real> vold(
nx_,0.0);
645 std::vector<Real> jnew(
nx_+2,0.0);
646 std::vector<Real> jold(
nx_+2,0.0);
647 for (
unsigned t =
nt_+1; t > 0; t--) {
648 for (
unsigned n = 0; n <
nx_; n++) {
650 vnew[n] = (*vp)[(t-2)*nx_+n];
657 for (
unsigned n = 0; n < nx_+2; n++) {
659 (*jvp)[(t-1)*(nx_+2)+n] = jnew[n] + jold[n];
661 jold.assign(jnew.begin(),jnew.end());
667 Teuchos::RCP<std::vector<Real> > ijvp =
668 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(ijv)).getVector());
669 Teuchos::RCP<const std::vector<Real> > vp =
671 Teuchos::RCP<const std::vector<Real> > up =
673 std::vector<Real> J(
nx_,0.0);
674 std::vector<Real> r(
nx_,0.0);
675 std::vector<Real> s(
nx_,0.0);
676 std::vector<Real> unew(
nx_,0.0);
677 std::vector<Real> d(
nx_,0.0);
678 std::vector<Real> dl(
nx_-1,0.0);
679 std::vector<Real> du(
nx_-1,0.0);
680 for (
unsigned t =
nt_; t > 0; t--) {
682 for (
unsigned n = 0; n <
nx_; n++) {
683 r[n] = (*vp)[(t-1)*nx_+n];
684 unew[n] = (*up)[(t-1)*nx_+n];
693 for (
unsigned n = 0; n <
nx_; n++) {
694 (*ijvp)[(t-1)*nx_+n] = s[n];
701 Teuchos::RCP<std::vector<Real> > hwvp =
702 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<
ROL::StdVector<Real> >(hwv)).getVector());
703 Teuchos::RCP<const std::vector<Real> > wp =
705 Teuchos::RCP<const std::vector<Real> > vp =
707 std::vector<Real> snew(
nx_,0.0);
708 std::vector<Real> wnew(
nx_,0.0);
709 std::vector<Real> wold(
nx_,0.0);
710 std::vector<Real> vnew(
nx_,0.0);
711 for (
unsigned t =
nt_; t > 0; t--) {
712 for (
unsigned n = 0; n <
nx_; n++) {
713 vnew[n] = (*vp)[(t-1)*nx_+n];
714 wnew[n] = (*wp)[(t-1)*nx_+n];
717 for (
unsigned n = 0; n <
nx_; n++) {
718 (*hwvp)[(t-1)*nx_+n] = snew[n];
720 wold.assign(wnew.begin(),wnew.end());
756 case 1: val = ((x<=0.5) ? 1.0 : 0.0);
break;
757 case 2: val = 1.0;
break;
758 case 3: val = std::abs(std::sin(8.0*M_PI*x));
break;
759 case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5));
break;
764 Real
dot(
const std::vector<Real> &x,
const std::vector<Real> &y) {
766 Real c = ((x.size()==
nx_) ? 4.0 : 2.0);
767 for (
unsigned i=0; i<x.size(); i++) {
769 ip +=
dx_/6.0*(c*x[i] + x[i+1])*y[i];
771 else if ( i == x.size()-1 ) {
772 ip +=
dx_/6.0*(x[i-1] + c*x[i])*y[i];
775 ip +=
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
781 void apply_mass(std::vector<Real> &Mu,
const std::vector<Real> &u ) {
782 Mu.resize(u.size(),0.0);
783 Real c = ((u.size()==
nx_) ? 4.0 : 2.0);
784 for (
unsigned i=0; i<u.size(); i++) {
786 Mu[i] =
dx_/6.0*(c*u[i] + u[i+1]);
788 else if ( i == u.size()-1 ) {
789 Mu[i] =
dx_/6.0*(u[i-1] + c*u[i]);
792 Mu[i] =
dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
804 dx_ = 1.0/((Real)nx+1.0);
809 Teuchos::RCP<const std::vector<Real> > up =
811 Teuchos::RCP<const std::vector<Real> > zp =
814 std::vector<Real> U(
nx_,0.0);
815 std::vector<Real> Z(
nx_+2,0.0);
816 for (
unsigned n = 0; n <
nx_+2; n++) {
821 for (
unsigned t = 0; t <
nt_; t++) {
822 ss = ((t < nt_-1) ?
dt_ : 0.5*
dt_);
823 for (
unsigned n = 0; n <
nx_; n++) {
826 val += 0.5*ss*
dot(U,U);
827 for (
unsigned n = 0; n < nx_+2; n++) {
828 Z[n] = (*zp)[(t+1)*(nx_+2)+n];
836 Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >(
838 Teuchos::RCP<const std::vector<Real> > up =
841 std::vector<Real> U(
nx_,0.0);
842 std::vector<Real> M(
nx_,0.0);
844 for (
unsigned t = 0; t <
nt_; t++) {
845 ss = ((t < nt_-1) ?
dt_ : 0.5*
dt_);
846 for (
unsigned n = 0; n <
nx_; n++) {
850 for (
unsigned n = 0; n <
nx_; n++) {
851 (*gp)[t*nx_+n] = ss*M[n];
857 Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >(
859 Teuchos::RCP<const std::vector<Real> > zp =
862 std::vector<Real> Z(
nx_+2,0.0);
863 std::vector<Real> M(
nx_+2,0.0);
865 for (
unsigned t = 0; t <
nt_+1; t++) {
866 ss = ((t < nt_ && t > 0) ?
dt_ : 0.5*
dt_);
867 for (
unsigned n = 0; n <
nx_+2; n++) {
868 Z[n] = (*zp)[t*(nx_+2)+n];
871 for (
unsigned n = 0; n < nx_+2; n++) {
872 (*gp)[t*(nx_+2)+n] = ss*
alpha_*M[n];
879 Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >(
881 Teuchos::RCP<const std::vector<Real> > vp =
884 std::vector<Real> V(
nx_,0.0);
885 std::vector<Real> M(
nx_,0.0);
887 for (
unsigned t = 0; t <
nt_; t++) {
888 ss = ((t < nt_-1) ?
dt_ : 0.5*
dt_);
889 for (
unsigned n = 0; n <
nx_; n++) {
890 V[n] = (*vp)[t*nx_+n];
893 for (
unsigned n = 0; n <
nx_; n++) {
894 (*hvp)[t*nx_+n] = ss*M[n];
911 Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >(
913 Teuchos::RCP<const std::vector<Real> > vp =
916 std::vector<Real> V(
nx_+2,0.0);
917 std::vector<Real> M(
nx_+2,0.0);
919 for (
unsigned t = 0; t <
nt_+1; t++) {
920 ss = ((t < nt_ && t > 0) ?
dt_ : 0.5*
dt_);
921 for (
unsigned n = 0; n <
nx_+2; n++) {
922 V[n] = (*vp)[t*(nx_+2)+n];
925 for (
unsigned n = 0; n < nx_+2; n++) {
926 (*hvp)[t*(nx_+2)+n] = ss*
alpha_*M[n];
934 int main(
int argc,
char *argv[]) {
936 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
939 int iprint = argc - 1;
940 Teuchos::RCP<std::ostream> outStream;
941 Teuchos::oblackholestream bhs;
943 outStream = Teuchos::rcp(&std::cout,
false);
945 outStream = Teuchos::rcp(&bhs,
false);
962 Teuchos::RCP<std::vector<RealT> > z_rcp = Teuchos::rcp(
new std::vector<RealT> ((nx+2)*(nt+1), 1.0) );
963 Teuchos::RCP<std::vector<RealT> > gz_rcp = Teuchos::rcp(
new std::vector<RealT> ((nx+2)*(nt+1), 1.0) );
964 Teuchos::RCP<std::vector<RealT> > yz_rcp = Teuchos::rcp(
new std::vector<RealT> ((nx+2)*(nt+1), 1.0) );
965 for (
int i=0; i<(nx+2)*(nt+1); i++) {
972 Teuchos::RCP<ROL::Vector<RealT> > zp = Teuchos::rcp(&z,
false);
973 Teuchos::RCP<ROL::Vector<RealT> > gzp = Teuchos::rcp(&gz,
false);
974 Teuchos::RCP<ROL::Vector<RealT> > yzp = Teuchos::rcp(&yz,
false);
976 Teuchos::RCP<std::vector<RealT> > u_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
977 Teuchos::RCP<std::vector<RealT> > gu_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
978 Teuchos::RCP<std::vector<RealT> > yu_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
979 for (
int i=0; i<nx*nt; i++) {
986 Teuchos::RCP<ROL::Vector<RealT> > up = Teuchos::rcp(&u,
false);
987 Teuchos::RCP<ROL::Vector<RealT> > gup = Teuchos::rcp(&gu,
false);
988 Teuchos::RCP<ROL::Vector<RealT> > yup = Teuchos::rcp(&yu,
false);
990 Teuchos::RCP<std::vector<RealT> > c_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
991 Teuchos::RCP<std::vector<RealT> > l_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
1013 Teuchos::RCP<std::vector<RealT> > p_rcp = Teuchos::rcp(
new std::vector<RealT> (nx*nt, 1.0) );
1015 Teuchos::RCP<ROL::Vector<RealT> > pp = Teuchos::rcp(&p,
false);
1016 Teuchos::RCP<ROL::Objective_SimOpt<RealT> > pobj = Teuchos::rcp(&obj,
false);
1017 Teuchos::RCP<ROL::EqualityConstraint_SimOpt<RealT> > pcon = Teuchos::rcp(&con,
false);
1023 std::string filename =
"input.xml";
1024 Teuchos::RCP<Teuchos::ParameterList> parlist_tr = Teuchos::rcp(
new Teuchos::ParameterList() );
1025 Teuchos::updateParametersFromXmlFile( filename, Teuchos::Ptr<Teuchos::ParameterList>(&*parlist_tr) );
1035 std::clock_t timer_tr = std::clock();
1036 algo_tr.
run(z,robj,
true,*outStream);
1037 *outStream <<
"Trust-Region Newton required " << (std::clock()-timer_tr)/(
RealT)CLOCKS_PER_SEC
1041 RealT ctol = 1.e-14;
1046 std::clock_t timer_sqp = std::clock();
1047 algo_sqp.
run(x,g,l,c,obj,con,
true,*outStream);
1048 *outStream <<
"Composite-Step SQP required " << (std::clock()-timer_sqp)/(
RealT)CLOCKS_PER_SEC
1051 std::ofstream control;
1052 control.open(
"control.txt");
1053 for (
int t = 0; t < nt+1; t++) {
1054 for (
int n = 0; n < nx+2; n++) {
1057 << (*z_rcp)[t*(nx+2)+n] <<
"\n";
1062 std::ofstream state;
1063 state.open(
"state.txt");
1064 for (
int t = 0; t < nt; t++) {
1065 for (
int n = 0; n < nx; n++) {
1068 << (*u_rcp)[t*nx+n] <<
"\n";
1073 catch (std::logic_error err) {
1074 *outStream << err.what() <<
"\n";
1079 std::cout <<
"End Result: TEST FAILED\n";
1081 std::cout <<
"End Result: TEST PASSED\n";
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)
Provides the interface to evaluate simulation-based objective functions.
void run_newton(std::vector< Real > &u, const std::vector< Real > &znew, const std::vector< Real > &uold, const std::vector< Real > &zold)
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
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)
Real evaluate_target(Real x)
Defines the linear algebra or vector space interface for simulation-based optimization.
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
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 .
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Contains definitions of custom data types in ROL.
virtual std::vector< std::vector< Real > > checkApplyAdjointJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS)
Finite-difference check for the application of the adjoint of constraint 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)
Real compute_norm(const std::vector< Real > &r)
EqualityConstraint_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)
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
void apply_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, bool adjoint=false)
int main(int argc, char *argv[])
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
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)
Defines the equality constraint operator interface for simulation-based optimization.
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 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.
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.
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.
Provides the interface to evaluate simulation-based reduced 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...
void apply_pde_jacobian_new(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
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...
Provides the std::vector implementation of the ROL::Vector interface.
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 .
virtual Real checkSolve(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, const ROL::Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Implements the computation of optimization steps with composite-step trust-region SQP methods...
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
std::vector< Real > u_init_
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Provides an interface to check status of optimization algorithms.
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 solve(ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
void scale(std::vector< Real > &u, const Real alpha=0.0)
Objective_BurgersControl(Real alpha=1.e-4, int nx=128, int nt=100, Real T=1.0)
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 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 value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
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 apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
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)
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_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)
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.
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)
virtual Real norm() const =0
Returns where .
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. ...
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
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.
void apply_pde_jacobian_old(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
Provides the interface to compute optimization steps with trust regions.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
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 adjoint of the partial constraint Hessian at , , to vector in direction ...