ROL
example_01.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 
70 #include <iostream>
71 
72 #include "Teuchos_oblackholestream.hpp"
73 #include "Teuchos_GlobalMPISession.hpp"
74 #include "Teuchos_XMLParameterListHelpers.hpp"
75 
76 #include "ROL_StdVector.hpp"
77 #include "ROL_Objective.hpp"
79 #include "ROL_CompositeStepSQP.hpp"
80 #include "ROL_Algorithm.hpp"
81 
82 
83 using namespace ROL;
84 
85 template<class Real>
86 class Objective_GrossPitaevskii : public Objective<Real> {
87 
88  // STL Vector
89  typedef std::vector<Real> svec;
90 
91  // ROL StdVector
93 
94  // Pointer to const STL vector
95  typedef Teuchos::RCP<const svec> pcsv;
96 
97  // Pointer to STL vector
98  typedef Teuchos::RCP<svec> psv;
99 
100 
101  private:
102 
104  Real g_;
105 
107  int nx_;
108 
110  Real dx_;
111 
114 
116 
118  void applyK(const Vector<Real> &v, Vector<Real> &kv) {
119 
120  // Pointer to direction vector
121  pcsv vp = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(v))).getVector();
122 
123  // Pointer to action of Hessian on direction vector
124  psv kvp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(kv)).getVector());
125 
126  Real dx2 = dx_*dx_;
127 
128  (*kvp)[0] = (2.0*(*vp)[0]-(*vp)[1])/dx2;
129 
130  for(int i=1;i<nx_-1;++i) {
131  (*kvp)[i] = (2.0*(*vp)[i]-(*vp)[i-1]-(*vp)[i+1])/dx2;
132  }
133 
134  (*kvp)[nx_-1] = (2.0*(*vp)[nx_-1]-(*vp)[nx_-2])/dx2;
135 
136  }
137 
138  public:
139 
140  Objective_GrossPitaevskii(const Real &g, const Vector<Real> &V) : g_(g),
141  Vp_((Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(V))).getVector()) {
142 
143  nx_ = Vp_->size();
144  dx_ = (1.0/(1.0+nx_));
145  }
146 
148 
152  Real value( const Vector<Real> &psi, Real &tol ) {
153 
154  // Pointer to opt vector
155  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
156 
157  // Pointer to K applied to opt vector
158  psv kpsip = Teuchos::rcp( new svec(nx_, 0.0) );
159  rvec kpsi(kpsip);
160 
161  Real J = 0;
162 
163  this->applyK(psi,kpsi);
164 
165  for(int i=0;i<nx_;++i) {
166  J += (*psip)[i]*(*kpsip)[i] + (*Vp_)[i]*pow((*psip)[i],2) + g_*pow((*psip)[i],4);
167  }
168 
169  // Rescale for numerical integration by trapezoidal rule
170  J *= 0.5*dx_;
171 
172  return J;
173  }
174 
176 
177  void gradient( Vector<Real> &g, const Vector<Real> &psi, Real &tol ) {
178 
179  // Pointer to opt vector
180  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
181 
182  // Pointer to gradient vector
183  psv gp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(g)).getVector());
184 
185  // Pointer to K applied to opt vector
186  psv kpsip = Teuchos::rcp( new svec(nx_, 0.0) );
187  rvec kpsi(kpsip);
188 
189  this->applyK(psi,kpsi);
190 
191  for(int i=0;i<nx_;++i) {
192  (*gp)[i] = ((*kpsip)[i] + (*Vp_)[i]*(*psip)[i] + 2.0*g_*pow((*psip)[i],3))*dx_;
193  }
194 
195  }
196 
197 
198 
200 
201  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol ) {
202 
203  // Pointer to opt vector
204  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
205 
206  // Pointer to direction vector
207  pcsv vp = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(v))).getVector();
208 
209  // Pointer to action of Hessian on direction vector
210  psv hvp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(hv)).getVector());
211 
212  this->applyK(v,hv);
213 
214  for(int i=0;i<nx_;++i) {
215  (*hvp)[i] *= dx_;
216  (*hvp)[i] += ( (*Vp_)[i] + 6.0*g_*pow((*psip)[i],2) )*(*vp)[i]*dx_;
217  }
218 
219  }
220 
221 };
222 
223 
224 template<class Real>
226 
227  // STL Vector
228  typedef std::vector<Real> svec;
229 
230  // ROL StdVector
232 
233  // Pointer to const STL vector
234  typedef Teuchos::RCP<const svec> pcsv;
235 
236  // Pointer to STL vector
237  typedef Teuchos::RCP<svec> psv;
238 
239  private:
240  int nx_;
241  Real dx_;
242 
243  public:
244  Normalization_Constraint(int n, Real dx) : nx_(n), dx_(dx) {}
245 
247 
251  void value(Vector<Real> &c, const Vector<Real> &psi, Real &tol){
252 
253  // Pointer to constraint vector (only one element)
254  psv cp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(c)).getVector());
255 
256  // Pointer to optimization vector
257  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
258 
259  (*cp)[0] = -1.0;
260  for(int i=0;i<nx_;++i) {
261  (*cp)[0] += dx_*pow((*psip)[i],2);
262  }
263  }
264 
266 
268  void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
269 
270  // Pointer to action of Jacobian of constraint on direction vector (yields scalar)
271  psv jvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(jv)).getVector());
272 
273  // Pointer to direction vector
274  pcsv vp = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(v))).getVector();
275 
276  // Pointer to optimization vector
277  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
278 
279  (*jvp)[0] = 0;
280  for(int i=0;i<nx_;++i) {
281  (*jvp)[0] += 2.0*dx_*(*psip)[i]*(*vp)[i];
282  }
283  }
284 
286 
288  void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &psi, Real &tol){
289 
290  // Pointer to action of adjoint of Jacobian of constraint on direction vector (yields vector)
291  psv ajvp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(ajv)).getVector());
292 
293  // Pointer to direction vector
294  pcsv vp = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(v))).getVector();
295 
296  // Pointer to optimization vector
297  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
298 
299  for(int i=0;i<nx_;++i) {
300  (*ajvp)[i] = 2.0*dx_*(*psip)[i]*(*vp)[0];
301  }
302  }
303 
305 
309  const Vector<Real> &psi, Real &tol){
310 
311  // The pointer to action of constraint Hessian in u,v inner product
312  psv ahuvp = Teuchos::rcp_const_cast<svec>((Teuchos::dyn_cast<rvec>(ahuv)).getVector());
313 
314  // Pointer to direction vector u
315  pcsv up = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(u))).getVector();
316 
317  // Pointer to direction vector v
318  pcsv vp = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(v))).getVector();
319 
320  // Pointer to optimization vector
321  pcsv psip = (Teuchos::dyn_cast<rvec>(const_cast<Vector<Real> &>(psi))).getVector();
322 
323 
324  for(int i=0;i<nx_;++i) {
325  (*ahuvp)[i] = 2.0*dx_*(*vp)[i]*(*up)[0];
326  }
327  }
328 };
329 
330 
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:152
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:288
void value(Vector< Real > &c, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:251
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:201
StdVector< Real > rvec
Definition: example_01.hpp:231
Teuchos::RCP< const svec > pcsv
Definition: example_01.hpp:234
std::vector< Real > svec
Definition: example_01.hpp:89
Normalization_Constraint(int n, Real dx)
Definition: example_01.hpp:244
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
Teuchos::RCP< svec > psv
Definition: example_01.hpp:98
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:268
Defines the equality constraint operator interface.
Provides the std::vector implementation of the ROL::Vector interface.
std::vector< Real > svec
Definition: example_01.hpp:228
void applyK(const Vector< Real > &v, Vector< Real > &kv)
Apply finite difference operator.
Definition: example_01.hpp:118
void gradient(Vector< Real > &g, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:177
Teuchos::RCP< svec > psv
Definition: example_01.hpp:237
StdVector< Real > rvec
Definition: example_01.hpp:92
Objective_GrossPitaevskii(const Real &g, const Vector< Real > &V)
Definition: example_01.hpp:140
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &psi, Real &tol)
Evaluate .
Definition: example_01.hpp:308
Teuchos::RCP< const svec > pcsv
Definition: example_01.hpp:95