84 #include "Teuchos_oblackholestream.hpp"
85 #include "Teuchos_GlobalMPISession.hpp"
86 #include "Teuchos_XMLParameterListHelpers.hpp"
99 template <
class Real,
class Element=Real>
102 template <
class Real,
class Element=Real>
105 template <
class Real,
class Element=Real>
108 template <
class Real,
class Element=Real>
115 template <
class Real,
class Element>
119 Teuchos::RCP<std::vector<Element> > std_vec_;
120 mutable Teuchos::RCP<OptDualStdVector<Real> > dual_vec_;
122 Teuchos::RCP<FiniteDifference<Real> >
fd_;
128 std_vec_(std_vec), dual_vec_(Teuchos::null), fd_(fd) {}
132 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
133 unsigned dimension = std_vec_->size();
134 for (
unsigned i=0; i<dimension; i++) {
135 (*std_vec_)[i] += (*xvalptr)[i];
140 unsigned dimension = std_vec_->size();
141 for (
unsigned i=0; i<dimension; i++) {
142 (*std_vec_)[i] *= alpha;
151 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
153 Teuchos::RCP<std::vector<Real> > kxvalptr = Teuchos::rcp(
new std::vector<Real> (std_vec_->size(), 0.0) );
155 this->fd_->apply(xvalptr,kxvalptr);
157 unsigned dimension = std_vec_->size();
158 for (
unsigned i=0; i<dimension; i++) {
159 val += (*std_vec_)[i]*(*kxvalptr)[i];
166 val = std::sqrt( dot(*
this) );
170 Teuchos::RCP<Vector<Real> >
clone()
const {
171 return Teuchos::rcp(
new OptStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size()) ),fd_ ) );
174 Teuchos::RCP<const std::vector<Element> >
getVector()
const {
178 Teuchos::RCP<Vector<Real> >
basis(
const int i )
const {
179 Teuchos::RCP<OptStdVector> e = Teuchos::rcp(
new OptStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size(), 0.0)), fd_ ) );
180 (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
189 Teuchos::RCP<std::vector<Element> > dual_vecp = Teuchos::rcp(
new std::vector<Element>(*std_vec_));
191 this->fd_->apply(dual_vecp);
199 template <
class Real,
class Element>
203 Teuchos::RCP<std::vector<Element> > std_vec_;
204 mutable Teuchos::RCP<OptStdVector<Real> > dual_vec_;
205 Teuchos::RCP<FiniteDifference<Real> >
fd_;
210 std_vec_(std_vec), dual_vec_(Teuchos::null), fd_(fd) {}
214 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
215 unsigned dimension = std_vec_->size();
216 for (
unsigned i=0; i<dimension; i++) {
217 (*std_vec_)[i] += (*xvalptr)[i];
222 unsigned dimension = std_vec_->size();
223 for (
unsigned i=0; i<dimension; i++) {
224 (*std_vec_)[i] *= alpha;
231 Teuchos::RCP<const std::vector<Element> > kxvalptr = ex.
getVector();
232 Teuchos::RCP<std::vector<Real> > xvalptr = Teuchos::rcp(
new std::vector<Real> (std_vec_->size(), 0.0) );
233 this->fd_->solve(kxvalptr,xvalptr);
235 unsigned dimension = std_vec_->size();
236 for (
unsigned i=0; i<dimension; i++) {
237 val += (*std_vec_)[i]*(*xvalptr)[i];
244 val = std::sqrt( dot(*
this) );
248 Teuchos::RCP<Vector<Real> >
clone()
const {
249 return Teuchos::rcp(
new OptDualStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size()) ), fd_ ) );
252 Teuchos::RCP<const std::vector<Element> >
getVector()
const {
256 Teuchos::RCP<Vector<Real> >
basis(
const int i )
const {
257 Teuchos::RCP<OptDualStdVector> e = Teuchos::rcp(
new OptDualStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size(), 0.0)),fd_ ) );
258 (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
265 Teuchos::RCP<std::vector<Element> > dual_vecp = Teuchos::rcp(
new std::vector<Element>(*std_vec_));
268 this->fd_->solve(dual_vecp);
278 template <
class Real,
class Element>
282 Teuchos::RCP<std::vector<Element> > std_vec_;
283 mutable Teuchos::RCP<ConDualStdVector<Real> > dual_vec_;
287 ConStdVector(
const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
291 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
292 unsigned dimension = std_vec_->size();
293 for (
unsigned i=0; i<dimension; i++) {
294 (*std_vec_)[i] += (*xvalptr)[i];
299 unsigned dimension = std_vec_->size();
300 for (
unsigned i=0; i<dimension; i++) {
301 (*std_vec_)[i] *= alpha;
308 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
312 unsigned dimension = std_vec_->size();
313 for (
unsigned i=0; i<dimension; i++) {
314 val += (*std_vec_)[i]*(*xvalptr)[i];
321 val = std::sqrt( dot(*
this) );
325 Teuchos::RCP<Vector<Real> >
clone()
const {
326 return Teuchos::rcp(
new ConStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size())) ) );
329 Teuchos::RCP<const std::vector<Element> >
getVector()
const {
333 Teuchos::RCP<Vector<Real> >
basis(
const int i )
const {
334 Teuchos::RCP<ConStdVector> e = Teuchos::rcp(
new ConStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size(), 0.0)) ) );
335 (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
342 dual_vec_ = Teuchos::rcp(
new ConDualStdVector<Real>( Teuchos::rcp(
new std::vector<Element>(*std_vec_) ) ) );
350 template <
class Real,
class Element>
354 Teuchos::RCP<std::vector<Element> > std_vec_;
355 mutable Teuchos::RCP<ConStdVector<Real> > dual_vec_;
359 ConDualStdVector(
const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
363 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
364 unsigned dimension = std_vec_->size();
365 for (
unsigned i=0; i<dimension; i++) {
366 (*std_vec_)[i] += (*xvalptr)[i];
371 unsigned dimension = std_vec_->size();
372 for (
unsigned i=0; i<dimension; i++) {
373 (*std_vec_)[i] *= alpha;
380 Teuchos::RCP<const std::vector<Element> > xvalptr = ex.
getVector();
381 unsigned dimension = std_vec_->size();
382 for (
unsigned i=0; i<dimension; i++) {
383 val += (*std_vec_)[i]*(*xvalptr)[i];
390 val = std::sqrt( dot(*
this) );
394 Teuchos::RCP<Vector<Real> >
clone()
const {
395 return Teuchos::rcp(
new ConDualStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size())) ) );
398 Teuchos::RCP<const std::vector<Element> >
getVector()
const {
402 Teuchos::RCP<Vector<Real> >
basis(
const int i )
const {
403 Teuchos::RCP<ConDualStdVector> e = Teuchos::rcp(
new ConDualStdVector( Teuchos::rcp(
new std::vector<Element>(std_vec_->size(), 0.0)) ) );
404 (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
411 dual_vec_ = Teuchos::rcp(
new ConStdVector<Real>( Teuchos::rcp(
new std::vector<Element>(*std_vec_) ) ) );
422 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real> >
437 Teuchos::RCP<const std::vector<Real> >
Vp_;
439 Teuchos::RCP<FiniteDifference<Real> >
fd_;
447 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(v))).getVector();
450 Teuchos::RCP<std::vector<Real> > kvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(kv)).getVector());
454 (*kvp)[0] = (2.0*(*vp)[0]-(*vp)[1])/dx2;
456 for(
int i=1;i<nx_-1;++i) {
457 (*kvp)[i] = (2.0*(*vp)[i]-(*vp)[i-1]-(*vp)[i+1])/dx2;
460 (*kvp)[nx_-1] = (2.0*(*vp)[nx_-1]-(*vp)[nx_-2])/dx2;
467 Vp_((Teuchos::dyn_cast<
StdVector<Real> >(const_cast<
Vector<Real> &>(V))).getVector()), fd_(fd) {
470 dx_ = (1.0/(1.0+nx_));
481 Teuchos::RCP<const std::vector<Real> > psip =
482 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
487 Teuchos::RCP<std::vector<Real> > kpsip = Teuchos::rcp(
new std::vector<Real> (nx_, 0.0) );
488 XDual kpsi(kpsip,this->fd_);
492 this->applyK(psi,kpsi);
494 for(
int i=0;i<nx_;++i) {
495 J += (*psip)[i]*(*kpsip)[i] + (*Vp_)[i]*pow((*psip)[i],2) + g_*pow((*psip)[i],4);
509 Teuchos::RCP<const std::vector<Real> > psip =
510 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
514 Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(g)).getVector());
517 Teuchos::RCP<std::vector<Real> > kpsip = Teuchos::rcp(
new std::vector<Real> (nx_, 0.0) );
518 XDual kpsi(kpsip,this->fd_);
520 this->applyK(psi,kpsi);
522 for(
int i=0;i<nx_;++i) {
523 (*gp)[i] = ((*kpsip)[i] + (*Vp_)[i]*(*psip)[i] + 2.0*g_*pow((*psip)[i],3))*dx_;
535 Teuchos::RCP<const std::vector<Real> > psip =
536 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
540 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(v))).getVector();
543 Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(hv)).getVector());
547 for(
int i=0;i<nx_;++i) {
549 (*hvp)[i] += ( (*Vp_)[i] + 6.0*g_*pow((*psip)[i],2) )*(*vp)[i]*dx_;
558 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real>,
class CPrim=StdVector<Real>,
class CDual=StdVector<Real> >
564 Teuchos::RCP<FiniteDifference<Real> >
fd_;
569 nx_(n), dx_(dx), fd_(fd), exactsolve_(exactsolve) {}
579 Teuchos::RCP<std::vector<Real> > cp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<CPrim>(c)).getVector());
582 Teuchos::RCP<const std::vector<Real> > psip =
583 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
588 for(
int i=0;i<nx_;++i) {
589 (*cp)[0] += dx_*pow((*psip)[i],2);
599 Teuchos::RCP<std::vector<Real> > jvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<CPrim>(jv)).getVector());
602 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(v))).getVector();
605 Teuchos::RCP<const std::vector<Real> > psip =
606 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
611 for(
int i=0;i<nx_;++i) {
612 (*jvp)[0] += 2.0*dx_*(*psip)[i]*(*vp)[i];
622 Teuchos::RCP<std::vector<Real> > ajvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(ajv)).getVector());
625 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<CDual>(
const_cast<Vector<Real> &
>(v))).getVector();
628 Teuchos::RCP<const std::vector<Real> > psip =
629 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
633 for(
int i=0;i<nx_;++i) {
634 (*ajvp)[i] = 2.0*dx_*(*psip)[i]*(*vp)[0];
646 Teuchos::RCP<std::vector<Real> > ahuvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(ahuv)).getVector());
649 Teuchos::RCP<const std::vector<Real> > up = (Teuchos::dyn_cast<CDual>(
const_cast<Vector<Real> &
>(u))).getVector();
652 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(v))).getVector();
655 Teuchos::RCP<const std::vector<Real> > psip =
656 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
659 for(
int i=0;i<nx_;++i) {
660 (*ahuvp)[i] = 2.0*dx_*(*vp)[i]*(*up)[0];
672 Teuchos::RCP<std::vector<Real> > v1p =
673 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XPrim>(v1)).getVector());
674 Teuchos::RCP<std::vector<Real> > v2p =
675 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<CDual>(v2)).getVector());
676 Teuchos::RCP<const std::vector<Real> > b1p =
677 (Teuchos::dyn_cast<XDual>(
const_cast<Vector<Real> &
>(b1))).getVector();
678 Teuchos::RCP<const std::vector<Real> > b2p =
679 (Teuchos::dyn_cast<CPrim>(
const_cast<Vector<Real> &
>(b2))).getVector();
680 Teuchos::RCP<const std::vector<Real> > psip =
681 (Teuchos::dyn_cast<XPrim>(
const_cast<Vector<Real> &
>(psi))).getVector();
683 Teuchos::RCP<std::vector<Real> > jacp = Teuchos::rcp(
new std::vector<Real> (nx_, 0.0) );
684 Teuchos::RCP<std::vector<Real> > b1dp = Teuchos::rcp(
new std::vector<Real> (nx_, 0.0) );
688 for(
int i=0;i<nx_;++i) {
689 (*jacp)[i] = (*psip)[i];
690 (*b1dp)[i] = (*b1p)[i];
701 Real d = 1.0/jac.dot(jac);
702 Real p = jac.dot(b1d);
704 (*v2p)[0] = d*(p-(*b2p)[0]);
707 v1.
scale(-(*v2p)[0]);
710 return std::vector<Real>(0);
Provides the interface to evaluate objective functions.
Teuchos::RCP< FiniteDifference< Real > > fd_
Real value(const Vector< Real > &psi, Real &tol)
Evaluate .
virtual void scale(const Real alpha)=0
Compute where .
void scale(const Real alpha)
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
ConStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec)
Objective_GrossPitaevskii(const Real &g, const Vector< Real > &V, Teuchos::RCP< FiniteDifference< Real > > fd)
Real norm() const
Returns where .
void value(Vector< Real > &c, const Vector< Real > &psi, Real &tol)
Evaluate .
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Real dot(const Vector< Real > &x) const
Compute where .
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Teuchos::RCP< const std::vector< Element > > getVector() const
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
OptStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec, Teuchos::RCP< FiniteDifference< Real > >fd)
int dimension() const
Return dimension of the vector space.
Teuchos::RCP< FiniteDifference< Real > > fd_
Teuchos::RCP< const std::vector< Element > > getVector() const
Defines the linear algebra or vector space interface.
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Teuchos::RCP< Vector< Real > > basis(const int i) const
Return i-th basis vector.
std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &psi, Real &tol)
Teuchos::RCP< Vector< Real > > basis(const int i) const
Return i-th basis vector.
virtual std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
Approximately solves the augmented system where , , , , is an identity operator, and is a zero operator.
Defines the equality constraint operator interface.
Provides the std::vector implementation of the ROL::Vector interface.
void scale(const Real alpha)
Compute where .
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Real dot(const Vector< Real > &x) const
Compute where .
void scale(const Real alpha)
Compute where .
Teuchos::RCP< const std::vector< Element > > getVector() const
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
int dimension() const
Return dimension of the vector space.
int dimension() const
Return dimension of the vector space.
OptDualStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec, Teuchos::RCP< FiniteDifference< Real > >fd)
Teuchos::RCP< Vector< Real > > basis(const int i) const
Return i-th basis vector.
void applyK(const Vector< Real > &v, Vector< Real > &kv)
Apply finite difference operator.
Real norm() const
Returns where .
Teuchos::RCP< const std::vector< Real > > Vp_
void plus(const Vector< Real > &x)
Compute , where .
int dimension() const
Return dimension of the vector space.
Real dot(const Vector< Real > &x) const
Compute where .
void plus(const Vector< Real > &x)
Compute , where .
const Vector< Real > & dual() const
Modify the dual of vector u to be .
Teuchos::RCP< FiniteDifference< Real > > fd_
Real norm() const
Returns where .
Real dot(const Vector< Real > &x) const
Modify the dot product between primal variables to be .
ConDualStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec)
void gradient(Vector< Real > &g, const Vector< Real > &psi, Real &tol)
Evaluate .
virtual void set(const Vector &x)
Set where .
Normalization_Constraint(int n, Real dx, Teuchos::RCP< FiniteDifference< Real > > fd, bool exactsolve)
Teuchos::RCP< const std::vector< Element > > getVector() const
void plus(const Vector< Real > &x)
Compute , where .
Teuchos::RCP< Vector< Real > > basis(const int i) const
Return i-th basis vector.
void plus(const Vector< Real > &x)
Compute , where .
Real norm() const
Returns where .
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
void scale(const Real alpha)
Compute where .
Teuchos::RCP< FiniteDifference< Real > > fd_