73 #include "Teuchos_GlobalMPISession.hpp"
75 #include "ROL_ParameterList.hpp"
108 ROL::Ptr<const vector>
Vp_;
112 return dynamic_cast<const SV&
>(x).getVector();
117 return dynamic_cast<SV&
>(x).getVector();
127 using namespace Teuchos;
130 ROL::Ptr<const vector> vp = getVector(v);
133 ROL::Ptr<vector> kvp = getVector(kv);
137 (*kvp)[0] = (2.0*(*vp)[0]-(*vp)[1])/dx2;
139 for(
uint i=1;i<nx_-1;++i) {
140 (*kvp)[i] = (2.0*(*vp)[i]-(*vp)[i-1]-(*vp)[i+1])/dx2;
143 (*kvp)[nx_-1] = (2.0*(*vp)[nx_-1]-(*vp)[nx_-2])/dx2;
151 dx_ = (1.0/(1.0+nx_));
164 ROL::Ptr<const vector> psip = getVector(psi);
167 ROL::Ptr<V> kpsi = psi.
clone();
168 ROL::Ptr<vector> kpsip = getVector(*kpsi);
172 this->applyK(psi,*kpsi);
174 for(
uint i=0;i<nx_;++i) {
175 J += (*psip)[i]*(*kpsip)[i] + (*Vp_)[i]*pow((*psip)[i],2) + g_*pow((*psip)[i],4);
191 ROL::Ptr<const vector> psip = getVector(psi);
194 ROL::Ptr<vector> gp = getVector(g);
197 ROL::Ptr<V> kpsi = psi.
clone();
198 ROL::Ptr<vector> kpsip = getVector(*kpsi);
202 for(
uint i=0;i<nx_;++i) {
203 (*gp)[i] = ((*kpsip)[i] + (*Vp_)[i]*(*psip)[i] + 2.0*g_*pow((*psip)[i],3))*dx_;
217 ROL::Ptr<const vector> psip = getVector(psi);
220 ROL::Ptr<const vector> vp = getVector(v);
223 ROL::Ptr<vector> hvp = getVector(hv);
227 for(
uint i=0;i<nx_;++i) {
229 (*hvp)[i] += ( (*Vp_)[i] + 6.0*g_*pow((*psip)[i],2) )*(*vp)[i]*dx_;
252 return dynamic_cast<const SV&
>(x).getVector();
257 return dynamic_cast<SV&
>(x).getVector();
273 ROL::Ptr<vector> cp = getVector(c);
276 ROL::Ptr<const vector> psip = getVector(psi);
279 for(
uint i=0;i<nx_;++i) {
280 (*cp)[0] += dx_*pow((*psip)[i],2);
292 ROL::Ptr<vector> jvp = getVector(jv);
295 ROL::Ptr<const vector> vp = getVector(v);
298 ROL::Ptr<const vector> psip = getVector(psi);
301 for(
uint i=0;i<nx_;++i) {
302 (*jvp)[0] += 2.0*dx_*(*psip)[i]*(*vp)[i];
314 ROL::Ptr<vector> ajvp = getVector(ajv);
317 ROL::Ptr<const vector> vp = getVector(v);
320 ROL::Ptr<const vector> psip = getVector(psi);
322 for(
uint i=0;i<nx_;++i) {
323 (*ajvp)[i] = 2.0*dx_*(*psip)[i]*(*vp)[0];
337 ROL::Ptr<vector> ahuvp = getVector(ahuv);
340 ROL::Ptr<const vector> up = getVector(u);
343 ROL::Ptr<const vector> vp = getVector(v);
346 ROL::Ptr<const vector> psip = getVector(psi);
348 for(
uint i=0;i<nx_;++i) {
349 (*ahuvp)[i] = 2.0*dx_*(*vp)[i]*(*up)[0];
Provides the interface to evaluate objective functions.
typename PV< Real >::size_type size_type
Real value(const Vector< Real > &psi, Real &tol)
Evaluate .
ROL::Ptr< const vector > getVector(const V &x)
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
std::vector< Real > vector
ROL::Ptr< vector > getVector(V &x)
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 .
Normalization_Constraint(int n, Real dx)
Defines the linear algebra or vector space interface.
std::vector< Real > vector
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
ROL::Ptr< const vector > Vp_
void applyK(const Vector< Real > &v, Vector< Real > &kv)
Apply finite difference operator.
ROL::Ptr< const vector > getVector(const V &x)
ROL::Ptr< vector > getVector(V &x)
void gradient(Vector< Real > &g, const Vector< Real > &psi, Real &tol)
Evaluate .
Objective_GrossPitaevskii(const Real &g, const Vector< Real > &V)
Defines the general constraint operator interface.
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .