62 #include "EpetraExt_VectorOut.h"
70 bool nonlinear_expansion =
false;
72 bool symmetric =
false;
74 double g_mean_exp = 0.172988;
75 double g_std_dev_exp = 0.0380007;
80 MPI_Init(&argc,&argv);
97 MyPID = globalComm->
MyPID();
101 for (
int i=0; i<num_KL; i++)
106 int sz = basis->size();
108 if (nonlinear_expansion)
109 Cijk = basis->computeTripleProductTensor();
111 Cijk = basis->computeLinearTripleProductTensor();
116 std::cout <<
"Stochastic Galerkin expansion size = " << sz << std::endl;
119 int num_spatial_procs = -1;
121 num_spatial_procs = std::atoi(argv[1]);
123 parallelParams.
set(
"Number of Spatial Processors", num_spatial_procs);
139 nonlinear_expansion, symmetric));
145 sgParams->
sublist(
"SG Operator");
147 sgParams->
sublist(
"SG Preconditioner");
148 if (!nonlinear_expansion) {
149 sgParams->
set(
"Parameter Expansion Type",
"Linear");
150 sgParams->
set(
"Jacobian Expansion Type",
"Linear");
152 sgOpParams.
set(
"Operator Method",
"Matrix Free");
153 sgPrecParams.
set(
"Preconditioner Method",
"Mean-based");
156 sgPrecParams.
set(
"Symmetric Gauss-Seidel", symmetric);
157 sgPrecParams.
set(
"Mean Preconditioner Type",
"ML");
159 sgPrecParams.
sublist(
"Mean Preconditioner Parameters");
160 precParams.
set(
"default values",
"SA");
161 precParams.
set(
"ML output", 0);
162 precParams.
set(
"max levels",5);
163 precParams.
set(
"increasing or decreasing",
"increasing");
164 precParams.
set(
"aggregation: type",
"Uncoupled");
165 precParams.
set(
"smoother: type",
"ML symmetric Gauss-Seidel");
166 precParams.
set(
"smoother: sweeps",2);
167 precParams.
set(
"smoother: pre or post",
"both");
168 precParams.
set(
"coarse: max size", 200);
169 #ifdef HAVE_ML_AMESOS
170 precParams.
set(
"coarse: type",
"Amesos-KLU");
172 precParams.
set(
"coarse: type",
"Jacobi");
178 expansion, sg_parallel_data,
186 basis->evaluateBases(point, basis_vals);
189 for (
int i=0; i<num_KL; i++) {
190 sg_p_poly->
term(i,0)[i] = 0.0;
191 sg_p_poly->
term(i,1)[i] = 1.0 / basis_vals[i+1];
193 std::cout << *sg_p_poly << std::endl;
199 sg_x->PutScalar(0.0);
209 EpetraExt::ModelEvaluator::InArgs sg_inArgs = sg_model->
createInArgs();
210 EpetraExt::ModelEvaluator::OutArgs sg_outArgs = sg_model->
createOutArgs();
211 sg_inArgs.set_p(1, sg_p);
212 sg_inArgs.set_x(sg_x);
213 sg_outArgs.set_f(sg_f);
214 sg_outArgs.set_W(sg_J);
215 sg_outArgs.set_WPrec(sg_M);
218 sg_model->
evalModel(sg_inArgs, sg_outArgs);
222 sg_f->Norm2(&norm_f);
224 std::cout <<
"\nInitial residual norm = " << norm_f << std::endl;
229 aztec.SetAztecOption(AZ_solver, AZ_cg);
231 aztec.SetAztecOption(AZ_solver, AZ_gmres);
232 aztec.SetAztecOption(AZ_precond, AZ_none);
233 aztec.SetAztecOption(AZ_kspace, 20);
234 aztec.SetAztecOption(AZ_conv, AZ_r0);
235 aztec.SetAztecOption(AZ_output, 1);
236 aztec.SetUserOperator(sg_J.get());
237 aztec.SetPrecOperator(sg_M.
get());
238 aztec.SetLHS(sg_dx.get());
239 aztec.SetRHS(sg_f.
get());
242 aztec.Iterate(1000, 1e-12);
245 sg_x->Update(-1.0, *sg_dx, 1.0);
248 EpetraExt::VectorToMatrixMarketFile(
"stochastic_solution.mm", *sg_x);
251 EpetraExt::VectorToMatrixMarketFile(
"stochastic_RHS.mm", *sg_f);
257 EpetraExt::RowMatrixToMatrixMarketFile(
"stochastic_operator.mm",
267 EpetraExt::VectorToMatrixMarketFile(
"mean_gal.mm", mean);
268 EpetraExt::VectorToMatrixMarketFile(
"std_dev_gal.mm", std_dev);
271 EpetraExt::ModelEvaluator::OutArgs sg_outArgs2 = sg_model->
createOutArgs();
274 sg_f->PutScalar(0.0);
275 sg_outArgs2.set_f(sg_f);
276 sg_outArgs2.set_g(0, sg_g);
277 sg_model->
evalModel(sg_inArgs, sg_outArgs2);
280 sg_f->Norm2(&norm_f);
282 std::cout <<
"\nFinal residual norm = " << norm_f << std::endl;
291 std::cout.precision(16);
295 std::cout << std::endl;
296 std::cout <<
"Response Mean = " << std::endl << g_mean << std::endl;
297 std::cout <<
"Response Std. Dev. = " << std::endl << g_std_dev << std::endl;
301 if (norm_f < 1.0e-10 &&
302 std::abs(g_mean[0]-g_mean_exp) < g_tol &&
303 std::abs(g_std_dev[0]-g_std_dev_exp) < g_tol)
307 std::cout <<
"Example Passed!" << std::endl;
309 std::cout <<
"Example Failed!" << std::endl;
319 catch (std::exception& e) {
320 std::cout << e.what() << std::endl;
322 catch (std::string& s) {
323 std::cout << s << std::endl;
326 std::cout << s << std::endl;
329 std::cout <<
"Caught unknown exception!" << std::endl;
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.
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< const Epetra_Map > get_x_map() const
Return solution vector map.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
void computeStandardDeviation(Epetra_Vector &v) const
Compute standard deviation.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< const EpetraExt::MultiComm > getMultiComm() const
Get global comm.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void computeMean(Epetra_Vector &v) const
Compute mean.
Teuchos::RCP< EpetraExt::BlockVector > getBlockVector()
Get block vector.
RCP< ParameterList > sublist(const RCP< ParameterList > ¶mList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
virtual int MyPID() const =0
Nonlinear, stochastic Galerkin ModelEvaluator.
OutArgs createOutArgs() const
Create OutArgs.
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.
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)
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.
InArgs createInArgs() const
Create 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< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
coeff_type & term(ordinal_type dimension, ordinal_type order)
Get term for dimension dimension and order order.
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner operator.
static void zeroOutTimers()