22 #include "NOX_Epetra.H"
23 #include "NOX_Epetra_LinearSystem_Stratimikos.H"
24 #include "BelosTypes.hpp"
34 #include "EpetraExt_VectorOut.h"
35 #include "EpetraExt_BlockUtility.h"
43 bool nonlinear_expansion =
false;
45 bool symmetric =
true;
46 bool use_solver =
true;
47 std::string solver_type =
"RGMRES";
51 MPI_Init(&argc,&argv);
58 Epetra_Object::SetTracebackMode(1);
70 MyPID = globalComm->
MyPID();
74 for (
int i=0; i<num_KL; i++)
82 Teuchos::rcp(
new Stokhos::SparseGridQuadrature<int,double>(basis));
84 quad->getQuadPoints();
87 quad->getBasisAtQuadPoints();
88 int sz = basis->size();
89 int num_mp = quad_weights.
size();
91 std::cout <<
"Stochastic Galerkin expansion size = " << sz << std::endl;
94 int num_spatial_procs = -1;
96 num_spatial_procs = std::atoi(argv[1]);
110 detPrecParams->
set(
"Preconditioner Type",
"ML");
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");
125 precParams.
set(
"coarse: type",
"Jacobi");
131 nonlinear_expansion, symmetric,
138 mpParams->
sublist(
"MP Preconditioner");
139 mpPrecParams.
set(
"Preconditioner Method",
"Block Diagonal");
140 mpPrecParams.
set(
"MP Preconditioner Type",
"ML");
142 mpPrecParams.
sublist(
"MP Preconditioner Parameters");
143 pointPrecParams = precParams;
148 mp_block_map, mpParams));
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];
164 mp_x_init->
init(0.0);
172 noxParams->
set(
"Nonlinear Solver",
"Line Search Based");
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 +
188 NOX::Utils utils(printParams);
192 searchParams.
set(
"Method",
"Full Step");
196 dirParams.
set(
"Method",
"Newton");
198 newtonParams.
set(
"Forcing Term Method",
"Constant");
203 stratLinSolParams.
sublist(
"Stratimikos");
206 stratParams.
set(
"Linear Solver Type",
"Belos");
208 stratParams.
sublist(
"Linear Solver Types").sublist(
"Belos");
210 if (solver_type ==
"GMRES") {
211 belosParams.
set(
"Solver Type",
"Block GMRES");
213 &(belosParams.
sublist(
"Solver Types").sublist(
"Block GMRES"));
215 else if (solver_type ==
"CG") {
216 belosParams.
set(
"Solver Type",
"Block CG");
218 &(belosParams.
sublist(
"Solver Types").sublist(
"Block CG"));
220 else if (solver_type ==
"RGMRES") {
221 belosParams.
set(
"Solver Type",
"GCRODR");
223 &(belosParams.
sublist(
"Solver Types").sublist(
"GCRODR"));
224 belosSolverParams->
set(
"Num Recycled Blocks", 20);
226 else if (solver_type ==
"RCG") {
227 belosParams.
set(
"Solver Type",
"RCG");
229 &(belosParams.
sublist(
"Solver Types").sublist(
"RCG"));
232 belosSolverParams->
set(
"Num Recycled Blocks", 10);
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);
245 belosSolverParams->
set(
"Verbosity",
248 Belos::StatusTestDetails);
249 stratLinSolParams.
set(
"Preconditioner",
"User Defined");
251 belosParams.
sublist(
"VerboseObject");
252 verboseParams.
set(
"Verbosity Level",
"medium");
256 statusParams.
set(
"Test Type",
"Combo");
257 statusParams.
set(
"Number of Tests", 2);
258 statusParams.
set(
"Combo Type",
"OR");
260 normF.
set(
"Test Type",
"NormF");
261 normF.
set(
"Tolerance", 1e-10);
262 normF.
set(
"Scale Type",
"Scaled");
264 maxIters.
set(
"Test Type",
"MaxIters");
265 maxIters.
set(
"Maximum Iterations", 1);
269 Teuchos::rcp(
new NOX::Epetra::ModelEvaluatorInterface(mp_model));
282 newtonParams.
sublist(
"Linear Solver");
283 outerSolParams.
sublist(
"Deterministic Solver Parameters") =
285 outerSolParams.
set(
"Preconditioner Strategy",
"Mean");
289 Teuchos::rcp(
new NOX::Epetra::ModelEvaluatorInterface(model));
301 inner_iJac, inner_A, inner_iPrec, inner_M,
304 Teuchos::rcp(
new NOX::Epetra::LinearSystemMPBD(printParams,
311 newtonParams.
sublist(
"Stratimikos Linear Solver") = stratLinSolParams;
313 Teuchos::rcp(
new NOX::Epetra::LinearSystemStratimikos(printParams,
321 Teuchos::rcp(
new NOX::Epetra::Group(printParams, iReq, *u, linsys));
325 NOX::StatusTest::buildStatusTests(statusParams, utils);
329 NOX::Solver::buildSolver(grp, statusTests, noxParams);
332 NOX::StatusTest::StatusType status = solver->solve();
335 const NOX::Epetra::Group& finalGroup =
336 dynamic_cast<const NOX::Epetra::Group&
>(solver->getSolutionGroup());
338 (
dynamic_cast<const NOX::Epetra::Vector&
>(finalGroup.getX())).getEpetraVector();
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);
358 *(model->
get_x_map()), *mp_local_map, *globalComm));
361 *(model->
get_g_map(0)), *mp_local_map, *globalComm));
366 mp_local_x.Import(*mp_x, mp_x_importer,
Add);
367 mp_local_g.Import(*mp_g, mp_g_importer,
Add);
374 multi_comm,
View, mp_local_x);
379 multi_comm, View, mp_local_g);
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);
396 EpetraExt::VectorToMatrixMarketFile(
"ni_stochastic_solution.mm",
397 *(sg_x.getBlockVector()));
402 sg_x.computeMean(mean);
403 sg_x.computeStandardDeviation(std_dev);
404 EpetraExt::VectorToMatrixMarketFile(
"mean_gal.mm", mean);
409 sg_g.computeMean(g_mean);
410 sg_g.computeStandardDeviation(g_std_dev);
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;
418 if (status == NOX::StatusTest::Converged && MyPID == 0)
419 utils.out() <<
"Example Passed!" << std::endl;
428 catch (std::exception& e) {
429 std::cout << e.what() << std::endl;
431 catch (std::string& s) {
432 std::cout << s << std::endl;
435 std::cout << s << std::endl;
438 std::cout <<
"Caught unknown exception!" <<std:: endl;
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 > ¶mList, 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's.
OutArgs createOutArgs() const
Create OutArgs.
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.
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.
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.