54 #include "ROL_TeuchosBatchManager.hpp"
56 #include "Teuchos_LAPACK.hpp"
84 void update(std::vector<Real> &u,
const std::vector<Real> &s,
const Real alpha=1.0)
const {
85 for (
unsigned i=0; i<u.size(); i++) {
90 void axpy(std::vector<Real> &out,
const Real a,
const std::vector<Real> &x,
const std::vector<Real> &y)
const {
91 for (
unsigned i=0; i < x.size(); i++) {
92 out[i] = a*x[i] + y[i];
96 void scale(std::vector<Real> &u,
const Real alpha=0.0)
const {
97 for (
unsigned i=0; i<u.size(); i++) {
102 void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
103 const std::vector<Real> &r,
const bool transpose =
false)
const {
104 if ( r.size() == 1 ) {
105 u.resize(1,r[0]/d[0]);
107 else if ( r.size() == 2 ) {
109 Real det = d[0]*d[1] - dl[0]*du[0];
110 u[0] = (d[1]*r[0] - du[0]*r[1])/det;
111 u[1] = (d[0]*r[1] - dl[0]*r[0])/det;
114 u.assign(r.begin(),r.end());
116 Teuchos::LAPACK<int,Real> lp;
117 std::vector<Real> du2(r.size()-2,0.0);
118 std::vector<int> ipiv(r.size(),0);
123 lp.GTTRF(dim,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
128 lp.GTTRS(trans,dim,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
133 BurgersFEM(
int nx = 128, Real nl = 1.0, Real cH1 = 1.0, Real cL2 = 1.0)
137 nu_ = std::pow(10.0,nu-2.0);
155 Real
compute_L2_dot(
const std::vector<Real> &x,
const std::vector<Real> &y)
const {
157 Real c = (((int)x.size()==
nx_) ? 4.0 : 2.0);
158 for (
unsigned i=0; i<x.size(); i++) {
160 ip +=
dx_/6.0*(c*x[i] + x[i+1])*y[i];
162 else if ( i == x.size()-1 ) {
163 ip +=
dx_/6.0*(x[i-1] + c*x[i])*y[i];
166 ip +=
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
178 void apply_mass(std::vector<Real> &Mu,
const std::vector<Real> &u )
const {
179 Mu.resize(u.size(),0.0);
180 Real c = (((int)u.size()==
nx_) ? 4.0 : 2.0);
181 for (
unsigned i=0; i<u.size(); i++) {
183 Mu[i] =
dx_/6.0*(c*u[i] + u[i+1]);
185 else if ( i == u.size()-1 ) {
186 Mu[i] =
dx_/6.0*(u[i-1] + c*u[i]);
189 Mu[i] =
dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
196 unsigned nx = u.size();
198 std::vector<Real> dl(nx-1,
dx_/6.0);
199 std::vector<Real> d(nx,2.0*
dx_/3.0);
200 std::vector<Real> du(nx-1,
dx_/6.0);
201 if ( (
int)nx !=
nx_ ) {
210 std::vector<Real> u(
nx_,0.0), Mu(
nx_,0.0), iMMu(
nx_,0.0), diff(
nx_,0.0);
211 for (
int i = 0; i <
nx_; i++) {
212 u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
216 axpy(diff,-1.0,iMMu,u);
219 outStream <<
"Test Inverse State Mass Matrix\n";
220 outStream <<
" ||u - inv(M)Mu|| = " << error <<
"\n";
221 outStream <<
" ||u|| = " << normu <<
"\n";
222 outStream <<
" Relative Error = " << error/normu <<
"\n";
225 u.resize(nx_+2,0.0); Mu.resize(nx_+2,0.0); iMMu.resize(nx_+2,0.0); diff.resize(nx_+2,0.0);
226 for (
int i = 0; i < nx_+2; i++) {
227 u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
231 axpy(diff,-1.0,iMMu,u);
234 outStream <<
"Test Inverse Control Mass Matrix\n";
235 outStream <<
" ||z - inv(M)Mz|| = " << error <<
"\n";
236 outStream <<
" ||z|| = " << normu <<
"\n";
237 outStream <<
" Relative Error = " << error/normu <<
"\n";
245 Real
compute_H1_dot(
const std::vector<Real> &x,
const std::vector<Real> &y)
const {
247 for (
int i=0; i<
nx_; i++) {
249 ip +=
cL2_*
dx_/6.0*(4.0*x[i] + x[i+1])*y[i];
250 ip +=
cH1_*(2.0*x[i] - x[i+1])/
dx_*y[i];
252 else if ( i == nx_-1 ) {
253 ip +=
cL2_*
dx_/6.0*(x[i-1] + 4.0*x[i])*y[i];
254 ip +=
cH1_*(2.0*x[i] - x[i-1])/
dx_*y[i];
257 ip +=
cL2_*
dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
258 ip +=
cH1_*(2.0*x[i] - x[i-1] - x[i+1])/
dx_*y[i];
270 void apply_H1(std::vector<Real> &Mu,
const std::vector<Real> &u )
const {
272 for (
int i=0; i<
nx_; i++) {
274 Mu[i] =
cL2_*
dx_/6.0*(4.0*u[i] + u[i+1])
275 +
cH1_*(2.0*u[i] - u[i+1])/
dx_;
277 else if ( i == nx_-1 ) {
278 Mu[i] =
cL2_*
dx_/6.0*(u[i-1] + 4.0*u[i])
279 +
cH1_*(2.0*u[i] - u[i-1])/
dx_;
282 Mu[i] =
cL2_*
dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1])
283 +
cH1_*(2.0*u[i] - u[i-1] - u[i+1])/
dx_;
298 std::vector<Real> u(
nx_,0.0), Mu(
nx_,0.0), iMMu(
nx_,0.0), diff(
nx_,0.0);
299 for (
int i = 0; i <
nx_; i++) {
300 u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
304 axpy(diff,-1.0,iMMu,u);
307 outStream <<
"Test Inverse State H1 Matrix\n";
308 outStream <<
" ||u - inv(M)Mu|| = " << error <<
"\n";
309 outStream <<
" ||u|| = " << normu <<
"\n";
310 outStream <<
" Relative Error = " << error/normu <<
"\n";
319 const std::vector<Real> &z)
const {
322 for (
int i=0; i<
nx_; i++) {
325 r[i] =
nu_/
dx_*(2.0*u[i]-u[i+1]);
328 r[i] =
nu_/
dx_*(2.0*u[i]-u[i-1]);
331 r[i] =
nu_/
dx_*(2.0*u[i]-u[i-1]-u[i+1]);
335 r[i] +=
nl_*u[i+1]*(u[i]+u[i+1])/6.0;
338 r[i] -=
nl_*u[i-1]*(u[i-1]+u[i])/6.0;
341 r[i] -=
dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
355 std::vector<Real> &d,
356 std::vector<Real> &du,
357 const std::vector<Real> &u)
const {
366 for (
int i=0; i<
nx_; i++) {
368 dl[i] +=
nl_*(-2.0*u[i]-u[i+1])/6.0;
369 d[i] +=
nl_*u[i+1]/6.0;
372 d[i] -=
nl_*u[i-1]/6.0;
373 du[i-1] +=
nl_*(u[i-1]+2.0*u[i])/6.0;
383 const std::vector<Real> &v,
384 const std::vector<Real> &u,
385 const std::vector<Real> &z)
const {
387 for (
int i = 0; i <
nx_; i++) {
390 jv[i] += -
nu_/
dx_*v[i-1]-
nl_*(u[i-1]/6.0*v[i]+(u[i]+2.0*u[i-1])/6.0*v[i-1]);
393 jv[i] += -
nu_/
dx_*v[i+1]+
nl_*(u[i+1]/6.0*v[i]+(u[i]+2.0*u[i+1])/6.0*v[i+1]);
396 jv[ 0] -=
nl_*
u0_/6.0*v[0];
397 jv[nx_-1] +=
nl_*
u1_/6.0*v[nx_-1];
402 const std::vector<Real> &v,
403 const std::vector<Real> &u,
404 const std::vector<Real> &z)
const {
406 std::vector<Real> d(
nx_,0.0);
407 std::vector<Real> dl(
nx_-1,0.0);
408 std::vector<Real> du(
nx_-1,0.0);
416 const std::vector<Real> &v,
417 const std::vector<Real> &u,
418 const std::vector<Real> &z)
const {
420 for (
int i = 0; i <
nx_; i++) {
421 ajv[i] =
nu_/
dx_*2.0*v[i];
423 ajv[i] += -
nu_/
dx_*v[i-1]-
nl_*(u[i-1]/6.0*v[i]
424 -(u[i-1]+2.0*u[i])/6.0*v[i-1]);
427 ajv[i] += -
nu_/
dx_*v[i+1]+
nl_*(u[i+1]/6.0*v[i]
428 -(u[i+1]+2.0*u[i])/6.0*v[i+1]);
431 ajv[ 0] -=
nl_*
u0_/6.0*v[0];
432 ajv[nx_-1] +=
nl_*
u1_/6.0*v[nx_-1];
437 const std::vector<Real> &v,
438 const std::vector<Real> &u,
439 const std::vector<Real> &z)
const {
441 std::vector<Real> d(
nx_,0.0);
442 std::vector<Real> du(
nx_-1,0.0);
443 std::vector<Real> dl(
nx_-1,0.0);
454 const std::vector<Real> &v,
455 const std::vector<Real> &u,
456 const std::vector<Real> &z)
const {
457 for (
int i=0; i<
nx_; i++) {
459 jv[i] = -
dx_/6.0*(v[i]+4.0*v[i+1]+v[i+2]);
465 const std::vector<Real> &v,
466 const std::vector<Real> &u,
467 const std::vector<Real> &z)
const {
468 for (
int i=0; i<
nx_+2; i++) {
470 jv[i] = -
dx_/6.0*v[i];
473 jv[i] = -
dx_/6.0*(4.0*v[i-1]+v[i]);
475 else if ( i == nx_ ) {
476 jv[i] = -
dx_/6.0*(4.0*v[i-1]+v[i-2]);
478 else if ( i == nx_+1 ) {
479 jv[i] = -
dx_/6.0*v[i-2];
482 jv[i] = -
dx_/6.0*(v[i-2]+4.0*v[i-1]+v[i]);
491 const std::vector<Real> &w,
492 const std::vector<Real> &v,
493 const std::vector<Real> &u,
494 const std::vector<Real> &z)
const {
495 for (
int i=0; i<
nx_; i++) {
499 ahwv[i] += (w[i]*v[i+1] - w[i+1]*(2.0*v[i]+v[i+1]))/6.0;
502 ahwv[i] += (w[i-1]*(v[i-1]+2.0*v[i]) - w[i]*v[i-1])/6.0;
508 const std::vector<Real> &w,
509 const std::vector<Real> &v,
510 const std::vector<Real> &u,
511 const std::vector<Real> &z) {
512 ahwv.assign(u.size(),0.0);
515 const std::vector<Real> &w,
516 const std::vector<Real> &v,
517 const std::vector<Real> &u,
518 const std::vector<Real> &z) {
519 ahwv.assign(z.size(),0.0);
522 const std::vector<Real> &w,
523 const std::vector<Real> &v,
524 const std::vector<Real> &u,
525 const std::vector<Real> &z) {
526 ahwv.assign(z.size(),0.0);
533 ROL::Ptr<std::vector<Real> >
vec_;
534 ROL::Ptr<BurgersFEM<Real> >
fem_;
536 mutable ROL::Ptr<L2VectorDual<Real> >
dual_vec_;
545 const std::vector<Real>& xval = *ex.
getVector();
546 std::copy(xval.begin(),xval.end(),
vec_->begin());
551 const std::vector<Real>& xval = *ex.
getVector();
554 (*vec_)[i] += xval[i];
567 const std::vector<Real>& xval = *ex.
getVector();
568 return fem_->compute_L2_dot(xval,*
vec_);
573 val = std::sqrt(
dot(*
this) );
577 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
578 return ROL::makePtr<L2VectorPrimal>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
589 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
590 ROL::Ptr<L2VectorPrimal> e
591 = ROL::makePtr<L2VectorPrimal>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
592 (*e->getVector())[i] = 1.0;
601 dual_vec_ = ROL::makePtr<L2VectorDual<Real>>(
602 ROL::makePtr<std::vector<Real>>(*vec_),
fem_);
604 fem_->apply_mass(*(ROL::constPtrCast<std::vector<Real> >(
dual_vec_->getVector())),*
vec_);
613 ROL::Ptr<std::vector<Real> >
vec_;
614 ROL::Ptr<BurgersFEM<Real> >
fem_;
616 mutable ROL::Ptr<L2VectorPrimal<Real> >
dual_vec_;
625 const std::vector<Real>& xval = *ex.
getVector();
626 std::copy(xval.begin(),xval.end(),
vec_->begin());
631 const std::vector<Real>& xval = *ex.
getVector();
634 (*vec_)[i] += xval[i];
647 const std::vector<Real>& xval = *ex.
getVector();
649 std::vector<Real> Mx(dimension,0.0);
650 fem_->apply_inverse_mass(Mx,xval);
652 for (
unsigned i = 0; i <
dimension; i++) {
653 val += Mx[i]*(*vec_)[i];
660 val = std::sqrt(
dot(*
this) );
664 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
665 return ROL::makePtr<L2VectorDual>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
676 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
677 ROL::Ptr<L2VectorDual> e
678 = ROL::makePtr<L2VectorDual>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
679 (*e->getVector())[i] = 1.0;
688 dual_vec_ = ROL::makePtr<L2VectorPrimal<Real>>(
689 ROL::makePtr<std::vector<Real>>(*vec_),
fem_);
691 fem_->apply_inverse_mass(*(ROL::constPtrCast<std::vector<Real> >(
dual_vec_->getVector())),*
vec_);
700 ROL::Ptr<std::vector<Real> >
vec_;
701 ROL::Ptr<BurgersFEM<Real> >
fem_;
703 mutable ROL::Ptr<H1VectorDual<Real> >
dual_vec_;
712 const std::vector<Real>& xval = *ex.
getVector();
713 std::copy(xval.begin(),xval.end(),
vec_->begin());
718 const std::vector<Real>& xval = *ex.
getVector();
721 (*vec_)[i] += xval[i];
734 const std::vector<Real>& xval = *ex.
getVector();
735 return fem_->compute_H1_dot(xval,*
vec_);
740 val = std::sqrt(
dot(*
this) );
744 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
745 return ROL::makePtr<H1VectorPrimal>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
756 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
757 ROL::Ptr<H1VectorPrimal> e
758 = ROL::makePtr<H1VectorPrimal>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
759 (*e->getVector())[i] = 1.0;
768 dual_vec_ = ROL::makePtr<H1VectorDual<Real>>(
769 ROL::makePtr<std::vector<Real>>(*vec_),
fem_);
771 fem_->apply_H1(*(ROL::constPtrCast<std::vector<Real> >(
dual_vec_->getVector())),*
vec_);
780 ROL::Ptr<std::vector<Real> >
vec_;
781 ROL::Ptr<BurgersFEM<Real> >
fem_;
783 mutable ROL::Ptr<H1VectorPrimal<Real> >
dual_vec_;
792 const std::vector<Real>& xval = *ex.
getVector();
793 std::copy(xval.begin(),xval.end(),
vec_->begin());
798 const std::vector<Real>& xval = *ex.
getVector();
801 (*vec_)[i] += xval[i];
814 const std::vector<Real>& xval = *ex.
getVector();
816 std::vector<Real> Mx(dimension,0.0);
817 fem_->apply_inverse_H1(Mx,xval);
819 for (
unsigned i = 0; i <
dimension; i++) {
820 val += Mx[i]*(*vec_)[i];
827 val = std::sqrt(
dot(*
this) );
831 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
832 return ROL::makePtr<H1VectorDual>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
843 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
844 ROL::Ptr<H1VectorDual> e
845 = ROL::makePtr<H1VectorDual>( ROL::makePtr<std::vector<Real>>(
vec_->size(),0.0),
fem_);
846 (*e->getVector())[i] = 1.0;
855 dual_vec_ = ROL::makePtr<H1VectorPrimal<Real>>(
856 ROL::makePtr<std::vector<Real>>(*vec_),
fem_);
858 fem_->apply_inverse_H1(*(ROL::constPtrCast<std::vector<Real> >(
dual_vec_->getVector())),*
vec_);
877 ROL::Ptr<BurgersFEM<Real> >
fem_;
886 ROL::Ptr<std::vector<Real> > cp =
888 ROL::Ptr<const std::vector<Real> > up =
890 ROL::Ptr<const std::vector<Real> > zp =
893 const std::vector<Real> param
895 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
897 fem_->compute_residual(*cp,*up,*zp);
902 ROL::Ptr<std::vector<Real> > jvp =
904 ROL::Ptr<const std::vector<Real> > vp =
906 ROL::Ptr<const std::vector<Real> > up =
908 ROL::Ptr<const std::vector<Real> > zp =
911 const std::vector<Real> param
913 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
915 fem_->apply_pde_jacobian(*jvp,*vp,*up,*zp);
920 ROL::Ptr<std::vector<Real> > jvp =
922 ROL::Ptr<const std::vector<Real> > vp =
924 ROL::Ptr<const std::vector<Real> > up =
926 ROL::Ptr<const std::vector<Real> > zp =
929 const std::vector<Real> param
931 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
933 fem_->apply_control_jacobian(*jvp,*vp,*up,*zp);
938 ROL::Ptr<std::vector<Real> > ijvp =
940 ROL::Ptr<const std::vector<Real> > vp =
942 ROL::Ptr<const std::vector<Real> > up =
944 ROL::Ptr<const std::vector<Real> > zp =
947 const std::vector<Real> param
949 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
951 fem_->apply_inverse_pde_jacobian(*ijvp,*vp,*up,*zp);
956 ROL::Ptr<std::vector<Real> > jvp =
958 ROL::Ptr<const std::vector<Real> > vp =
960 ROL::Ptr<const std::vector<Real> > up =
962 ROL::Ptr<const std::vector<Real> > zp =
965 const std::vector<Real> param
967 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
969 fem_->apply_adjoint_pde_jacobian(*jvp,*vp,*up,*zp);
974 ROL::Ptr<std::vector<Real> > jvp =
976 ROL::Ptr<const std::vector<Real> > vp =
978 ROL::Ptr<const std::vector<Real> > up =
980 ROL::Ptr<const std::vector<Real> > zp =
983 const std::vector<Real> param
985 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
987 fem_->apply_adjoint_control_jacobian(*jvp,*vp,*up,*zp);
992 ROL::Ptr<std::vector<Real> > iajvp =
994 ROL::Ptr<const std::vector<Real> > vp =
996 ROL::Ptr<const std::vector<Real> > up =
998 ROL::Ptr<const std::vector<Real> > zp =
1001 const std::vector<Real> param
1003 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1005 fem_->apply_inverse_adjoint_pde_jacobian(*iajvp,*vp,*up,*zp);
1011 ROL::Ptr<std::vector<Real> > ahwvp =
1013 ROL::Ptr<const std::vector<Real> > wp =
1015 ROL::Ptr<const std::vector<Real> > vp =
1017 ROL::Ptr<const std::vector<Real> > up =
1019 ROL::Ptr<const std::vector<Real> > zp =
1022 const std::vector<Real> param
1024 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1026 fem_->apply_adjoint_pde_hessian(*ahwvp,*wp,*vp,*up,*zp);
1036 ROL::Ptr<std::vector<Real> > ahwvp =
1038 ROL::Ptr<const std::vector<Real> > wp =
1040 ROL::Ptr<const std::vector<Real> > vp =
1042 ROL::Ptr<const std::vector<Real> > up =
1044 ROL::Ptr<const std::vector<Real> > zp =
1047 const std::vector<Real> param
1049 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1051 fem_->apply_adjoint_control_pde_hessian(*ahwvp,*wp,*vp,*up,*zp);
1060 ROL::Ptr<std::vector<Real> > ahwvp =
1062 ROL::Ptr<const std::vector<Real> > wp =
1064 ROL::Ptr<const std::vector<Real> > vp =
1066 ROL::Ptr<const std::vector<Real> > up =
1068 ROL::Ptr<const std::vector<Real> > zp =
1071 const std::vector<Real> param
1073 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1075 fem_->apply_adjoint_pde_control_hessian(*ahwvp,*wp,*vp,*up,*zp);
1084 ROL::Ptr<std::vector<Real> > ahwvp =
1086 ROL::Ptr<const std::vector<Real> > wp =
1088 ROL::Ptr<const std::vector<Real> > vp =
1090 ROL::Ptr<const std::vector<Real> > up =
1092 ROL::Ptr<const std::vector<Real> > zp =
1095 const std::vector<Real> param
1097 fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1099 fem_->apply_adjoint_control_hessian(*ahwvp,*wp,*vp,*up,*zp);
1107 template<
class Real>
1118 ROL::Ptr<BurgersFEM<Real> >
fem_;
1119 ROL::Ptr<ROL::Vector<Real> >
ud_;
1120 ROL::Ptr<ROL::Vector<Real> >
diff_;
1132 ROL::Ptr<const std::vector<Real> > up =
1134 ROL::Ptr<const std::vector<Real> > zp =
1136 ROL::Ptr<const std::vector<Real> > udp =
1139 std::vector<Real> diff(udp->size(),0.0);
1140 for (
unsigned i = 0; i < udp->size(); i++) {
1141 diff[i] = (*up)[i] - (*udp)[i];
1143 return 0.5*(
fem_->compute_L2_dot(diff,diff) +
alpha_*
fem_->compute_L2_dot(*zp,*zp));
1147 ROL::Ptr<std::vector<Real> > gp =
1149 ROL::Ptr<const std::vector<Real> > up =
1151 ROL::Ptr<const std::vector<Real> > udp =
1154 std::vector<Real> diff(udp->size(),0.0);
1155 for (
unsigned i = 0; i < udp->size(); i++) {
1156 diff[i] = (*up)[i] - (*udp)[i];
1158 fem_->apply_mass(*gp,diff);
1162 ROL::Ptr<std::vector<Real> > gp =
1164 ROL::Ptr<const std::vector<Real> > zp =
1167 fem_->apply_mass(*gp,*zp);
1168 for (
unsigned i = 0; i < zp->size(); i++) {
1175 ROL::Ptr<std::vector<Real> > hvp =
1177 ROL::Ptr<const std::vector<Real> > vp =
1180 fem_->apply_mass(*hvp,*vp);
1195 ROL::Ptr<std::vector<Real> > hvp =
1197 ROL::Ptr<const std::vector<Real> > vp =
1200 fem_->apply_mass(*hvp,*vp);
1201 for (
unsigned i = 0; i < vp->size(); i++) {
1207 template<
class Real,
class Ordinal>
1215 catch (std::exception &e) {
1222 : ROL::TeuchosBatchManager<Real,Ordinal>(comm) {}
1224 ROL::Ptr<std::vector<Real> > input_ptr;
1226 int dim_i = input_ptr->size();
1227 ROL::Ptr<std::vector<Real> > output_ptr;
1229 int dim_o = output_ptr->size();
1230 if ( dim_i != dim_o ) {
1231 std::cout <<
"L2VectorBatchManager: DIMENSION MISMATCH ON RANK "
1232 << ROL::TeuchosBatchManager<Real,Ordinal>::batchID() <<
"\n";
1235 ROL::TeuchosBatchManager<Real,Ordinal>::sumAll(&(*input_ptr)[0],&(*output_ptr)[0],dim_i);
1240 template<
class Real,
class Ordinal>
1248 catch (std::exception &e) {
1255 : ROL::TeuchosBatchManager<Real,Ordinal>(comm) {}
1257 ROL::Ptr<std::vector<Real> > input_ptr;
1259 int dim_i = input_ptr->size();
1260 ROL::Ptr<std::vector<Real> > output_ptr;
1262 int dim_o = output_ptr->size();
1263 if ( dim_i != dim_o ) {
1264 std::cout <<
"H1VectorBatchManager: DIMENSION MISMATCH ON RANK "
1265 << ROL::TeuchosBatchManager<Real,Ordinal>::batchID() <<
"\n";
1268 ROL::TeuchosBatchManager<Real,Ordinal>::sumAll(&(*input_ptr)[0],&(*output_ptr)[0],dim_i);
1273 template<
class Real>
1274 Real
random(
const ROL::Ptr<
const Teuchos::Comm<int> > &comm) {
1276 if ( Teuchos::rank<int>(*comm)==0 ) {
1277 val = (Real)rand()/(Real)RAND_MAX;
1279 Teuchos::broadcast<int,Real>(*comm,0,1,&val);
H1VectorPrimal< Real > DualConstraintVector
BurgersFEM(int nx=128, Real nl=1.0, Real cH1=1.0, Real cL2=1.0)
L2VectorDual(const ROL::Ptr< std::vector< Real > > &vec, const ROL::Ptr< BurgersFEM< Real > > &fem)
ROL::Ptr< std::vector< Real > > vec_
Provides the interface to evaluate simulation-based objective functions.
ROL::Ptr< ROL::Vector< Real > > diff_
Real compute_H1_dot(const std::vector< Real > &x, const std::vector< Real > &y) const
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 norm() const
Returns where .
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u) const
void axpy(std::vector< Real > &out, const Real a, const std::vector< Real > &x, const std::vector< Real > &y) const
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
void apply_adjoint_pde_control_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
ROL::Ptr< std::vector< Real > > vec_
int dimension() const
Return dimension of the vector space.
ROL::Ptr< L2VectorDual< Real > > dual_vec_
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 plus(const ROL::Vector< Real > &x)
Compute , where .
void apply_adjoint_pde_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
ROL::Ptr< BurgersFEM< Real > > fem_
H1VectorPrimal(const ROL::Ptr< std::vector< Real > > &vec, const ROL::Ptr< BurgersFEM< Real > > &fem)
Real norm() const
Returns where .
Contains definitions of custom data types in ROL.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
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 ...
L2VectorBatchManager(const ROL::Ptr< const Teuchos::Comm< Ordinal > > &comm)
ROL::Ptr< const std::vector< Real > > getVector() const
void plus(const ROL::Vector< Real > &x)
Compute , where .
ROL::Ptr< BurgersFEM< Real > > fem_
const std::vector< Real > getParameter(void) const
Real compute_H1_norm(const std::vector< Real > &r) const
Real norm() const
Returns where .
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 ...
ROL::Ptr< std::vector< Real > > vec_
H1VectorDual(const ROL::Ptr< std::vector< Real > > &vec, const ROL::Ptr< BurgersFEM< Real > > &fem)
virtual void zero()
Set to zero vector.
ROL::Ptr< std::vector< Real > > getVector()
Defines the linear algebra or vector space interface.
void cast_vector(ROL::Ptr< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
ROL::Ptr< BurgersFEM< Real > > fem_
H1VectorDual< Real > DualStateVector
H1VectorDual< Real > PrimalConstraintVector
void scale(const Real alpha)
Compute where .
L2VectorPrimal(const ROL::Ptr< std::vector< Real > > &vec, const ROL::Ptr< BurgersFEM< Real > > &fem)
Real dot(const ROL::Vector< Real > &x) const
Compute where .
H1VectorPrimal< Real > PrimalStateVector
void hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Constraint_BurgersControl(ROL::Ptr< BurgersFEM< Real > > &fem, bool useHessian=true)
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.
void plus(const ROL::Vector< Real > &x)
Compute , where .
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 mesh_spacing(void) const
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
int dimension() const
Return dimension of the vector space.
ROL::Ptr< std::vector< Real > > vec_
void test_inverse_mass(std::ostream &outStream=std::cout)
void apply_pde_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
ROL::Ptr< BurgersFEM< Real > > fem_
ROL::Ptr< ROL::Vector< Real > > ud_
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
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...
ROL::Ptr< const std::vector< Real > > getVector() const
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
ROL::Ptr< const std::vector< Real > > getVector() const
void set(const ROL::Vector< Real > &x)
Set where .
void scale(const Real alpha)
Compute where .
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u) const
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
H1VectorDual< Real > DualStateVector
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 .
Real compute_L2_dot(const std::vector< Real > &x, const std::vector< Real > &y) const
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 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) const
ROL::Ptr< std::vector< Real > > getVector()
L2VectorPrimal< Real > PrimalControlVector
void apply_inverse_mass(std::vector< Real > &Mu, const std::vector< Real > &u) const
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
ROL::Ptr< std::vector< Real > > getVector()
ROL::Ptr< BurgersFEM< Real > > fem_
ROL::Ptr< const std::vector< Real > > getVector() const
ROL::Ptr< BurgersFEM< Real > > fem_
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0) const
void scale(const Real alpha)
Compute where .
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 random(const ROL::Ptr< const Teuchos::Comm< int > > &comm)
void set(const ROL::Vector< Real > &x)
Set where .
void apply_adjoint_control_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
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) const
H1VectorBatchManager(const ROL::Ptr< const Teuchos::Comm< Ordinal > > &comm)
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void set_problem_data(const Real nu, const Real f, const Real u0, const Real u1)
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_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
void test_inverse_H1(std::ostream &outStream=std::cout)
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void apply_adjoint_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Real dot(const ROL::Vector< Real > &x) const
Compute where .
void cast_vector(ROL::Ptr< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
L2VectorPrimal< Real > PrimalControlVector
H1VectorPrimal< Real > PrimalStateVector
void apply_inverse_pde_jacobian(std::vector< Real > &ijv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
L2VectorDual< Real > DualControlVector
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
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 apply_H1(std::vector< Real > &Mu, const std::vector< Real > &u) const
void plus(const ROL::Vector< Real > &x)
Compute , where .
void sumAll(ROL::Vector< Real > &input, ROL::Vector< Real > &output)
void apply_inverse_H1(std::vector< Real > &Mu, const std::vector< Real > &u) const
void apply_adjoint_pde_jacobian(std::vector< Real > &ajv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
void set(const ROL::Vector< Real > &x)
Set where .
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void apply_inverse_adjoint_pde_jacobian(std::vector< Real > &iajv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
void apply_adjoint_control_pde_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
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...
int dimension() const
Return dimension of the vector space.
Defines the constraint operator interface for simulation-based optimization.
L2VectorDual< Real > DualControlVector
ROL::Ptr< H1VectorPrimal< Real > > dual_vec_
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
void sumAll(ROL::Vector< Real > &input, ROL::Vector< Real > &output)
void set(const ROL::Vector< Real > &x)
Set where .
Real norm() const
Returns where .
void scale(std::vector< Real > &u, const Real alpha=0.0) const
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ROL::Ptr< H1VectorDual< Real > > dual_vec_
Objective_BurgersControl(const ROL::Ptr< BurgersFEM< Real > > &fem, const ROL::Ptr< ROL::Vector< Real > > &ud, Real alpha=1.e-4)
ROL::Ptr< std::vector< Real > > getVector()
void scale(const Real alpha)
Compute where .
int dimension() const
Return dimension of the vector space.
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Real compute_L2_norm(const std::vector< Real > &r) const
ROL::Ptr< L2VectorPrimal< Real > > dual_vec_