ROL
gross-pitaevskii/example_02.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
82 #include <iostream>
83 
84 #include "Teuchos_oblackholestream.hpp"
85 #include "Teuchos_GlobalMPISession.hpp"
86 #include "Teuchos_XMLParameterListHelpers.hpp"
87 
88 #include "ROL_StdVector.hpp"
89 #include "ROL_Objective.hpp"
91 #include "ROL_CompositeStepSQP.hpp"
92 #include "ROL_Algorithm.hpp"
93 
95 
96 
97 using namespace ROL;
98 
99 template <class Real, class Element=Real>
100 class OptStdVector; // Optimization space.
101 
102 template <class Real, class Element=Real>
103 class OptDualStdVector; // Dual optimization space.
104 
105 template <class Real, class Element=Real>
106 class ConStdVector; // Constraint space.
107 
108 template <class Real, class Element=Real>
109 class ConDualStdVector; // Dual constraint space.
110 
111 // Vector space definitions:
112 
113 
114 // Optimization space.
115 template <class Real, class Element>
116 class OptStdVector : public Vector<Real> {
117 
118 private:
119 Teuchos::RCP<std::vector<Element> > std_vec_;
120 mutable Teuchos::RCP<OptDualStdVector<Real> > dual_vec_;
121 
122 Teuchos::RCP<FiniteDifference<Real> > fd_;
123 
124 
125 public:
126 
127 OptStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec, Teuchos::RCP<FiniteDifference<Real> >fd) :
128  std_vec_(std_vec), dual_vec_(Teuchos::null), fd_(fd) {}
129 
130 void plus( const Vector<Real> &x ) {
131  OptStdVector &ex = Teuchos::dyn_cast<OptStdVector>(const_cast <Vector<Real> &>(x));
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];
136  }
137 }
138 
139 void scale( const Real alpha ) {
140  unsigned dimension = std_vec_->size();
141  for (unsigned i=0; i<dimension; i++) {
142  (*std_vec_)[i] *= alpha;
143  }
144 }
145 
146 
148 Real dot( const Vector<Real> &x ) const {
149  Real val = 0;
150  OptStdVector<Real, Element> & ex = Teuchos::dyn_cast<OptStdVector<Real, Element> >(const_cast <Vector<Real> &>(x));
151  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
152 
153  Teuchos::RCP<std::vector<Real> > kxvalptr = Teuchos::rcp( new std::vector<Real> (std_vec_->size(), 0.0) );
154 
155  this->fd_->apply(xvalptr,kxvalptr);
156 
157  unsigned dimension = std_vec_->size();
158  for (unsigned i=0; i<dimension; i++) {
159  val += (*std_vec_)[i]*(*kxvalptr)[i];
160  }
161  return val;
162 }
163 
164 Real norm() const {
165  Real val = 0;
166  val = std::sqrt( dot(*this) );
167  return val;
168 }
169 
170 Teuchos::RCP<Vector<Real> > clone() const {
171  return Teuchos::rcp( new OptStdVector( Teuchos::rcp( new std::vector<Element>(std_vec_->size()) ),fd_ ) );
172 }
173 
174 Teuchos::RCP<const std::vector<Element> > getVector() const {
175  return std_vec_;
176 }
177 
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;
181  return e;
182 }
183 
184 int dimension() const {return std_vec_->size();}
185 
186 
188 const Vector<Real> & dual() const {
189  Teuchos::RCP<std::vector<Element> > dual_vecp = Teuchos::rcp(new std::vector<Element>(*std_vec_));
190  dual_vec_ = Teuchos::rcp( new OptDualStdVector<Real>( dual_vecp, fd_ ) );
191  this->fd_->apply(dual_vecp);
192  return *dual_vec_;
193 }
194 
195 }; // class OptStdVector
196 
197 
198 // Dual optimization space.
199 template <class Real, class Element>
200 class OptDualStdVector : public Vector<Real> {
201 
202 private:
203 Teuchos::RCP<std::vector<Element> > std_vec_;
204 mutable Teuchos::RCP<OptStdVector<Real> > dual_vec_;
205 Teuchos::RCP<FiniteDifference<Real> > fd_;
206 
207 public:
208 
209 OptDualStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec, Teuchos::RCP<FiniteDifference<Real> >fd) :
210  std_vec_(std_vec), dual_vec_(Teuchos::null), fd_(fd) {}
211 
212 void plus( const Vector<Real> &x ) {
213  OptDualStdVector &ex = Teuchos::dyn_cast<OptDualStdVector>(const_cast <Vector<Real> &>(x));
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];
218  }
219 }
220 
221 void scale( const Real alpha ) {
222  unsigned dimension = std_vec_->size();
223  for (unsigned i=0; i<dimension; i++) {
224  (*std_vec_)[i] *= alpha;
225  }
226 }
227 
228 Real dot( const Vector<Real> &x ) const {
229  Real val = 0;
230  OptDualStdVector<Real, Element> & ex = Teuchos::dyn_cast<OptDualStdVector<Real, Element> >(const_cast <Vector<Real> &>(x));
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);
234 
235  unsigned dimension = std_vec_->size();
236  for (unsigned i=0; i<dimension; i++) {
237  val += (*std_vec_)[i]*(*xvalptr)[i];
238  }
239  return val;
240 }
241 
242 Real norm() const {
243  Real val = 0;
244  val = std::sqrt( dot(*this) );
245  return val;
246 }
247 
248 Teuchos::RCP<Vector<Real> > clone() const {
249  return Teuchos::rcp( new OptDualStdVector( Teuchos::rcp( new std::vector<Element>(std_vec_->size()) ), fd_ ) );
250 }
251 
252 Teuchos::RCP<const std::vector<Element> > getVector() const {
253  return std_vec_;
254 }
255 
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;
259  return e;
260 }
261 
262 int dimension() const {return std_vec_->size();}
263 
264 const Vector<Real> & dual() const {
265  Teuchos::RCP<std::vector<Element> > dual_vecp = Teuchos::rcp(new std::vector<Element>(*std_vec_));// = new std::vector<Element>(*std_vec_);
266  dual_vec_ = Teuchos::rcp( new OptStdVector<Real>( dual_vecp, fd_ ) );
267 
268  this->fd_->solve(dual_vecp);
269  return *dual_vec_;
270 }
271 
272 }; // class OptDualStdVector
273 
274 
275 
276 
277 // Constraint space.
278 template <class Real, class Element>
279 class ConStdVector : public Vector<Real> {
280 
281 private:
282 Teuchos::RCP<std::vector<Element> > std_vec_;
283 mutable Teuchos::RCP<ConDualStdVector<Real> > dual_vec_;
284 
285 public:
286 
287 ConStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
288 
289 void plus( const Vector<Real> &x ) {
290  ConStdVector &ex = Teuchos::dyn_cast<ConStdVector>(const_cast <Vector<Real> &>(x));
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];
295  }
296 }
297 
298 void scale( const Real alpha ) {
299  unsigned dimension = std_vec_->size();
300  for (unsigned i=0; i<dimension; i++) {
301  (*std_vec_)[i] *= alpha;
302  }
303 }
304 
305 Real dot( const Vector<Real> &x ) const {
306  Real val = 0;
307  ConStdVector<Real, Element> & ex = Teuchos::dyn_cast<ConStdVector<Real, Element> >(const_cast <Vector<Real> &>(x));
308  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
309 
310 
311 
312  unsigned dimension = std_vec_->size();
313  for (unsigned i=0; i<dimension; i++) {
314  val += (*std_vec_)[i]*(*xvalptr)[i];
315  }
316  return val;
317 }
318 
319 Real norm() const {
320  Real val = 0;
321  val = std::sqrt( dot(*this) );
322  return val;
323 }
324 
325 Teuchos::RCP<Vector<Real> > clone() const {
326  return Teuchos::rcp( new ConStdVector( Teuchos::rcp(new std::vector<Element>(std_vec_->size())) ) );
327 }
328 
329 Teuchos::RCP<const std::vector<Element> > getVector() const {
330  return std_vec_;
331 }
332 
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;
336  return e;
337 }
338 
339 int dimension() const {return std_vec_->size();}
340 
341 const Vector<Real> & dual() const {
342  dual_vec_ = Teuchos::rcp( new ConDualStdVector<Real>( Teuchos::rcp( new std::vector<Element>(*std_vec_) ) ) );
343  return *dual_vec_;
344 }
345 
346 }; // class ConStdVector
347 
348 
349 // Dual constraint space.
350 template <class Real, class Element>
351 class ConDualStdVector : public Vector<Real> {
352 private:
353 
354 Teuchos::RCP<std::vector<Element> > std_vec_;
355 mutable Teuchos::RCP<ConStdVector<Real> > dual_vec_;
356 
357 public:
358 
359 ConDualStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
360 
361 void plus( const Vector<Real> &x ) {
362  ConDualStdVector &ex = Teuchos::dyn_cast<ConDualStdVector>(const_cast <Vector<Real> &>(x));
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];
367  }
368 }
369 
370 void scale( const Real alpha ) {
371  unsigned dimension = std_vec_->size();
372  for (unsigned i=0; i<dimension; i++) {
373  (*std_vec_)[i] *= alpha;
374  }
375 }
376 
377 Real dot( const Vector<Real> &x ) const {
378  Real val = 0;
379  ConDualStdVector<Real, Element> & ex = Teuchos::dyn_cast<ConDualStdVector<Real, Element> >(const_cast <Vector<Real> &>(x));
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];
384  }
385  return val;
386 }
387 
388 Real norm() const {
389  Real val = 0;
390  val = std::sqrt( dot(*this) );
391  return val;
392 }
393 
394 Teuchos::RCP<Vector<Real> > clone() const {
395  return Teuchos::rcp( new ConDualStdVector( Teuchos::rcp(new std::vector<Element>(std_vec_->size())) ) );
396 }
397 
398 Teuchos::RCP<const std::vector<Element> > getVector() const {
399  return std_vec_;
400 }
401 
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;
405  return e;
406 }
407 
408 int dimension() const {return std_vec_->size();}
409 
410 const Vector<Real> & dual() const {
411  dual_vec_ = Teuchos::rcp( new ConStdVector<Real>( Teuchos::rcp( new std::vector<Element>(*std_vec_) ) ) );
412  return *dual_vec_;
413 }
414 
415 }; // class ConDualStdVector
416 
417 /*** End of declaration of four vector spaces. ***/
418 
419 
420 
422 template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real> >
423 class Objective_GrossPitaevskii : public Objective<Real> {
424 
425  private:
426 
428  Real g_;
429 
431  int nx_;
432 
434  Real dx_;
435 
437  Teuchos::RCP<const std::vector<Real> > Vp_;
438 
439  Teuchos::RCP<FiniteDifference<Real> > fd_;
440 
442 
444  void applyK(const Vector<Real> &v, Vector<Real> &kv) {
445 
446  // Pointer to direction vector
447  Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(v))).getVector();
448 
449  // Pointer to action of Hessian on direction vector
450  Teuchos::RCP<std::vector<Real> > kvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(kv)).getVector());
451 
452  Real dx2 = dx_*dx_;
453 
454  (*kvp)[0] = (2.0*(*vp)[0]-(*vp)[1])/dx2;
455 
456  for(int i=1;i<nx_-1;++i) {
457  (*kvp)[i] = (2.0*(*vp)[i]-(*vp)[i-1]-(*vp)[i+1])/dx2;
458  }
459 
460  (*kvp)[nx_-1] = (2.0*(*vp)[nx_-1]-(*vp)[nx_-2])/dx2;
461 
462  }
463 
464  public:
465 
466  Objective_GrossPitaevskii(const Real &g, const Vector<Real> &V, Teuchos::RCP<FiniteDifference<Real> > fd) : g_(g),
467  Vp_((Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(V))).getVector()), fd_(fd) {
468 
469  nx_ = Vp_->size();
470  dx_ = (1.0/(1.0+nx_));
471  }
472 
474 
478  Real value( const Vector<Real> &psi, Real &tol ) {
479 
480  // Pointer to opt vector
481  Teuchos::RCP<const std::vector<Real> > psip =
482  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
483 
484 
485 
486  // Pointer to K applied to opt vector
487  Teuchos::RCP<std::vector<Real> > kpsip = Teuchos::rcp( new std::vector<Real> (nx_, 0.0) );
488  XDual kpsi(kpsip,this->fd_);
489 
490  Real J = 0;
491 
492  this->applyK(psi,kpsi);
493 
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);
496  }
497 
498  // Rescale for numerical integration by trapezoidal rule
499  J *= 0.5*dx_;
500 
501  return J;
502  }
503 
505 
506  void gradient( Vector<Real> &g, const Vector<Real> &psi, Real &tol ) {
507 
508  // Pointer to opt vector
509  Teuchos::RCP<const std::vector<Real> > psip =
510  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
511 
512 
513  // Pointer to gradient vector
514  Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(g)).getVector());
515 
516  // Pointer to K applied to opt vector
517  Teuchos::RCP<std::vector<Real> > kpsip = Teuchos::rcp( new std::vector<Real> (nx_, 0.0) );
518  XDual kpsi(kpsip,this->fd_);
519 
520  this->applyK(psi,kpsi);
521 
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_;
524  }
525 
526  }
527 
528 
529 
531 
532  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol ) {
533 
534  // Pointer to opt vector
535  Teuchos::RCP<const std::vector<Real> > psip =
536  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
537 
538 
539  // Pointer to direction vector
540  Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(v))).getVector();
541 
542  // Pointer to action of Hessian on direction vector
543  Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(hv)).getVector());
544 
545  this->applyK(v,hv);
546 
547  for(int i=0;i<nx_;++i) {
548  (*hvp)[i] *= dx_;
549  (*hvp)[i] += ( (*Vp_)[i] + 6.0*g_*pow((*psip)[i],2) )*(*vp)[i]*dx_;
550  }
551 
552  }
553 
554 };
555 
556 
558 template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real>, class CPrim=StdVector<Real>, class CDual=StdVector<Real> >
559 class Normalization_Constraint : public EqualityConstraint<Real> {
560 
561  private:
562  int nx_;
563  Real dx_;
564  Teuchos::RCP<FiniteDifference<Real> > fd_;
565  bool exactsolve_;
566 
567  public:
568  Normalization_Constraint(int n, Real dx, Teuchos::RCP<FiniteDifference<Real> > fd, bool exactsolve) :
569  nx_(n), dx_(dx), fd_(fd), exactsolve_(exactsolve) {}
570 
572 
576  void value(Vector<Real> &c, const Vector<Real> &psi, Real &tol){
577 
578  // Pointer to constraint vector (only one element)
579  Teuchos::RCP<std::vector<Real> > cp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<CPrim>(c)).getVector());
580 
581  // Pointer to opt vector
582  Teuchos::RCP<const std::vector<Real> > psip =
583  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
584 
585 
586 
587  (*cp)[0] = -1.0;
588  for(int i=0;i<nx_;++i) {
589  (*cp)[0] += dx_*pow((*psip)[i],2);
590  }
591  }
592 
594 
596  void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
597 
598  // Pointer to action of Jacobian of constraint on direction vector (yields scalar)
599  Teuchos::RCP<std::vector<Real> > jvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<CPrim>(jv)).getVector());
600 
601  // Pointer to direction vector
602  Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(v))).getVector();
603 
604  // Pointer to opt vector
605  Teuchos::RCP<const std::vector<Real> > psip =
606  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
607 
608 
609 
610  (*jvp)[0] = 0;
611  for(int i=0;i<nx_;++i) {
612  (*jvp)[0] += 2.0*dx_*(*psip)[i]*(*vp)[i];
613  }
614  }
615 
617 
619  void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
620 
621  // Pointer to action of adjoint of Jacobian of constraint on direction vector (yields vector)
622  Teuchos::RCP<std::vector<Real> > ajvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(ajv)).getVector());
623 
624  // Pointer to direction vector
625  Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<CDual>(const_cast<Vector<Real> &>(v))).getVector();
626 
627  // Pointer to opt vector
628  Teuchos::RCP<const std::vector<Real> > psip =
629  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
630 
631 
632 
633  for(int i=0;i<nx_;++i) {
634  (*ajvp)[i] = 2.0*dx_*(*psip)[i]*(*vp)[0];
635  }
636  }
637 
639 
643  const Vector<Real> &psi, Real &tol){
644 
645  // The pointer to action of constraint Hessian in u,v inner product
646  Teuchos::RCP<std::vector<Real> > ahuvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<XDual>(ahuv)).getVector());
647 
648  // Pointer to direction vector u
649  Teuchos::RCP<const std::vector<Real> > up = (Teuchos::dyn_cast<CDual>(const_cast<Vector<Real> &>(u))).getVector();
650 
651  // Pointer to direction vector v
652  Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(v))).getVector();
653 
654  // Pointer to opt vector
655  Teuchos::RCP<const std::vector<Real> > psip =
656  (Teuchos::dyn_cast<XPrim>(const_cast<Vector<Real> &>(psi))).getVector();
657 
658 
659  for(int i=0;i<nx_;++i) {
660  (*ahuvp)[i] = 2.0*dx_*(*vp)[i]*(*up)[0];
661  }
662  }
663 
669  std::vector<Real> solveAugmentedSystem(Vector<Real> &v1, Vector<Real> &v2, const Vector<Real> &b1,
670  const Vector<Real> &b2, const Vector<Real> &psi, Real &tol) {
671  if(exactsolve_) {
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();
682 
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) );
685 
686 
687 
688  for(int i=0;i<nx_;++i) {
689  (*jacp)[i] = (*psip)[i];
690  (*b1dp)[i] = (*b1p)[i];
691  }
692 
693  // The Jacobian of the constraint is \f$c'(\psi)=2dx\psi\f$
694  XDual jac(jacp,fd_);
695  jac.scale(2.0*dx_);
696 
697  // A Dual-in-name-only version of b1, so we can compute the desired inner products involving inv(K)
698  XDual b1d(b1dp,fd_);
699 
700  // \f$ (c'K^{-1}*c'^\ast)^{-1} \f$
701  Real d = 1.0/jac.dot(jac);
702  Real p = jac.dot(b1d);
703 
704  (*v2p)[0] = d*(p-(*b2p)[0]);
705 
706  v1.set(jac.dual());
707  v1.scale(-(*v2p)[0]);
708  v1.plus(b1d.dual());
709 
710  return std::vector<Real>(0);
711  }
712  else{
713  return EqualityConstraint<Real>::solveAugmentedSystem(v1,v2,b1,b2,psi,tol);
714  }
715  }
716 };
717 
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.
Definition: ROL_Vector.hpp:72
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 .
Definition: ROL_Vector.hpp:194
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_