Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
linear2d_diffusion_pce_mpni.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 // ModelEvaluator implementing our problem
11 #include "twoD_diffusion_ME.hpp"
12 
13 // Epetra communicator
14 #ifdef HAVE_MPI
15 #include "Epetra_MpiComm.h"
16 #else
17 #include "Epetra_SerialComm.h"
18 #endif
19 
20 // NOX
21 #include "NOX.H"
22 #include "NOX_Epetra.H"
23 #include "NOX_Epetra_LinearSystem_Stratimikos.H"
24 #include "BelosTypes.hpp"
25 
26 // Stokhos Stochastic Galerkin
27 #include "Stokhos_Epetra.hpp"
29 
30 // Timing utilities
31 #include "Teuchos_TimeMonitor.hpp"
32 
33 // I/O utilities
34 #include "EpetraExt_VectorOut.h"
35 #include "EpetraExt_BlockUtility.h"
36 
37 int main(int argc, char *argv[]) {
38  int n = 32; // spatial discretization (per dimension)
39  int num_KL = 2; // number of KL terms
40  int p = 3; // polynomial order
41  double mu = 0.1; // mean of exponential random field
42  double s = 0.2; // std. dev. of exponential r.f.
43  bool nonlinear_expansion = false; // nonlinear expansion of diffusion coeff
44  // (e.g., log-normal)
45  bool symmetric = true; // use symmetric formulation
46  bool use_solver = true;
47  std::string solver_type = "RGMRES";
48 
49 // Initialize MPI
50 #ifdef HAVE_MPI
51  MPI_Init(&argc,&argv);
52 #endif
53 
54  int MyPID;
55 
56  try {
57 
58  Epetra_Object::SetTracebackMode(1);
59 
60  {
61  TEUCHOS_FUNC_TIME_MONITOR("Total PCE Calculation Time");
62 
63  // Create a communicator for Epetra objects
65 #ifdef HAVE_MPI
66  globalComm = Teuchos::rcp(new Epetra_MpiComm(MPI_COMM_WORLD));
67 #else
68  globalComm = Teuchos::rcp(new Epetra_SerialComm);
69 #endif
70  MyPID = globalComm->MyPID();
71 
72  // Create Stochastic Galerkin basis and quadrature
74  for (int i=0; i<num_KL; i++)
75  bases[i] = Teuchos::rcp(new Stokhos::LegendreBasis<int,double>(p, true));
78  1e-12));
79  // Teuchos::RCP<const Stokhos::Quadrature<int,double> > quad =
80  // Teuchos::rcp(new Stokhos::TensorProductQuadrature<int,double>(basis));
82  Teuchos::rcp(new Stokhos::SparseGridQuadrature<int,double>(basis));
83  const Teuchos::Array< Teuchos::Array<double> >& quad_points =
84  quad->getQuadPoints();
85  const Teuchos::Array<double>& quad_weights = quad->getQuadWeights();
86  const Teuchos::Array< Teuchos::Array<double> >& quad_values =
87  quad->getBasisAtQuadPoints();
88  int sz = basis->size();
89  int num_mp = quad_weights.size();
90  if (MyPID == 0)
91  std::cout << "Stochastic Galerkin expansion size = " << sz << std::endl;
92 
93  // Create multi-point parallel distribution
94  int num_spatial_procs = -1;
95  if (argc > 1)
96  num_spatial_procs = std::atoi(argv[1]);
98  Stokhos::buildMultiComm(*globalComm, num_mp, num_spatial_procs);
100  Stokhos::getStochasticComm(multi_comm);
102  Stokhos::getSpatialComm(multi_comm);
103  Teuchos::RCP<const Epetra_Map> mp_block_map =
104  Teuchos::rcp(new Epetra_Map(num_mp, 0, *mp_comm));
105  int num_my_mp = mp_block_map->NumMyElements();
106  //mp_block_map->Print(std::cout);
107 
108  Teuchos::RCP<Teuchos::ParameterList> detPrecParams =
110  detPrecParams->set("Preconditioner Type", "ML");
111  Teuchos::ParameterList& precParams =
112  detPrecParams->sublist("Preconditioner Parameters");
113  precParams.set("default values", "SA");
114  precParams.set("ML output", 0);
115  precParams.set("max levels",5);
116  precParams.set("increasing or decreasing","increasing");
117  precParams.set("aggregation: type", "Uncoupled");
118  precParams.set("smoother: type","ML symmetric Gauss-Seidel");
119  precParams.set("smoother: sweeps",2);
120  precParams.set("smoother: pre or post", "both");
121  precParams.set("coarse: max size", 200);
122 #ifdef HAVE_ML_AMESOS
123  precParams.set("coarse: type","Amesos-KLU");
124 #else
125  precParams.set("coarse: type","Jacobi");
126 #endif
127 
128  // Create application
130  Teuchos::rcp(new twoD_diffusion_ME(app_comm, n, num_KL, mu, s, basis,
131  nonlinear_expansion, symmetric,
132  detPrecParams));
133 
134  // Setup multi-point algorithmic parameters
137  Teuchos::ParameterList& mpPrecParams =
138  mpParams->sublist("MP Preconditioner");
139  mpPrecParams.set("Preconditioner Method", "Block Diagonal");
140  mpPrecParams.set("MP Preconditioner Type", "ML");
141  Teuchos::ParameterList& pointPrecParams =
142  mpPrecParams.sublist("MP Preconditioner Parameters");
143  pointPrecParams = precParams;
144 
145  // Create stochastic Galerkin model evaluator
147  Teuchos::rcp(new Stokhos::MPModelEvaluator(model, multi_comm,
148  mp_block_map, mpParams));
149 
150  // Set up multi-point parameters
152  mp_model->create_p_mp(0);
153  int my_mp_begin = mp_block_map->MinMyGID();
154  for (int j=0; j<num_my_mp; j++) {
155  for (int i=0; i<num_KL; i++) {
156  (*mp_p_init)[j][i] = quad_points[j+my_mp_begin][i];
157  }
158  }
159  mp_model->set_p_mp_init(0, *mp_p_init);
160 
161  // Setup multi-point initial guess
163  mp_model->create_x_mp();
164  mp_x_init->init(0.0);
165  mp_model->set_x_mp_init(*mp_x_init);
166 
167  // Set up NOX parameters
170 
171  // Set the nonlinear solver method
172  noxParams->set("Nonlinear Solver", "Line Search Based");
173 
174  // Set the printing parameters in the "Printing" sublist
175  Teuchos::ParameterList& printParams = noxParams->sublist("Printing");
176  printParams.set("MyPID", MyPID);
177  printParams.set("Output Precision", 3);
178  printParams.set("Output Processor", 0);
179  printParams.set("Output Information",
180  NOX::Utils::OuterIteration +
181  NOX::Utils::OuterIterationStatusTest +
182  NOX::Utils::InnerIteration +
183  NOX::Utils::LinearSolverDetails +
184  NOX::Utils::Warning +
185  NOX::Utils::Error);
186 
187  // Create printing utilities
188  NOX::Utils utils(printParams);
189 
190  // Sublist for line search
191  Teuchos::ParameterList& searchParams = noxParams->sublist("Line Search");
192  searchParams.set("Method", "Full Step");
193 
194  // Sublist for direction
195  Teuchos::ParameterList& dirParams = noxParams->sublist("Direction");
196  dirParams.set("Method", "Newton");
197  Teuchos::ParameterList& newtonParams = dirParams.sublist("Newton");
198  newtonParams.set("Forcing Term Method", "Constant");
199 
200  // Alternative linear solver list for Stratimikos
201  Teuchos::ParameterList stratLinSolParams;
202  Teuchos::ParameterList& stratParams =
203  stratLinSolParams.sublist("Stratimikos");
204 
205  // Sublist for linear solver for the Newton method
206  stratParams.set("Linear Solver Type", "Belos");
207  Teuchos::ParameterList& belosParams =
208  stratParams.sublist("Linear Solver Types").sublist("Belos");
209  Teuchos::ParameterList* belosSolverParams = NULL;
210  if (solver_type == "GMRES") {
211  belosParams.set("Solver Type","Block GMRES");
212  belosSolverParams =
213  &(belosParams.sublist("Solver Types").sublist("Block GMRES"));
214  }
215  else if (solver_type == "CG") {
216  belosParams.set("Solver Type","Block CG");
217  belosSolverParams =
218  &(belosParams.sublist("Solver Types").sublist("Block CG"));
219  }
220  else if (solver_type == "RGMRES") {
221  belosParams.set("Solver Type","GCRODR");
222  belosSolverParams =
223  &(belosParams.sublist("Solver Types").sublist("GCRODR"));
224  belosSolverParams->set("Num Recycled Blocks", 20);
225  }
226  else if (solver_type == "RCG") {
227  belosParams.set("Solver Type","RCG");
228  belosSolverParams =
229  &(belosParams.sublist("Solver Types").sublist("RCG"));
232  belosSolverParams->set("Num Recycled Blocks", 10);
233  }
234  else
235  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
236  "Unknown solver type " << solver_type);
237  belosSolverParams->set("Convergence Tolerance", 1e-12);
238  belosSolverParams->set("Maximum Iterations", 1000);
239  if (solver_type != "CG")
240  belosSolverParams->set("Num Blocks", 100);
241  belosSolverParams->set("Block Size", 1);
242  belosSolverParams->set("Output Frequency",1);
243  belosSolverParams->set("Output Style",1);
244  //belosSolverParams->set("Verbosity",33);
245  belosSolverParams->set("Verbosity",
246  Belos::Errors +
247  Belos::Warnings +
248  Belos::StatusTestDetails);
249  stratLinSolParams.set("Preconditioner", "User Defined");
250  Teuchos::ParameterList& verboseParams =
251  belosParams.sublist("VerboseObject");
252  verboseParams.set("Verbosity Level", "medium");
253 
254  // Sublist for convergence tests
255  Teuchos::ParameterList& statusParams = noxParams->sublist("Status Tests");
256  statusParams.set("Test Type", "Combo");
257  statusParams.set("Number of Tests", 2);
258  statusParams.set("Combo Type", "OR");
259  Teuchos::ParameterList& normF = statusParams.sublist("Test 0");
260  normF.set("Test Type", "NormF");
261  normF.set("Tolerance", 1e-10);
262  normF.set("Scale Type", "Scaled");
263  Teuchos::ParameterList& maxIters = statusParams.sublist("Test 1");
264  maxIters.set("Test Type", "MaxIters");
265  maxIters.set("Maximum Iterations", 1);
266 
267  // Create NOX interface
269  Teuchos::rcp(new NOX::Epetra::ModelEvaluatorInterface(mp_model));
270 
271  // Create NOX linear system object
273 
276  Teuchos::RCP<Epetra_Operator> M = mp_model->create_WPrec()->PrecOp;
280  if (use_solver) {
281  Teuchos::ParameterList& outerSolParams =
282  newtonParams.sublist("Linear Solver");
283  outerSolParams.sublist("Deterministic Solver Parameters") =
284  stratLinSolParams;
285  outerSolParams.set("Preconditioner Strategy", "Mean");
286  Teuchos::RCP<Epetra_Operator> inner_A = model->create_W();
287  Teuchos::RCP<Epetra_Operator> inner_M = model->create_WPrec()->PrecOp;
289  Teuchos::rcp(new NOX::Epetra::ModelEvaluatorInterface(model));
291  inner_nox_interface;
293  inner_nox_interface;
295  inner_nox_interface;
298  Teuchos::rcp(new NOX::Epetra::LinearSystemStratimikos(
299  printParams,
300  stratLinSolParams,
301  inner_iJac, inner_A, inner_iPrec, inner_M,
302  *inner_u, true));
303  linsys =
304  Teuchos::rcp(new NOX::Epetra::LinearSystemMPBD(printParams,
305  outerSolParams,
306  inner_linsys,
307  iReq, iJac, A,
308  model->get_x_map()));
309  }
310  else {
311  newtonParams.sublist("Stratimikos Linear Solver") = stratLinSolParams;
312  linsys =
313  Teuchos::rcp(new NOX::Epetra::LinearSystemStratimikos(printParams,
314  stratLinSolParams,
315  iJac, A, iPrec, M,
316  *u, true));
317  }
318 
319  // Build NOX group
321  Teuchos::rcp(new NOX::Epetra::Group(printParams, iReq, *u, linsys));
322 
323  // Create the Solver convergence test
325  NOX::StatusTest::buildStatusTests(statusParams, utils);
326 
327  // Create the solver
329  NOX::Solver::buildSolver(grp, statusTests, noxParams);
330 
331  // Solve the system
332  NOX::StatusTest::StatusType status = solver->solve();
333 
334  // Get final solution
335  const NOX::Epetra::Group& finalGroup =
336  dynamic_cast<const NOX::Epetra::Group&>(solver->getSolutionGroup());
337  const Epetra_Vector& finalSolution =
338  (dynamic_cast<const NOX::Epetra::Vector&>(finalGroup.getX())).getEpetraVector();
340  Teuchos::rcp(&finalSolution, false);
341 
342  // Evaluate final responses
343  Teuchos::RCP<const Epetra_Vector> mp_p = mp_model->get_p_init(1);
345  Teuchos::rcp(new Epetra_Vector(*(mp_model->get_g_map(0))));
346  EpetraExt::ModelEvaluator::InArgs mp_inArgs = mp_model->createInArgs();
347  EpetraExt::ModelEvaluator::OutArgs mp_outArgs = mp_model->createOutArgs();
348  mp_inArgs.set_x(mp_x);
349  mp_inArgs.set_p(1, mp_p);
350  mp_outArgs.set_g(0, mp_g);
351  mp_model->evalModel(mp_inArgs, mp_outArgs);
352 
353  // Import x and g
354  Teuchos::RCP<Epetra_LocalMap> mp_local_map =
355  Teuchos::rcp(new Epetra_LocalMap(num_mp, 0, *globalComm));
356  Teuchos::RCP<const Epetra_Map> mp_local_x_map =
357  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
358  *(model->get_x_map()), *mp_local_map, *globalComm));
359  Teuchos::RCP<const Epetra_Map> mp_local_g_map =
360  Teuchos::rcp(EpetraExt::BlockUtility::GenerateBlockMap(
361  *(model->get_g_map(0)), *mp_local_map, *globalComm));
362  Epetra_Import mp_x_importer(*mp_local_x_map, mp_x->Map());
363  Epetra_Import mp_g_importer(*mp_local_g_map, mp_g->Map());
364  Epetra_Vector mp_local_x(*mp_local_x_map);
365  Epetra_Vector mp_local_g(*mp_local_g_map);
366  mp_local_x.Import(*mp_x, mp_x_importer, Add);
367  mp_local_g.Import(*mp_g, mp_g_importer, Add);
368 
369  // Compute PC expansions
371  mp_local_map,
372  model->get_x_map(),
373  mp_local_x_map,
374  multi_comm, View, mp_local_x);
376  mp_local_map,
377  model->get_g_map(0),
378  mp_local_g_map,
379  multi_comm, View, mp_local_g);
380  Teuchos::RCP<Epetra_LocalMap> sg_block_map =
381  Teuchos::rcp(new Epetra_LocalMap(sz, 0, *globalComm));
382  Stokhos::EpetraVectorOrthogPoly sg_x(basis, sg_block_map,
383  model->get_x_map(), multi_comm);
384  Stokhos::EpetraVectorOrthogPoly sg_g(basis, sg_block_map,
385  model->get_g_map(0), multi_comm);
386  for (int i=0; i<sz; i++) {
387  sg_x[i].PutScalar(0.0);
388  sg_g[i].PutScalar(0.0);
389  for (int j=0; j<num_mp; j++) {
390  sg_x[i].Update(quad_weights[j]*quad_values[j][i], mp_x_vec[j], 1.0);
391  sg_g[i].Update(quad_weights[j]*quad_values[j][i], mp_g_vec[j], 1.0);
392  }
393  }
394 
395  // Save SG solution to file
396  EpetraExt::VectorToMatrixMarketFile("ni_stochastic_solution.mm",
397  *(sg_x.getBlockVector()));
398 
399  // Save mean and variance to file
400  Epetra_Vector mean(*(model->get_x_map()));
401  Epetra_Vector std_dev(*(model->get_x_map()));
402  sg_x.computeMean(mean);
403  sg_x.computeStandardDeviation(std_dev);
404  EpetraExt::VectorToMatrixMarketFile("mean_gal.mm", mean);
405 
406  // Print mean and standard deviation of responses
407  Epetra_Vector g_mean(*(model->get_g_map(0)));
408  Epetra_Vector g_std_dev(*(model->get_g_map(0)));
409  sg_g.computeMean(g_mean);
410  sg_g.computeStandardDeviation(g_std_dev);
411  // std::cout << "\nResponse Expansion = " << std::endl;
412  // std::cout.precision(12);
413  // sg_g.print(std::cout);
414  std::cout << std::endl;
415  std::cout << "Response Mean = " << std::endl << g_mean << std::endl;
416  std::cout << "Response Std. Dev. = " << std::endl << g_std_dev << std::endl;
417 
418  if (status == NOX::StatusTest::Converged && MyPID == 0)
419  utils.out() << "Example Passed!" << std::endl;
420 
421  }
422 
425 
426  }
427 
428  catch (std::exception& e) {
429  std::cout << e.what() << std::endl;
430  }
431  catch (std::string& s) {
432  std::cout << s << std::endl;
433  }
434  catch (char *s) {
435  std::cout << s << std::endl;
436  }
437  catch (...) {
438  std::cout << "Caught unknown exception!" <<std:: endl;
439  }
440 
441 #ifdef HAVE_MPI
442  MPI_Finalize() ;
443 #endif
444 
445 }
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::RCP< const Epetra_Comm > getStochasticComm(const Teuchos::RCP< const EpetraExt::MultiComm > &globalMultiComm)
void set_x_mp_init(const Stokhos::ProductEpetraVector &x_mp_in)
Set initial multi-point solution.
void init(const value_type &val)
Initialize coefficients.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Stokhos::ProductEpetraVector > create_x_mp(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create multi-point vector using x map and owned mp map.
Multi-point model evaluator.
RCP< ParameterList > sublist(const RCP< ParameterList > &paramList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< const EpetraExt::MultiComm > buildMultiComm(const Epetra_Comm &globalComm, int num_global_stochastic_blocks, int num_spatial_procs=-1)
virtual int MyPID() const =0
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner operator.
int NumMyElements() const
ModelEvaluator for a linear 2-D diffusion problem.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
virtual const Epetra_BlockMap & Map() const =0
A container class for products of Epetra_Vector&#39;s.
OutArgs createOutArgs() const
Create OutArgs.
int MinMyGID() const
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
Return response function map.
Legendre polynomial basis.
int main(int argc, char **argv)
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
size_type size() const
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner for W.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
InArgs createInArgs() const
Create InArgs.
Teuchos::RCP< const Epetra_Comm > getSpatialComm(const Teuchos::RCP< const EpetraExt::MultiComm > &globalMultiComm)
void set_p_mp_init(int i, const Stokhos::ProductEpetraVector &p_mp_in)
Set initial multi-point parameter.
int n
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
static void zeroOutTimers()
Teuchos::RCP< Stokhos::ProductEpetraVector > create_p_mp(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create multi-point vector using p map.