ROL
gross-pitaevskii/example_02.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
48#include <iostream>
49
50#include "ROL_Stream.hpp"
51#include "Teuchos_GlobalMPISession.hpp"
52
53#include "ROL_ParameterList.hpp"
54#include "ROL_StdVector.hpp"
55#include "ROL_Objective.hpp"
56#include "ROL_Constraint.hpp"
57#include "ROL_Algorithm.hpp"
58#include "ROL_CompositeStep.hpp"
60
62
63
64using namespace ROL;
65
66template <class Real, class Element=Real>
67class OptStdVector; // Optimization space.
68
69template <class Real, class Element=Real>
70class OptDualStdVector; // Dual optimization space.
71
72template <class Real, class Element=Real>
73class ConStdVector; // Constraint space.
74
75template <class Real, class Element=Real>
76class ConDualStdVector; // Dual constraint space.
77
78// Vector space definitions:
79
80
81// Optimization space.
82template <class Real, class Element>
83class OptStdVector : public Vector<Real> {
84
85 typedef std::vector<Element> vector;
86 typedef typename vector::size_type uint;
87
88private:
89ROL::Ptr<std::vector<Element> > std_vec_;
90mutable ROL::Ptr<OptDualStdVector<Real> > dual_vec_;
91
92ROL::Ptr<FiniteDifference<Real> > fd_;
93
94
95public:
96
97OptStdVector(const ROL::Ptr<std::vector<Element> > & std_vec, ROL::Ptr<FiniteDifference<Real> >fd) :
98 std_vec_(std_vec), dual_vec_(ROL::nullPtr), fd_(fd) {}
99
100void plus( const Vector<Real> &x ) {
101 const OptStdVector &ex = dynamic_cast<const OptStdVector&>(x);
102 ROL::Ptr<const vector> xvalptr = ex.getVector();
103 uint dimension = std_vec_->size();
104 for (uint i=0; i<dimension; i++) {
105 (*std_vec_)[i] += (*xvalptr)[i];
106 }
107}
108
109void scale( const Real alpha ) {
110 uint dimension = std_vec_->size();
111 for (uint i=0; i<dimension; i++) {
112 (*std_vec_)[i] *= alpha;
113 }
114}
115
116
118Real dot( const Vector<Real> &x ) const {
119 Real val = 0;
120 const OptStdVector<Real, Element> & ex = dynamic_cast<const OptStdVector&>(x);
121 ROL::Ptr<const vector> xvalptr = ex.getVector();
122
123 ROL::Ptr<vector> kxvalptr = ROL::makePtr<vector>(std_vec_->size(), 0.0);
124
125 fd_->apply(xvalptr,kxvalptr);
126
127 uint dimension = std_vec_->size();
128 for (uint i=0; i<dimension; i++) {
129 val += (*std_vec_)[i]*(*kxvalptr)[i];
130 }
131 return val;
132}
133
134Real norm() const {
135 Real val = 0;
136 val = std::sqrt( dot(*this) );
137 return val;
138}
139
140ROL::Ptr<Vector<Real> > clone() const {
141 return ROL::makePtr<OptStdVector>( ROL::makePtr<vector>(std_vec_->size()), fd_ );
142}
143
144ROL::Ptr<const vector> getVector() const {
145 return std_vec_;
146}
147
148ROL::Ptr<vector> getVector() {
149 return std_vec_;
150}
151
152ROL::Ptr<Vector<Real> > basis( const int i ) const {
153 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(),0.0);
154 ROL::Ptr<OptStdVector> e = ROL::makePtr<OptStdVector>( e_ptr, fd_ );
155 (*e_ptr)[i]= 1.0;
156 return e;
157}
158
159int dimension() const {return static_cast<int>(std_vec_->size());}
160
161
163const Vector<Real> & dual() const {
164 ROL::Ptr<vector> dual_vecp = ROL::makePtr<vector>(*std_vec_);
165 dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( dual_vecp, fd_ );
166 fd_->apply(dual_vecp);
167 return *dual_vec_;
168}
169
170}; // class OptStdVector
171
172
173// Dual optimization space.
174template <class Real, class Element>
175class OptDualStdVector : public Vector<Real> {
176
177 typedef std::vector<Element> vector;
178 typedef typename vector::size_type uint;
179
180private:
181ROL::Ptr<std::vector<Element> > std_vec_;
182mutable ROL::Ptr<OptStdVector<Real> > dual_vec_;
183ROL::Ptr<FiniteDifference<Real> > fd_;
184
185public:
186
187OptDualStdVector(const ROL::Ptr<std::vector<Element> > & std_vec, ROL::Ptr<FiniteDifference<Real> >fd) :
188 std_vec_(std_vec), dual_vec_(ROL::nullPtr), fd_(fd) {}
189
190void plus( const Vector<Real> &x ) {
191 const OptDualStdVector &ex = dynamic_cast<const OptDualStdVector&>(x);
192 ROL::Ptr<const vector> xvalptr = ex.getVector();
193 uint dimension = std_vec_->size();
194 for (uint i=0; i<dimension; i++) {
195 (*std_vec_)[i] += (*xvalptr)[i];
196 }
197}
198
199void scale( const Real alpha ) {
200 uint dimension = std_vec_->size();
201 for (uint i=0; i<dimension; i++) {
202 (*std_vec_)[i] *= alpha;
203 }
204}
205
206Real dot( const Vector<Real> &x ) const {
207 Real val = 0;
208 const OptDualStdVector<Real, Element> & ex = dynamic_cast<const OptDualStdVector<Real, Element>&>(x);
209 ROL::Ptr<const vector> kxvalptr = ex.getVector();
210 ROL::Ptr<vector> xvalptr = ROL::makePtr<vector>(std_vec_->size(), 0.0);
211 fd_->solve(kxvalptr,xvalptr);
212
213 uint dimension = std_vec_->size();
214 for (unsigned i=0; i<dimension; i++) {
215 val += (*std_vec_)[i]*(*xvalptr)[i];
216 }
217 return val;
218}
219
220Real norm() const {
221 Real val = 0;
222 val = std::sqrt( dot(*this) );
223 return val;
224}
225
226ROL::Ptr<Vector<Real> > clone() const {
227 return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()), fd_ );
228}
229
230ROL::Ptr<const std::vector<Element> > getVector() const {
231 return std_vec_;
232}
233
234ROL::Ptr<std::vector<Element> > getVector() {
235 return std_vec_;
236}
237
238ROL::Ptr<Vector<Real> > basis( const int i ) const {
239 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(), 0.0 );
240 ROL::Ptr<OptDualStdVector> e = ROL::makePtr<OptDualStdVector>( e_ptr,fd_ );
241 (*e_ptr)[i] = 1.0;
242 return e;
243}
244
245int dimension() const {return static_cast<int>(std_vec_->size());}
246
247const Vector<Real> & dual() const {
248 ROL::Ptr<vector> dual_vecp = ROL::makePtr<vector>(*std_vec_);
249 dual_vec_ = ROL::makePtr<OptStdVector<Real>>( dual_vecp, fd_ );
250
251 fd_->solve(dual_vecp);
252 return *dual_vec_;
253}
254
255}; // class OptDualStdVector
256
257
258
259
260// Constraint space.
261template <class Real, class Element>
262class ConStdVector : public Vector<Real> {
263
264 typedef std::vector<Element> vector;
265 typedef typename vector::size_type uint;
266
267private:
268ROL::Ptr<std::vector<Element> > std_vec_;
269mutable ROL::Ptr<ConDualStdVector<Real> > dual_vec_;
270
271public:
272
273ConStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
274
275void plus( const Vector<Real> &x ) {
276 const ConStdVector &ex = dynamic_cast<const ConStdVector&>(x);
277 ROL::Ptr<const vector> xvalptr = ex.getVector();
278 uint dimension = std_vec_->size();
279 for (uint i=0; i<dimension; i++) {
280 (*std_vec_)[i] += (*xvalptr)[i];
281 }
282}
283
284void scale( const Real alpha ) {
285 uint dimension = std_vec_->size();
286 for (uint i=0; i<dimension; i++) {
287 (*std_vec_)[i] *= alpha;
288 }
289}
290
291Real dot( const Vector<Real> &x ) const {
292 Real val = 0;
293 const ConStdVector<Real, Element> & ex = dynamic_cast<const ConStdVector<Real, Element>&>(x);
294 ROL::Ptr<const vector> xvalptr = ex.getVector();
295
296 uint dimension = std_vec_->size();
297 for (uint i=0; i<dimension; i++) {
298 val += (*std_vec_)[i]*(*xvalptr)[i];
299 }
300 return val;
301}
302
303Real norm() const {
304 Real val = 0;
305 val = std::sqrt( dot(*this) );
306 return val;
307}
308
309ROL::Ptr<Vector<Real> > clone() const {
310 return ROL::makePtr<ConStdVector>( ROL::makePtr<vector>(std_vec_->size()));
311}
312
313ROL::Ptr<const std::vector<Element> > getVector() const {
314 return std_vec_;
315}
316
317ROL::Ptr<std::vector<Element> > getVector() {
318 return std_vec_;
319}
320
321ROL::Ptr<Vector<Real> > basis( const int i ) const {
322 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(),0.0);
323 ROL::Ptr<ConStdVector> e = ROL::makePtr<ConStdVector>( e_ptr);
324 (*e_ptr)[i] = 1.0;
325 return e;
326}
327
328int dimension() const {return static_cast<int>(std_vec_->size());}
329
330const Vector<Real> & dual() const {
331 dual_vec_ = ROL::makePtr<ConDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
332 return *dual_vec_;
333}
334
335}; // class ConStdVector
336
337
338// Dual constraint space.
339template <class Real, class Element>
340class ConDualStdVector : public Vector<Real> {
341
342 typedef std::vector<Element> vector;
343 typedef typename vector::size_type uint;
344
345private:
346
347ROL::Ptr<std::vector<Element> > std_vec_;
348mutable ROL::Ptr<ConStdVector<Real> > dual_vec_;
349
350public:
351
352ConDualStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
353
354void plus( const Vector<Real> &x ) {
355 const ConDualStdVector &ex = dynamic_cast<const ConDualStdVector&>(x);
356 ROL::Ptr<const vector> xvalptr = ex.getVector();
357 uint dimension = std_vec_->size();
358 for (uint i=0; i<dimension; i++) {
359 (*std_vec_)[i] += (*xvalptr)[i];
360 }
361}
362
363void scale( const Real alpha ) {
364 uint dimension = std_vec_->size();
365 for (uint i=0; i<dimension; i++) {
366 (*std_vec_)[i] *= alpha;
367 }
368}
369
370Real dot( const Vector<Real> &x ) const {
371 Real val = 0;
372 const ConDualStdVector<Real, Element> & ex = dynamic_cast<const ConDualStdVector<Real, Element>&>(x);
373 ROL::Ptr<const vector> xvalptr = ex.getVector();
374 uint dimension = std_vec_->size();
375 for (uint i=0; i<dimension; i++) {
376 val += (*std_vec_)[i]*(*xvalptr)[i];
377 }
378 return val;
379}
380
381Real norm() const {
382 Real val = 0;
383 val = std::sqrt( dot(*this) );
384 return val;
385}
386
387ROL::Ptr<Vector<Real> > clone() const {
388 return ROL::makePtr<ConDualStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()));
389}
390
391ROL::Ptr<const std::vector<Element> > getVector() const {
392 return std_vec_;
393}
394
395ROL::Ptr<std::vector<Element> > getVector() {
396 return std_vec_;
397}
398
399ROL::Ptr<Vector<Real> > basis( const int i ) const {
400 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(),0.0);
401 ROL::Ptr<ConDualStdVector> e = ROL::makePtr<ConDualStdVector>( e_ptr );
402 (*e_ptr)[i] = 1.0;
403 return e;
404}
405
406int dimension() const {return static_cast<int>(std_vec_->size());}
407
408const Vector<Real> & dual() const {
409 dual_vec_ = ROL::makePtr<ConStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
410 return *dual_vec_;
411}
412
413}; // class ConDualStdVector
414
415/*** End of declaration of four vector spaces. ***/
416
417
418
420template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real> >
421class Objective_GrossPitaevskii : public Objective<Real> {
422
423 typedef std::vector<Real> vector;
424 typedef typename vector::size_type uint;
425
426 private:
427
429 Real g_;
430
432 uint nx_;
433
435 Real dx_;
436
438 ROL::Ptr<const std::vector<Real> > Vp_;
439
440 ROL::Ptr<FiniteDifference<Real> > fd_;
441
443
445 void applyK(const Vector<Real> &v, Vector<Real> &kv) {
446
447
448
449 // Pointer to direction vector
450 ROL::Ptr<const vector> vp = dynamic_cast<const XPrim&>(v).getVector();
451
452 // Pointer to action of Hessian on direction vector
453 ROL::Ptr<vector> kvp = dynamic_cast<XDual&>(kv).getVector();
454
455 Real dx2 = dx_*dx_;
456
457 (*kvp)[0] = (2.0*(*vp)[0]-(*vp)[1])/dx2;
458
459 for(uint i=1;i<nx_-1;++i) {
460 (*kvp)[i] = (2.0*(*vp)[i]-(*vp)[i-1]-(*vp)[i+1])/dx2;
461 }
462
463 (*kvp)[nx_-1] = (2.0*(*vp)[nx_-1]-(*vp)[nx_-2])/dx2;
464
465 }
466
467 public:
468
469 Objective_GrossPitaevskii(const Real &g, const Vector<Real> &V, ROL::Ptr<FiniteDifference<Real> > fd) : g_(g),
470 Vp_((dynamic_cast<const StdVector<Real>&>(V)).getVector()), fd_(fd) {
471
472 nx_ = Vp_->size();
473 dx_ = (1.0/(1.0+nx_));
474 }
475
477
481 Real value( const Vector<Real> &psi, Real &tol ) {
482
483
484
485 // Pointer to opt vector
486 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
487
488 // Pointer to K applied to opt vector
489 ROL::Ptr<vector> kpsip = ROL::makePtr<vector>(nx_, 0.0);
490 XDual kpsi(kpsip,fd_);
491
492 Real J = 0;
493
494 applyK(psi,kpsi);
495
496 for(uint i=0;i<nx_;++i) {
497 J += (*psip)[i]*(*kpsip)[i] + (*Vp_)[i]*pow((*psip)[i],2) + g_*pow((*psip)[i],4);
498 }
499
500 // Rescale for numerical integration by trapezoidal rule
501 J *= 0.5*dx_;
502
503 return J;
504 }
505
507
508 void gradient( Vector<Real> &g, const Vector<Real> &psi, Real &tol ) {
509
510
511
512 // Pointer to opt vector
513 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
514
515 // Pointer to gradient vector
516 ROL::Ptr<vector> gp = dynamic_cast<XDual&>(g).getVector();
517
518 // Pointer to K applied to opt vector
519 ROL::Ptr<vector> kpsip = ROL::makePtr<vector>(nx_, 0.0);
520 XDual kpsi(kpsip,fd_);
521
522 applyK(psi,kpsi);
523
524 for(uint i=0;i<nx_;++i) {
525 (*gp)[i] = ((*kpsip)[i] + (*Vp_)[i]*(*psip)[i] + 2.0*g_*pow((*psip)[i],3))*dx_;
526 }
527
528 }
529
530
531
533
534 void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol ) {
535
536
537
538 // Pointer to opt vector
539 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
540
541 // Pointer to direction vector
542 ROL::Ptr<const vector> vp = dynamic_cast<const XPrim&>(v).getVector();
543
544 // Pointer to action of Hessian on direction vector
545 ROL::Ptr<vector> hvp = dynamic_cast<XDual&>(hv).getVector();
546
547 applyK(v,hv);
548
549 for(uint i=0;i<nx_;++i) {
550 (*hvp)[i] *= dx_;
551 (*hvp)[i] += ( (*Vp_)[i] + 6.0*g_*pow((*psip)[i],2) )*(*vp)[i]*dx_;
552 }
553
554 }
555
556};
557
558
560template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real>, class CPrim=StdVector<Real>, class CDual=StdVector<Real> >
561class Normalization_Constraint : public Constraint<Real> {
562
563 typedef std::vector<Real> vector;
564 typedef typename vector::size_type uint;
565
566 private:
567 uint nx_;
568 Real dx_;
569 ROL::Ptr<FiniteDifference<Real> > fd_;
571
572 public:
573 Normalization_Constraint(int n, Real dx, ROL::Ptr<FiniteDifference<Real> > fd, bool exactsolve) :
574 nx_(n), dx_(dx), fd_(fd), exactsolve_(exactsolve) {}
575
577
581 void value(Vector<Real> &c, const Vector<Real> &psi, Real &tol){
582
583
584
585 // Pointer to constraint vector (only one element)
586 ROL::Ptr<vector> cp = dynamic_cast<CPrim&>(c).getVector();
587
588 // Pointer to opt vector
589 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
590
591 (*cp)[0] = -1.0;
592 for(uint i=0;i<nx_;++i) {
593 (*cp)[0] += dx_*pow((*psip)[i],2);
594 }
595 }
596
598
600 void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
601
602
603
604 // Pointer to action of Jacobian of constraint on direction vector (yields scalar)
605 ROL::Ptr<vector> jvp = dynamic_cast<CPrim&>(jv).getVector();
606
607 // Pointer to direction vector
608 ROL::Ptr<const vector> vp = dynamic_cast<const XPrim&>(v).getVector();
609
610 // Pointer to opt vector
611 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
612
613 (*jvp)[0] = 0;
614 for(uint i=0;i<nx_;++i) {
615 (*jvp)[0] += 2.0*dx_*(*psip)[i]*(*vp)[i];
616 }
617 }
618
620
622 void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
623
624
625
626 // Pointer to action of adjoint of Jacobian of constraint on direction vector (yields vector)
627 ROL::Ptr<vector> ajvp = dynamic_cast<XDual&>(ajv).getVector();
628
629 // Pointer to direction vector
630 ROL::Ptr<const vector> vp = dynamic_cast<const CDual&>(v).getVector();
631
632 // Pointer to opt vector
633 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
634
635 for(uint i=0;i<nx_;++i) {
636 (*ajvp)[i] = 2.0*dx_*(*psip)[i]*(*vp)[0];
637 }
638 }
639
641
645 const Vector<Real> &psi, Real &tol){
646
647
648
649
650 // The pointer to action of constraint Hessian in u,v inner product
651 ROL::Ptr<vector> ahuvp = dynamic_cast<XDual&>(ahuv).getVector();
652
653 // Pointer to direction vector u
654 ROL::Ptr<const vector> up = dynamic_cast<const CDual&>(u).getVector();
655
656 // Pointer to direction vector v
657 ROL::Ptr<const vector> vp = dynamic_cast<const XPrim&>(v).getVector();
658
659 // Pointer to opt vector
660 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
661
662 for(uint i=0;i<nx_;++i) {
663 (*ahuvp)[i] = 2.0*dx_*(*vp)[i]*(*up)[0];
664 }
665 }
666
672 std::vector<Real> solveAugmentedSystem(Vector<Real> &v1, Vector<Real> &v2, const Vector<Real> &b1,
673 const Vector<Real> &b2, const Vector<Real> &psi, Real &tol) {
674
675
676
677 if(exactsolve_) {
678 ROL::Ptr<vector> v1p = dynamic_cast<XPrim&>(v1).getVector();
679 ROL::Ptr<vector> v2p = dynamic_cast<CDual&>(v2).getVector();
680 ROL::Ptr<const vector> b1p = dynamic_cast<const XDual&>(b1).getVector();
681 ROL::Ptr<const vector> b2p = dynamic_cast<const CPrim&>(b2).getVector();
682 ROL::Ptr<const vector> psip = dynamic_cast<const XPrim&>(psi).getVector();
683
684 ROL::Ptr<vector> jacp = ROL::makePtr<vector>(nx_, 0.0);
685 ROL::Ptr<vector> b1dp = ROL::makePtr<vector>(nx_, 0.0);
686
687 for(uint i=0;i<nx_;++i) {
688 (*jacp)[i] = (*psip)[i];
689 (*b1dp)[i] = (*b1p)[i];
690 }
691
692 // The Jacobian of the constraint is \f$c'(\psi)=2dx\psi\f$
693 XDual jac(jacp,fd_);
694 jac.scale(2.0*dx_);
695
696 // A Dual-in-name-only version of b1, so we can compute the desired inner products involving inv(K)
697 XDual b1d(b1dp,fd_);
698
699 // \f$ (c'K^{-1}*c'^\ast)^{-1} \f$
700 Real d = 1.0/jac.dot(jac);
701 Real p = jac.dot(b1d);
702
703 (*v2p)[0] = d*(p-(*b2p)[0]);
704
705 v1.set(jac.dual());
706 v1.scale(-(*v2p)[0]);
707 v1.plus(b1d.dual());
708
709 return std::vector<Real>(0);
710 }
711 else{
712 return Constraint<Real>::solveAugmentedSystem(v1,v2,b1,b2,psi,tol);
713 }
714 }
715};
716
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
ROL::Ptr< const std::vector< Element > > getVector() const
ROL::Ptr< ConStdVector< Real > > dual_vec_
Real dot(const Vector< Real > &x) const
Compute where .
ROL::Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
std::vector< Element > vector
int dimension() const
Return dimension of the vector space.
void plus(const Vector< Real > &x)
Compute , where .
Real norm() const
Returns where .
ConDualStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
void scale(const Real alpha)
Compute where .
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
ROL::Ptr< std::vector< Element > > getVector()
ROL::Ptr< std::vector< Element > > std_vec_
Real dot(const ROL::Vector< Real > &x) const
Compute where .
ROL::Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Real dot(const Vector< Real > &x) const
Compute where .
Real norm() const
Returns where .
ROL::Ptr< ConDualStdVector< Real > > dual_vec_
ROL::Ptr< std::vector< Element > > std_vec_
Real dot(const ROL::Vector< Real > &x) const
Compute where .
ROL::Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
void scale(const Real alpha)
Compute where .
std::vector< Element > vector
ConStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
ROL::Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
ROL::Ptr< std::vector< Element > > getVector()
ROL::Ptr< const std::vector< Element > > getVector() const
void plus(const Vector< Real > &x)
Compute , where .
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.
ROL::Ptr< FiniteDifference< Real > > fd_
std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &psi, Real &tol)
Normalization_Constraint(int n, Real dx, ROL::Ptr< FiniteDifference< Real > > fd, bool exactsolve)
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
void value(Vector< Real > &c, const Vector< Real > &psi, Real &tol)
Evaluate .
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
ROL::Ptr< const vector > getVector(const V &x)
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Objective_GrossPitaevskii(const Real &g, const Vector< Real > &V, ROL::Ptr< FiniteDifference< Real > > fd)
ROL::Ptr< const std::vector< Real > > Vp_
ROL::Ptr< FiniteDifference< Real > > fd_
ROL::Ptr< const vector > getVector(const V &x)
void applyK(const Vector< Real > &v, Vector< Real > &kv)
Apply finite difference operator.
void gradient(Vector< Real > &g, const Vector< Real > &psi, Real &tol)
Evaluate .
Real value(const Vector< Real > &psi, Real &tol)
Evaluate .
Real norm() const
Returns where .
int dimension() const
Return dimension of the vector space.
Real dot(const Vector< Real > &x) const
Compute where .
ROL::Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void scale(const Real alpha)
Compute where .
Real dot(const ROL::Vector< Real > &x) const
Compute where .
ROL::Ptr< std::vector< Element > > std_vec_
ROL::Ptr< const std::vector< Element > > getVector() const
ROL::Ptr< std::vector< Element > > getVector()
ROL::Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
void plus(const Vector< Real > &x)
Compute , where .
OptDualStdVector(const ROL::Ptr< std::vector< Element > > &std_vec, ROL::Ptr< FiniteDifference< Real > >fd)
std::vector< Element > vector
ROL::Ptr< OptStdVector< Real > > dual_vec_
ROL::Ptr< FiniteDifference< Real > > fd_
OptStdVector(const ROL::Ptr< std::vector< Element > > &std_vec, ROL::Ptr< FiniteDifference< Real > >fd)
void plus(const Vector< Real > &x)
Compute , where .
ROL::Ptr< vector > getVector()
ROL::Ptr< const std::vector< Element > > getVector() const
ROL::Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
int dimension() const
Return dimension of the vector space.
ROL::Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
const Vector< Real > & dual() const
Modify the dual of vector u to be .
Real dot(const Vector< Real > &x) const
Modify the dot product between primal variables to be .
ROL::Ptr< FiniteDifference< Real > > fd_
Real norm() const
Returns where .
ROL::Ptr< const vector > getVector() const
ROL::Ptr< OptDualStdVector< Real > > dual_vec_
ROL::Ptr< std::vector< Element > > std_vec_
void scale(const Real alpha)
Compute where .
std::vector< Element > vector
Defines the general constraint operator interface.
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
Provides the interface to evaluate objective functions.
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual void scale(const Real alpha)=0
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .