14 #include "NOX_Epetra.H"
40 noxParams->
set(
"Nonlinear Solver",
"Line Search Based");
44 printParams.
set(
"MyPID", MyPID);
45 printParams.
set(
"Output Precision", 3);
46 printParams.
set(
"Output Processor", 0);
47 printParams.
set(
"Output Information",
48 NOX::Utils::OuterIteration +
49 NOX::Utils::OuterIterationStatusTest +
50 NOX::Utils::InnerIteration +
53 NOX::Utils::LinearSolverDetails +
59 searchParams.
set(
"Method",
"Full Step");
63 dirParams.
set(
"Method",
"Newton");
65 newtonParams.
set(
"Forcing Term Method",
"Constant");
69 lsParams.
set(
"Aztec Solver",
"GMRES");
70 lsParams.
set(
"Max Iterations", 100);
71 lsParams.
set(
"Size of Krylov Subspace", 100);
72 lsParams.
set(
"Tolerance", 1e-4);
73 lsParams.
set(
"Output Frequency", 10);
74 lsParams.
set(
"Preconditioner",
"Ifpack");
78 statusParams.
set(
"Test Type",
"Combo");
79 statusParams.
set(
"Number of Tests", 2);
80 statusParams.
set(
"Combo Type",
"OR");
82 normF.
set(
"Test Type",
"NormF");
83 normF.
set(
"Tolerance", 1e-10);
84 normF.
set(
"Scale Type",
"Scaled");
86 maxIters.
set(
"Test Type",
"MaxIters");
87 maxIters.
set(
"Maximum Iterations", 10);
103 ParameterList& printParams = noxParams->sublist(
"Printing");
104 NOX::Utils utils(printParams);
107 RCP<NOX::Epetra::ModelEvaluatorInterface> nox_interface =
108 rcp(
new NOX::Epetra::ModelEvaluatorInterface(sg_model));
111 RCP<const Epetra_Vector> u = sg_model->get_x_init();
112 RCP<Epetra_Operator>
A = sg_model->create_W();
113 RCP<NOX::Epetra::Interface::Required> iReq = nox_interface;
114 RCP<NOX::Epetra::Interface::Jacobian> iJac = nox_interface;
115 RCP<NOX::Epetra::LinearSystemAztecOO> linsys;
116 RCP<Epetra_Operator> M = sg_model->create_WPrec()->PrecOp;
117 RCP<NOX::Epetra::Interface::Preconditioner> iPrec = nox_interface;
118 ParameterList& dirParams = noxParams->sublist(
"Direction");
119 ParameterList& newtonParams = dirParams.sublist(
"Newton");
120 ParameterList& lsParams = newtonParams.sublist(
"Linear Solver");
121 lsParams.set(
"Preconditioner",
"User Defined");
122 linsys =
rcp(
new NOX::Epetra::LinearSystemAztecOO(printParams, lsParams,
123 iJac, A, iPrec, M, *u));
126 RCP<NOX::Epetra::Group> grp =
127 rcp(
new NOX::Epetra::Group(printParams, iReq, *u, linsys));
130 ParameterList& statusParams = noxParams->sublist(
"Status Tests");
131 RCP<NOX::StatusTest::Generic> statusTests =
132 NOX::StatusTest::buildStatusTests(statusParams, utils);
135 RCP<NOX::Solver::Generic> solver =
136 NOX::Solver::buildSolver(grp, statusTests, noxParams);
142 const NOX::Epetra::Group& finalGroup =
143 dynamic_cast<const NOX::Epetra::Group&
>(solver.getSolutionGroup());
145 (
dynamic_cast<const NOX::Epetra::Vector&
>(finalGroup.getX())).getEpetraVector();
146 return finalSolution;
171 RCP<const Epetra_Comm> globalComm;
177 MyPID = globalComm->MyPID();
181 Array< RCP<const OneDOrthogPolyBasis<int,double> > > bases(1);
182 bases[0] =
rcp(
new LegendreBasis<int,double>(p));
183 RCP<const CompletePolynomialBasis<int,double> > basis =
184 rcp(
new CompletePolynomialBasis<int,double>(bases));
185 RCP<Stokhos::Sparse3Tensor<int,double> > Cijk =
186 basis->computeTripleProductTensor();
187 RCP<OrthogPolyExpansion<int,double> > expansion =
188 rcp(
new AlgebraicOrthogPolyExpansion<int,double>(basis, Cijk));
191 int num_spatial_procs = -1;
192 ParameterList parallelParams;
193 parallelParams.set(
"Number of Spatial Processors", num_spatial_procs);
194 RCP<ParallelData> sg_parallel_data =
195 rcp(
new ParallelData(basis, Cijk, globalComm, parallelParams));
196 RCP<const Epetra_Comm> app_comm = sg_parallel_data->getSpatialComm();
199 RCP<EpetraExt::ModelEvaluator> model =
rcp(
new SimpleME(app_comm));
202 RCP<ParameterList> sgParams =
rcp(
new ParameterList);
203 ParameterList& sgPrecParams = sgParams->sublist(
"SG Preconditioner");
204 sgPrecParams.set(
"Preconditioner Method",
"Mean-based");
207 sgPrecParams.set(
"Mean Preconditioner Type",
"Ifpack");
210 RCP<SGModelEvaluator> sg_model =
212 sg_parallel_data, sgParams));
217 RCP<EpetraVectorOrthogPoly> x_init_sg = sg_model->create_x_sg();
218 x_init_sg->init(0.0);
219 (*x_init_sg)[0] = *(model->get_x_init());
220 sg_model->set_x_sg_init(*x_init_sg);
226 RCP<EpetraVectorOrthogPoly> p_init_sg = sg_model->create_p_sg(0);
227 p_init_sg->init(0.0);
228 (*p_init_sg)[0] = *(model->get_p_init(0));
229 for (
int i=0; i<model->get_p_map(0)->NumMyElements(); i++)
230 (*p_init_sg)[i+1][i] = 1.0;
231 sg_model->set_p_sg_init(0, *p_init_sg);
232 std::cout <<
"Stochatic Galerkin parameter expansion = " << std::endl
233 << *p_init_sg << std::endl;
239 NOX::StatusTest::StatusType status = solver->solve();
245 RCP<Stokhos::EpetraVectorOrthogPoly> x_sg =
246 sg_model->create_x_sg(
View, &finalSolution);
249 std::cout <<
"Final Solution = " << std::endl;
250 std::cout << *x_sg << std::endl;
252 if (status != NOX::StatusTest::Converged)
257 if (success && MyPID == 0)
258 std::cout <<
"Example Passed!" << std::endl;
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)
Orthogonal polynomial expansions limited to algebraic operations.
Nonlinear, stochastic Galerkin ModelEvaluator.
Teuchos::RCP< Teuchos::ParameterList > create_nox_parameter_list(int MyPID)
Abstract base class for orthogonal polynomial-based expansions.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const Epetra_Vector & get_final_solution(const NOX::Solver::Generic &solver)
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
#define TEUCHOS_STANDARD_CATCH_STATEMENTS(VERBOSE, ERR_STREAM, SUCCESS_FLAG)
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
Legendre polynomial basis.
int main(int argc, char **argv)
Abstract base class for 1-D orthogonal polynomials.
Teuchos::RCP< NOX::Solver::Generic > create_nox_solver(int MyPID, const Teuchos::RCP< EpetraExt::ModelEvaluator > &sg_model)