5 #include "Teuchos_LAPACK.hpp" 
    6 #include "Teuchos_SerialDenseMatrix.hpp" 
   25   FEM(
int nx = 128, Real kl = 0.1, Real kr = 10.0) : 
nx_(nx), 
kl_(kl), 
kr_(kr) {
 
   29     pts_[0] = -0.9061798459386640; 
 
   30     pts_[1] = -0.5384693101056831;
 
   32     pts_[3] =  0.5384693101056831;
 
   33     pts_[4] =  0.9061798459386640;
 
   36     wts_[0] = 0.2369268850561891;
 
   37     wts_[1] = 0.4786286704993665;
 
   38     wts_[2] = 0.5688888888888889;
 
   39     wts_[3] = 0.4786286704993665;
 
   40     wts_[4] = 0.2369268850561891;
 
   43     for (
int i = 0; i < 5; i++) {
 
   47     for (
int i = 0; i < 5; i++) {
 
   55   void build_mesh(std::vector<Real> &x, 
const std::vector<Real> ¶m) {
 
   59     Real  xp = 0.1*param[0];
 
   75     Real dx1 = (xp+1.0)/(Real)nx1;
 
   76     Real dx2 = (1.0-xp)/(Real)nx2;
 
   79     x.resize(n1+n2+3,0.0);
 
   80     for (
int k = 0; k < n1+n2+3; k++) {
 
   81       x[k] = ((k < nx1) ? (Real)k*dx1-1.0 : ((k==nx1) ? xp : (Real)(k-nx1)*dx2+xp));
 
   93     for ( 
int i = 0; i < rem; i++ ) {
 
   94       if ( i%4 == 0 )      { nz3++; }
 
   95       else if ( i%4 == 1 ) { nz1++; }
 
   96       else if ( i%4 == 2 ) { nz4++; } 
 
   97       else if ( i%4 == 3 ) { nz2++; }
 
  101     for (
int k = 0; k < 
nz(); k++) {
 
  103         x[k] = 0.25*(Real)k/(Real)(nz1-1) - 1.0;
 
  105       if ( k >= nz1 && k < nz1+nz2 ) {
 
  106         x[k] = 0.5*(Real)(k-nz1+1)/(Real)nz2 - 0.75;
 
  108       if ( k >= nz1+nz2 && k < nz1+nz2+nz3 ) {
 
  109         x[k] = 0.5*(Real)(k-nz1-nz2+1)/(Real)nz3 - 0.25;
 
  111       if ( k >= nz1+nz2+nz3-1 ) {
 
  112         x[k] = 0.75*(Real)(k-nz1-nz2-nz3+1)/(Real)nz4 + 0.25;
 
  125   void build_force(std::vector<Real> &F, 
const std::vector<Real> ¶m) {
 
  127     std::vector<Real> x(
nu()+2,0.0);
 
  131     F.resize(x.size()-2,0.0);
 
  133     int size = 
static_cast<int>(x.size());
 
  134     for (
int i = 0; i < size-2; i++) {
 
  135       for (
int j = 0; j < 5; j++) {
 
  137         pt = (x[i+1]-x[i])*
pts_[j] + x[i];
 
  138         F[i] += 
wts_[j]*(pt-x[i])*std::exp(-std::pow(pt-0.5*param[1],2.0));
 
  140         pt = (x[i+2]-x[i+1])*
pts_[j] + x[i+1];
 
  141         F[i] += 
wts_[j]*(x[i+2]-pt)*std::exp(-std::pow(pt-0.5*param[1],2.0));
 
  148                   const std::vector<Real> ¶m) {
 
  150     std::vector<Real> x(
nu()+2,0.0);
 
  153     int xsize = 
static_cast<int>(x.size());
 
  155     d.resize(xsize-2,0.0);
 
  156     for ( 
int i = 0; i < xsize-2; i++ ) {
 
  157       if ( x[i+1] < 0.1*param[0] ) {
 
  158         d[i] = 
kl_/(x[i+1]-x[i]) + 
kl_/(x[i+2]-x[i+1]);
 
  160       else if ( x[i+1] > 0.1*param[0] ) {
 
  161         d[i] = 
kr_/(x[i+1]-x[i]) + 
kr_/(x[i+2]-x[i+1]);
 
  164         d[i] = 
kl_/(x[i+1]-x[i]) + 
kr_/(x[i+2]-x[i+1]);
 
  170     dl.resize(xsize-3,0.0);
 
  171     for ( 
int i = 0; i < xsize-3; i++ ) {
 
  172       if ( x[i+2] <= 0.1*param[0] ) {
 
  173         dl[i] = -
kl_/(x[i+2]-x[i+1]);
 
  175       else if ( x[i+2] > 0.1*param[0] ) {
 
  176         dl[i] = -
kr_/(x[i+2]-x[i+1]);
 
  181     du.assign(dl.begin(),dl.end());
 
  185   void apply_jacobian_1(std::vector<Real> &jv, 
const std::vector<Real> &v, 
const std::vector<Real> ¶m) {
 
  187     std::vector<Real> dl(
nu()-1,0.0);
 
  188     std::vector<Real> d(
nu(),0.0);
 
  189     std::vector<Real> du(
nu()-1,0.0);
 
  193     jv.resize(d.size(),0.0);
 
  194     int size = 
static_cast<int>(d.size());
 
  195     for ( 
int i = 0; i < size; ++i ) {
 
  197       jv[i] += ((i>0) ? dl[i-1]*v[i-1] : 0.0);
 
  198       jv[i] += ((i<size-1) ? du[i]*v[i+1] : 0.0);
 
  202   void apply_jacobian_2(std::vector<Real> &jv, 
const std::vector<Real> &v, 
const std::vector<Real> ¶m) {
 
  204     std::vector<Real> xu(
nu()+2,0.0);
 
  207     std::vector<Real> xz(
nz(),0.0);
 
  210     std::vector<Real> x(xu.size()+xz.size(),0.0);
 
  211     typename std::vector<Real>::iterator it;
 
  212     it = std::set_union(xu.begin(),xu.end(),xz.begin(),xz.end(),x.begin());
 
  213     x.resize(it-x.begin());
 
  215     std::vector<Real> vi;
 
  216     int xzsize = 
static_cast<int>(xz.size());
 
  217     for (it=x.begin(); it!=x.end(); it++) {
 
  218       for (
int i = 0; i < xzsize-1; i++) {
 
  219         if ( *it >= xz[i] && *it <= xz[i+1] ) {
 
  220           vi.push_back(v[i+1]*(*it-xz[i])/(xz[i+1]-xz[i]) + v[i]*(xz[i+1]-*it)/(xz[i+1]-xz[i]));
 
  225     int xsize = 
static_cast<int>(x.size());
 
  227     std::vector<Real> Mv(xsize,0.0);
 
  228     for (
int i = 0; i < xsize; i++) {
 
  229       Mv[i] += ((i>0) ? (x[i]-x[i-1])/6.0*vi[i-1] + (x[i]-x[i-1])/3.0*vi[i] : 0.0);
 
  230       Mv[i] += ((i<xsize-1) ? (x[i+1]-x[i])/3.0*vi[i] + (x[i+1]-x[i])/6.0*vi[i+1] : 0.0);
 
  233     int xusize = 
static_cast<int>(xu.size());
 
  235     jv.resize(xusize-2,0.0);
 
  236     for (
int i = 0; i < xusize-2; i++) {
 
  237       for (
int j = 0; j < xsize-1; j++) {
 
  238         jv[i] += (((x[j]>=xu[i  ]) && (x[j]< xu[i+1])) ? Mv[j]*(x[j]-xu[i  ])/(xu[i+1]-xu[i  ]) : 0.0);
 
  239         jv[i] += (((x[j]>=xu[i+1]) && (x[j]<=xu[i+2])) ? Mv[j]*(xu[i+2]-x[j])/(xu[i+2]-xu[i+1]) : 0.0);
 
  240         if ( x[j] > xu[i+2] ) { 
break; }
 
  246                           const std::vector<Real> ¶m){
 
  248     std::vector<Real> xu(
nu()+2,0.0);
 
  251     std::vector<Real> xz(
nz(),0.0);
 
  254     std::vector<Real> x(xu.size()+xz.size(),0.0);
 
  255     typename std::vector<Real>::iterator it;
 
  256     it = std::set_union(xu.begin(),xu.end(),xz.begin(),xz.end(),x.begin());
 
  257     x.resize(it-x.begin());
 
  259     std::vector<Real> vi;
 
  261     int xusize = 
static_cast<int>(xu.size());
 
  262     for (it=x.begin(); it!=x.end(); it++) {
 
  263       for (
int i = 0; i < xusize-1; i++) {
 
  264         if ( *it >= xu[i] && *it <= xu[i+1] ) {
 
  266           val += ((i!=0       ) ? v[i-1]*(xu[i+1]-*it)/(xu[i+1]-xu[i]) : 0.0);
 
  267           val += ((i!=xusize-2) ? v[i  ]*(*it-xu[i  ])/(xu[i+1]-xu[i]) : 0.0);
 
  274     int xsize = 
static_cast<int>(x.size());
 
  275     std::vector<Real> Mv(xsize,0.0);
 
  276     for (
int i = 0; i < xsize; i++) {
 
  277       Mv[i] += ((i>0) ? (x[i]-x[i-1])/6.0*vi[i-1] + (x[i]-x[i-1])/3.0*vi[i] : 0.0);
 
  278       Mv[i] += ((i<xsize-1) ? (x[i+1]-x[i])/3.0*vi[i] + (x[i+1]-x[i])/6.0*vi[i+1] : 0.0);
 
  283     int xzsize = 
static_cast<int>(xz.size());
 
  284     for (
int i = 0; i < xzsize; i++) {
 
  285       for (
int j = 0; j < xsize; j++) {
 
  287           jv[i] += (((x[j]>=xz[i  ]) && (x[j]<=xz[i+1])) ? Mv[j]*(xz[i+1]-x[j])/(xz[i+1]-xz[i  ]) : 0.0);
 
  288           if ( x[j] > xz[i+1] ) { 
break; }
 
  290         else if ( i == xzsize-1 ) {
 
  291           jv[i] += (((x[j]>=xz[i-1]) && (x[j]<=xz[i  ])) ? Mv[j]*(x[j]-xz[i-1])/(xz[i  ]-xz[i-1]) : 0.0);
 
  294           jv[i] += (((x[j]>=xz[i-1]) && (x[j]< xz[i  ])) ? Mv[j]*(x[j]-xz[i-1])/(xz[i  ]-xz[i-1]) : 0.0);
 
  295           jv[i] += (((x[j]>=xz[i  ]) && (x[j]<=xz[i+1])) ? Mv[j]*(xz[i+1]-x[j])/(xz[i+1]-xz[i  ]) : 0.0);
 
  296           if ( x[j] > xz[i+1] ) { 
break; }
 
  302   void apply_mass_state(std::vector<Real> &Mv, 
const std::vector<Real> &v, 
const std::vector<Real> ¶m) {
 
  307     int size = 
static_cast<int>(x.size());
 
  309     Mv.resize(size-2,0.0);
 
  310     for (
int i = 0; i < size-2; i++) {
 
  311       Mv[i] = (x[i+2]-x[i])/3.0*v[i];
 
  313         Mv[i] += (x[i+1]-x[i])/6.0*v[i-1];
 
  316         Mv[i] += (x[i+2]-x[i+1])/6.0*v[i+1];
 
  323     std::vector<Real> x(
nz(),0.0);
 
  326     int xsize = 
static_cast<Real
>(x.size());
 
  328     Mv.resize(xsize,0.0);
 
  329     for (
int i = 0; i < xsize; i++) {
 
  331         Mv[i] += (x[i]-x[i-1])/6.0*v[i-1] + (x[i]-x[i-1])/3.0*v[i];
 
  334         Mv[i] += (x[i+1]-x[i])/3.0*v[i] + (x[i+1]-x[i])/6.0*v[i+1];
 
  341     std::vector<Real> x(
nz(),0.0);
 
  344     std::vector<Real> d(
nz(),0.0);
 
  345     std::vector<Real> dl(
nz()-1,0.0);
 
  346     std::vector<Real> du(
nz()-1,0.0);
 
  347     for (
int i = 0; i < 
nz(); i++) {
 
  349         dl[i-1] = (x[i]-x[i-1])/6.0;
 
  350         d[i]   += (x[i]-x[i-1])/3.0;
 
  353         d[i] += (x[i+1]-x[i])/3.0;
 
  354         du[i] = (x[i+1]-x[i])/6.0;
 
  365   void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, 
 
  366                     std::vector<Real> &du, 
const std::vector<Real> &r, 
const bool transpose = 
false) {
 
  368     u.assign(r.begin(),r.end());
 
  370     Teuchos::LAPACK<int,Real> lp;
 
  371     std::vector<Real> du2(r.size()-2,0.0);
 
  372     std::vector<int> ipiv(r.size(),0);
 
  377     lp.GTTRF(n,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
 
  382     lp.GTTRS(trans,n,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
 
  393       case 1:  val = ((x[i]<0.5) ? 1.0 : 0.0);             
break;
 
  394       case 2:  val = 1.0;                                  
break;
 
  395       case 3:  val = std::abs(std::sin(8.0*M_PI*x[i]));    
break;
 
  396       case 4:  val = std::exp(-0.5*(x[i]-0.5)*(x[i]-0.5)); 
break;
 
  405   const ROL::Ptr<FEM<Real> > 
FEM_;
 
  413   void plus(std::vector<Real> &u, 
const std::vector<Real> &s, 
const Real alpha=1.0) {
 
  414     int size = 
static_cast<int>(u.size());
 
  415     for (
int i=0; i<size; i++) {
 
  421   void scale(std::vector<Real> &u, 
const Real alpha=0.0) {
 
  422     int size = 
static_cast<int>(u.size());
 
  423     for (
int i=0; i<size; i++) {
 
  440     ROL::Ptr<const std::vector<Real> > up = 
dynamic_cast<const ROL::StdVector<Real>&
>(u).getVector();
 
  441     ROL::Ptr<const std::vector<Real> > zp = 
dynamic_cast<const ROL::StdVector<Real>&
>(z).getVector();
 
  445     std::vector<Real> F(
FEM_->nu(),0.0);
 
  447     std::vector<Real> Bz(
FEM_->nu(),0.0);
 
  457     ROL::Ptr<const std::vector<Real> > zp = 
dynamic_cast<const ROL::StdVector<Real>&
>(z).getVector();
 
  459     std::vector<Real> dl(
FEM_->nu()-1,0.0);
 
  460     std::vector<Real> d(
FEM_->nu(),0.0);
 
  461     std::vector<Real> du(
FEM_->nu()-1,0.0);
 
  464     std::vector<Real> F(
FEM_->nu(),0.0);
 
  466     std::vector<Real> Bz(
FEM_->nu(),0.0);
 
  470     FEM_->linear_solve(*up,dl,d,du,F);
 
  479     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  486     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  494     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  496     std::vector<Real> dl(
FEM_->nu()-1,0.0);
 
  497     std::vector<Real> d(
FEM_->nu(),0.0);
 
  498     std::vector<Real> du(
FEM_->nu()-1,0.0);
 
  501     FEM_->linear_solve(*jvp,dl,d,du,*vp);
 
  513     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  549   const ROL::Ptr<FEM<Real> > 
FEM_;
 
  556   Real 
dot(
const std::vector<Real> &x, 
const std::vector<Real> &y) {
 
  558     int size = 
static_cast<int>(x.size());
 
  559     for (
int i=0; i<size; i++) {
 
  566   void scale(std::vector<Real> &u, 
const Real alpha=0.0) {
 
  567     int size = 
static_cast<int>(u.size());
 
  568     for (
int i=0; i<size; i++) {
 
  581     ROL::Ptr<const std::vector<Real> > up = 
dynamic_cast<const ROL::StdVector<Real>&
>(u).getVector();
 
  582     ROL::Ptr<const std::vector<Real> > zp = 
dynamic_cast<const ROL::StdVector<Real>&
>(z).getVector();
 
  583     int nz = 
FEM_->nz(), nu = 
FEM_->nu();
 
  585     std::vector<Real> Mz(nz);
 
  586     FEM_->apply_mass_control(Mz,*zp);
 
  589     std::vector<Real> diff(nu);
 
  590     for (
int i=0; i<nu; i++) {
 
  593     std::vector<Real> Mu(nu);
 
  595     Real uval = 0.5*
dot(Mu,diff);  
 
  601     ROL::Ptr<const std::vector<Real> > up = 
dynamic_cast<const ROL::StdVector<Real>&
>(u).getVector();
 
  604     std::vector<Real> diff(nu);
 
  605     for (
int i=0; i<nu; i++) {
 
  613     ROL::Ptr<const std::vector<Real> > zp = 
dynamic_cast<const ROL::StdVector<Real>&
>(z).getVector();
 
  615     FEM_->apply_mass_control(*gp,*zp);
 
  622     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  640     ROL::Ptr<const std::vector<Real> > vp = 
dynamic_cast<const ROL::StdVector<Real>&
>(v).getVector();
 
  642     FEM_->apply_mass_control(*hvp,*vp);
 
Provides the interface to evaluate simulation-based objective functions. 
void scale(std::vector< Real > &u, const Real alpha=0.0)
void plus(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.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. 
void build_mesh(std::vector< Real > &x)
void apply_adjoint_jacobian_2(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > ¶m)
Contains definitions of custom data types in ROL. 
DiffusionConstraint(const ROL::Ptr< FEM< Real > > &FEM)
const std::vector< Real > getParameter(void) const 
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator  at . 
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 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_jacobian_2(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > ¶m)
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
virtual void zero()
Set to zero vector. 
Defines the linear algebra or vector space interface. 
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 build_mesh(std::vector< Real > &x, const std::vector< Real > ¶m)
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 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. 
int getNumSolves(void) const 
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...
FEM(int nx=128, Real kl=0.1, Real kr=10.0)
void applyInverseAdjointJacobian_1(ROL::Vector< Real > &jv, 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 ...
DiffusionObjective(const ROL::Ptr< FEM< Real > > fem, const Real alpha=1.e-4)
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 applyInverseJacobian_1(ROL::Vector< Real > &jv, 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 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_control(std::vector< Real > &Mv, const std::vector< Real > &v)
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 scale(std::vector< Real > &u, const Real alpha=0.0)
void apply_jacobian_1(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > ¶m)
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
const std::vector< Real > getParameter(void) const 
void solve(ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve  for . 
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. 
Defines the constraint operator interface for simulation-based optimization. 
const ROL::Ptr< FEM< Real > > FEM_
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 build_jacobian_1(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > ¶m)
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 apply_mass_state(std::vector< Real > &Mv, const std::vector< Real > &v, const std::vector< Real > ¶m)
void apply_inverse_mass_control(std::vector< Real > &Mv, const std::vector< Real > &v)
void applyAdjointJacobian_1(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 the vector . This is the primary inter...
const ROL::Ptr< FEM< Real > > FEM_
void build_force(std::vector< Real > &F, const std::vector< Real > ¶m)
Real evaluate_target(int i, const std::vector< Real > ¶m)