ROL
gross-pitaevskii/example_03.cpp
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 
68 #include<fstream>
69 #include "example_03.hpp"
70 
71 typedef double RealT;
72 
73 int main(int argc, char* argv[]) {
74 
75  // Set up MPI
76  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
77 
78  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
79  int iprint = argc - 1;
80  Teuchos::RCP<std::ostream> outStream;
81  Teuchos::oblackholestream bhs; // outputs nothing
82  if (iprint > 0)
83  outStream = Teuchos::rcp(&std::cout, false);
84  else
85  outStream = Teuchos::rcp(&bhs, false);
86 
87  int errorFlag = 0;
88 
89 
90  Teuchos::ParameterList parlist;
91  std::string paramfile = "parameters.xml";
92  Teuchos::updateParametersFromXmlFile(paramfile,Teuchos::Ptr<Teuchos::ParameterList>(&parlist));
93 
94  int ni = parlist.get("Interior Grid Points",50);
95  double gnl = parlist.get("Nonlinearity Coefficient g",10.00);
96  bool useRiesz = parlist.get("Use Riesz Map",true);
97  std::string filename = parlist.get("Output File Name","results.csv");
98 
99  std::ofstream outfile(filename.c_str());
100 
101  // Number of quadrature points
102  int nq = 2*ni;
103 
104  Teuchos::RCP<Teuchos::LAPACK<int,RealT> >lapack =
105  Teuchos::rcp( new Teuchos::LAPACK<int,RealT>() );
106  Teuchos::RCP<NodalBasis<RealT> > nb = Teuchos::rcp(new NodalBasis<RealT>(lapack,ni+2,nq));
107 
108  // Exclude end point interpolants to impose Dirichlet condition
109  std::vector<RealT> L(ni*nq,0);
110  std::copy(nb->L_.begin()+nq,nb->L_.end()-nq,L.begin());
111  std::vector<RealT> Lp(ni*nq,0);
112  std::copy(nb->Lp_.begin()+nq,nb->Lp_.end()-nq,Lp.begin());
113 
114  // Quadrature points
115  std::vector<RealT> x(nb->xq_);
116 
117  // Quadrature weights
118  std::vector<RealT> w(nb->wq_);
119 
120  // Mass matrix
121  Teuchos::RCP<InnerProductMatrix<RealT> > mass =
122  Teuchos::rcp( new InnerProductMatrixSolver<RealT>(lapack,L,L,w,1) );
123 
124  // Kinetic energy matrix
125  Teuchos::RCP<InnerProductMatrix<RealT> > kinetic =
126  Teuchos::rcp( new InnerProductMatrixSolver<RealT>(lapack,Lp,Lp,w,1) );
127 
128  // Nonlinear (quartic) term. Set \f$\psi\f$ later
129  Teuchos::RCP<InnerProductMatrix<RealT> > nonlinear =
130  Teuchos::rcp( new InnerProductMatrix<RealT>(L,L,w,1) );
131 
132  // Confinement Potential
133  std::vector<RealT> v(nq,0);
134  for(int i=0;i<nq;++i){
135  v[i] = 100*x[i]*x[i];
136  }
137 
138  Teuchos::RCP<InnerProductMatrix<RealT> > potential =
139  Teuchos::rcp( new InnerProductMatrix<RealT>(L,L,w,v) );
140 
141  // Iteration Vector (pointer to optimzation vector)
142  Teuchos::RCP<std::vector<RealT> > psi_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
143 
144  // Normalized initial guess
145  for(int i=0;i<ni;++i){
146  (*psi_rcp)[i]=(1+nb->xi_[i+1])*(1-nb->xi_[i+1])*sqrt(15.0/16.0);
147  }
148  OptStdVector<RealT> psi(psi_rcp,useRiesz,kinetic);
149 
150  // Equality constraint value (scalar)
151  Teuchos::RCP<std::vector<RealT> > c_rcp = Teuchos::rcp( new std::vector<RealT> (1, 0.0) );
152  ConStdVector<RealT> c(c_rcp);
153 
154  // Lagrange multiplier value (scalar)
155  Teuchos::RCP<std::vector<RealT> > lam_rcp = Teuchos::rcp( new std::vector<RealT> (1, 0.0) );
156  ConDualStdVector<RealT> lam(lam_rcp);
157 
158 /*
159  Teuchos::RCP<std::vector<RealT> > aa_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 1.0) );
160  OptDualStdVector<RealT> av(aa_rcp,useRiesz,mass);
161  Teuchos::RCP<std::vector<RealT> > bb_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 2.0) );
162  OptDualStdVector<RealT> bv(bb_rcp,useRiesz,mass);
163  Teuchos::RCP<std::vector<RealT> > cc_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 3.0) );
164  OptDualStdVector<RealT> cv(cc_rcp,useRiesz,mass);
165  av.checkVector(bv,cv);
166 
167  Teuchos::RCP<std::vector<RealT> > dd_rcp = Teuchos::rcp( new std::vector<RealT> (1, 1.0) );
168  ConDualStdVector<RealT> dv(dd_rcp);
169  Teuchos::RCP<std::vector<RealT> > ee_rcp = Teuchos::rcp( new std::vector<RealT> (1, 2.0) );
170  ConDualStdVector<RealT> ev(ee_rcp);
171  Teuchos::RCP<std::vector<RealT> > ff_rcp = Teuchos::rcp( new std::vector<RealT> (1, 3.0) );
172  ConDualStdVector<RealT> fv(ff_rcp);
173  dv.checkVector(ev,fv);
174 */
175 
176  // Gradient
177  Teuchos::RCP<std::vector<RealT> > g_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
178  OptDualStdVector<RealT> g(g_rcp,useRiesz,kinetic);
179 
180  // Instantiate objective function
181  Teuchos::RCP<Objective<RealT> > obj =
183  (ni,gnl,nb,kinetic,potential,nonlinear));
184 
185  // Instantiate normalization constraint
186  Teuchos::RCP<EqualityConstraint<RealT> > constr =
189 
190  RealT left = -1e0, right = 1e0;
191  Teuchos::RCP<std::vector<RealT> > xtest_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
192  Teuchos::RCP<std::vector<RealT> > d_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
193  Teuchos::RCP<std::vector<RealT> > gd_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
194  Teuchos::RCP<std::vector<RealT> > vd_rcp = Teuchos::rcp( new std::vector<RealT> (ni, 0.0) );
195  Teuchos::RCP<std::vector<RealT> > vc_rcp = Teuchos::rcp( new std::vector<RealT> (1, 0.0) );
196  Teuchos::RCP<std::vector<RealT> > vl_rcp = Teuchos::rcp( new std::vector<RealT> (1, 0.0) );
197 
198  OptStdVector<RealT> xtest(xtest_rcp,useRiesz,kinetic);
199  OptStdVector<RealT> d(d_rcp,useRiesz,kinetic);
200  OptDualStdVector<RealT> gd(gd_rcp,useRiesz,kinetic);
201  OptStdVector<RealT> vd(vd_rcp,useRiesz,kinetic);
202  ConStdVector<RealT> vc(vc_rcp);
203  ConDualStdVector<RealT> vl(vl_rcp);
204 
205  for (int i=0; i<ni; i++) {
206  (*xtest_rcp)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
207  (*d_rcp)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
208  (*gd_rcp)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
209  (*vd_rcp)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
210  }
211 
212  // set vc, vl
213  (*vc_rcp)[0] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
214  (*vl_rcp)[0] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
215 
216  obj->checkGradient(xtest, g, d, true, *outStream); *outStream << "\n";
217  obj->checkHessVec(xtest, g, vd, true, *outStream); *outStream << "\n";
218  obj->checkHessSym(xtest, g, d, vd, true, *outStream); *outStream << "\n";
219  constr->checkApplyJacobian(xtest, vd, vc, true, *outStream); *outStream << "\n";
220  constr->checkApplyAdjointJacobian(xtest, vl, vc, g, true, *outStream); *outStream << "\n";
221  constr->checkApplyAdjointHessian(xtest, vl, d, g, true, *outStream); *outStream << "\n";
222 
223 
224 
225  // Define Step
226  parlist.set("Nominal SQP Optimality Solver Tolerance", 1.e-4);
227  parlist.set("Maximum Number of Krylov Iterations",80);
228  parlist.set("Absolute Krylov Tolerance",1e-4);
229 
230  CompositeStepSQP<RealT> step(parlist);
231 
232 
233  // Define Status Test
234  RealT gtol = 1e-12; // norm of gradient tolerance
235  RealT ctol = 1e-12; // norm of constraint tolerance
236  RealT stol = 1e-14; // norm of step tolerance
237  int maxit = 100; // maximum number of iterations
238  StatusTestSQP<RealT> status(gtol, ctol, stol, maxit);
239 
240  // Define Algorithm
241  DefaultAlgorithm<RealT> algo(step,status,false);
242 
243  // Run Algorithm
244  std::vector<std::string> output = algo.run(psi, g, lam, c, *obj, *constr, false);
245 
246  for ( unsigned i = 0; i < output.size(); i++ ) {
247  *outStream << output[i];
248  }
249 
250  if(algo.getState()->gnorm>1e-6) {
251  errorFlag += 1;
252  }
253 
254  // Print grid points and minimizer to file
255  outfile << "-1,0" << std::endl;
256  for(int i=0;i<ni;++i) {
257  outfile << nb->xi_[i+1] << "," << (*psi_rcp)[i] << std::endl;
258  }
259  outfile << "1,0" << std::endl;
260 
261  if (errorFlag != 0)
262  std::cout << "End Result: TEST FAILED\n";
263  else
264  std::cout << "End Result: TEST PASSED\n";
265 
266  return 0;
267 }
Teuchos::RCP< const AlgorithmState< Real > > getState(void) const
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
Implements the computation of optimization steps with composite-step trust-region SQP methods...
This class adds a solve method.
int main(int argc, char *argv[])
double RealT