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_interlaced.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 // AztecOO solver
21 #include "AztecOO.h"
22 
23 // Stokhos Stochastic Galerkin
24 #include "Stokhos_Epetra.hpp"
26 
27 // Timing utilities
28 #include "Teuchos_TimeMonitor.hpp"
29 
30 // I/O utilities
31 #include "EpetraExt_VectorOut.h"
32 #include "EpetraExt_RowMatrixOut.h"
33 
34 int main(int argc, char *argv[]) {
35  int n = 32; // spatial discretization (per dimension)
36  int num_KL = 2; // number of KL terms
37  int p = 3; // polynomial order
38  double mu = 0.2; // mean of exponential random field
39  double s = 0.1; // std. dev. of exponential r.f.
40  bool nonlinear_expansion = false; // nonlinear expansion of diffusion coeff
41  // (e.g., log-normal)
42  bool symmetric = false; // use symmetric formulation
43 
44  double g_mean_exp = 0.172988; // expected response mean
45  double g_std_dev_exp = 0.0380007; // expected response std. dev.
46  double g_tol = 1e-6; // tolerance on determining success
47 
48 // Initialize MPI
49 #ifdef HAVE_MPI
50  MPI_Init(&argc,&argv);
51 #endif
52 
53  int MyPID;
54 
55  try {
56 
57  {
58  TEUCHOS_FUNC_TIME_MONITOR("Total PCE Calculation Time");
59 
60  // Create a communicator for Epetra objects
62 #ifdef HAVE_MPI
63  globalComm = Teuchos::rcp(new Epetra_MpiComm(MPI_COMM_WORLD));
64 #else
65  globalComm = Teuchos::rcp(new Epetra_SerialComm);
66 #endif
67  MyPID = globalComm->MyPID();
68 
69  // Create Stochastic Galerkin basis and expansion
71  for (int i=0; i<num_KL; i++)
72  bases[i] = Teuchos::rcp(new Stokhos::LegendreBasis<int,double>(p,true));
75  1e-12));
76  int sz = basis->size();
78  if (nonlinear_expansion)
79  Cijk = basis->computeTripleProductTensor();
80  else
81  Cijk = basis->computeLinearTripleProductTensor();
84  Cijk));
85  if (MyPID == 0)
86  std::cout << "Stochastic Galerkin expansion size = " << sz << std::endl;
87 
88  // Create stochastic parallel distribution
89  int num_spatial_procs = -1;
90  Teuchos::ParameterList parallelParams;
91  parallelParams.set("Number of Spatial Processors", num_spatial_procs);
92  // parallelParams.set("Rebalance Stochastic Graph", true);
93  // Teuchos::ParameterList& isorropia_params =
94  // parallelParams.sublist("Isorropia");
95  // isorropia_params.set("Balance objective", "nonzeros");
96  Teuchos::RCP<Stokhos::ParallelData> sg_parallel_data =
97  Teuchos::rcp(new Stokhos::ParallelData(basis, Cijk, globalComm,
98  parallelParams));
100  sg_parallel_data->getMultiComm();
102  sg_parallel_data->getSpatialComm();
103 
104  // Create application
106  Teuchos::rcp(new twoD_diffusion_ME(app_comm, n, num_KL, s, mu, basis,
107  nonlinear_expansion, symmetric));
108 
109  // Setup stochastic Galerkin algorithmic parameters
112  if (!nonlinear_expansion) {
113  sgParams->set("Parameter Expansion Type", "Linear");
114  sgParams->set("Jacobian Expansion Type", "Linear");
115  }
116 
117  Teuchos::ParameterList precParams;
118  precParams.set("default values", "SA");
119  precParams.set("ML output", 0);
120  precParams.set("max levels",5);
121  precParams.set("increasing or decreasing","increasing");
122  precParams.set("aggregation: type", "Uncoupled");
123  precParams.set("smoother: type","ML symmetric Gauss-Seidel");
124  precParams.set("smoother: sweeps",2);
125  precParams.set("smoother: pre or post", "both");
126  precParams.set("coarse: max size", 200);
127  precParams.set("PDE equations",sz);
128 #ifdef HAVE_ML_AMESOS
129  precParams.set("coarse: type","Amesos-KLU");
130 #else
131  precParams.set("coarse: type","Jacobi");
132 #endif
133 
134  // Create stochastic Galerkin model evaluator
137  model, basis, Teuchos::null,
138  expansion, sg_parallel_data,
139  sgParams));
140 
141  // Set up stochastic parameters
142  // The current implementation of the model doesn't actually use these
143  // values, but is hard-coded to certain uncertainty models
144  Teuchos::Array<double> point(num_KL, 1.0);
145  Teuchos::Array<double> basis_vals(sz);
146  basis->evaluateBases(point, basis_vals);
148  sg_model->create_p_sg(0);
149  for (int i=0; i<num_KL; i++) {
150  sg_p_poly->term(i,0)[i] = 0.0;
151  sg_p_poly->term(i,1)[i] = 1.0 / basis_vals[i+1];
152  }
153 
154  // Create vectors and operators
157  Teuchos::rcp(new Epetra_Vector(*(sg_model->get_x_map())));
158  sg_x->PutScalar(0.0);
160  Teuchos::rcp(new Epetra_Vector(*(sg_model->get_f_map())));
162  Teuchos::rcp(new Epetra_Vector(*(sg_model->get_x_map())));
164  Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(sg_model->create_W());
166  Teuchos::rcp(new ML_Epetra::MultiLevelPreconditioner(*sg_J, precParams,
167  false));
168 
169  // Setup InArgs and OutArgs
170  EpetraExt::ModelEvaluator::InArgs sg_inArgs = sg_model->createInArgs();
171  EpetraExt::ModelEvaluator::OutArgs sg_outArgs = sg_model->createOutArgs();
172  sg_inArgs.set_p(1, sg_p);
173  sg_inArgs.set_x(sg_x);
174  sg_outArgs.set_f(sg_f);
175  sg_outArgs.set_W(sg_J);
176 
177  // Evaluate model
178  sg_model->evalModel(sg_inArgs, sg_outArgs);
179  sg_M->ComputePreconditioner();
180 
181  // Print initial residual norm
182  double norm_f;
183  sg_f->Norm2(&norm_f);
184  if (MyPID == 0)
185  std::cout << "\nInitial residual norm = " << norm_f << std::endl;
186 
187  // Setup AztecOO solver
188  AztecOO aztec;
189  if (symmetric)
190  aztec.SetAztecOption(AZ_solver, AZ_cg);
191  else
192  aztec.SetAztecOption(AZ_solver, AZ_gmres);
193  aztec.SetAztecOption(AZ_precond, AZ_none);
194  aztec.SetAztecOption(AZ_kspace, 20);
195  aztec.SetAztecOption(AZ_conv, AZ_r0);
196  aztec.SetAztecOption(AZ_output, 1);
197  aztec.SetUserOperator(sg_J.get());
198  aztec.SetPrecOperator(sg_M.get());
199  aztec.SetLHS(sg_dx.get());
200  aztec.SetRHS(sg_f.get());
201 
202  // Solve linear system
203  aztec.Iterate(1000, 1e-12);
204 
205  // Update x
206  sg_x->Update(-1.0, *sg_dx, 1.0);
207 
208  // Save solution to file
209  EpetraExt::VectorToMatrixMarketFile("stochastic_solution_interlaced.mm",
210  *sg_x);
211 
212  // Save RHS to file
213  EpetraExt::VectorToMatrixMarketFile("stochastic_RHS_interlaced.mm",
214  *sg_f);
215 
216  // Save operator to file
217  EpetraExt::RowMatrixToMatrixMarketFile("stochastic_operator_interlaced.mm",
218  *sg_J);
219 
220  // Save mean and variance to file
222  sg_model->create_x_sg(View, sg_x.get());
223  Epetra_Vector mean(*(model->get_x_map()));
224  Epetra_Vector std_dev(*(model->get_x_map()));
225  sg_x_poly->computeMean(mean);
226  sg_x_poly->computeStandardDeviation(std_dev);
227  EpetraExt::VectorToMatrixMarketFile("mean_gal_interlaced.mm", mean);
228  EpetraExt::VectorToMatrixMarketFile("std_dev_gal_interlaced.mm", std_dev);
229 
230  // Compute new residual & response function
231  EpetraExt::ModelEvaluator::OutArgs sg_outArgs2 = sg_model->createOutArgs();
233  Teuchos::rcp(new Epetra_Vector(*(sg_model->get_g_map(0))));
234  sg_f->PutScalar(0.0);
235  sg_outArgs2.set_f(sg_f);
236  sg_outArgs2.set_g(0, sg_g);
237  sg_model->evalModel(sg_inArgs, sg_outArgs2);
238 
239  // Print initial residual norm
240  sg_f->Norm2(&norm_f);
241  if (MyPID == 0)
242  std::cout << "\nFinal residual norm = " << norm_f << std::endl;
243 
244  // Print mean and standard deviation of responses
246  sg_model->create_g_sg(0, View, sg_g.get());
247  Epetra_Vector g_mean(*(model->get_g_map(0)));
248  Epetra_Vector g_std_dev(*(model->get_g_map(0)));
249  sg_g_poly->computeMean(g_mean);
250  sg_g_poly->computeStandardDeviation(g_std_dev);
251  std::cout.precision(16);
252  // std::cout << "\nResponse Expansion = " << std::endl;
253  // std::cout.precision(12);
254  // sg_g_poly->print(std::cout);
255  std::cout << std::endl;
256  std::cout << "Response Mean = " << std::endl << g_mean << std::endl;
257  std::cout << "Response Std. Dev. = " << std::endl << g_std_dev << std::endl;
258 
259  // Determine if example passed
260  bool passed = false;
261  if (norm_f < 1.0e-10 &&
262  std::abs(g_mean[0]-g_mean_exp) < g_tol &&
263  std::abs(g_std_dev[0]-g_std_dev_exp) < g_tol)
264  passed = true;
265  if (MyPID == 0) {
266  if (passed)
267  std::cout << "Example Passed!" << std::endl;
268  else
269  std::cout << "Example Failed!" << std::endl;
270  }
271 
272  }
273 
276 
277  }
278 
279  catch (std::exception& e) {
280  std::cout << e.what() << std::endl;
281  }
282  catch (std::string& s) {
283  std::cout << s << std::endl;
284  }
285  catch (char *s) {
286  std::cout << s << std::endl;
287  }
288  catch (...) {
289  std::cout << "Caught unknown exception!" << std::endl;
290  }
291 
292 #ifdef HAVE_MPI
293  MPI_Finalize() ;
294 #endif
295 
296 }
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
void computeStandardDeviation(Epetra_Vector &v) const
Compute standard deviation.
Teuchos::RCP< const EpetraExt::MultiComm > getMultiComm() const
Get global comm.
void computeMean(Epetra_Vector &v) const
Compute mean.
Teuchos::RCP< EpetraExt::BlockVector > getBlockVector()
Get block vector.
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual int MyPID() const =0
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_g_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using g map.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
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)
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
Return response function map.
Legendre polynomial basis.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
int main(int argc, char **argv)
Teuchos::RCP< const Epetra_Comm > getSpatialComm() const
Get spatial comm.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_p_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using p map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using x map and owned sg map.
coeff_type & term(ordinal_type dimension, ordinal_type order)
Get term for dimension dimension and order order.
int n
static void zeroOutTimers()
Nonlinear, stochastic Galerkin ModelEvaluator that constructs a interlaced Jacobian.