45 #include "EpetraExt_MatrixMatrix.h"
52 struct KL_Diffusion_Func {
57 KL_Diffusion_Func(
double xyLeft,
double xyRight,
58 double mean_,
double std_dev,
59 double L,
int num_KL) : mean(mean_), point(2)
62 rfParams.
set(
"Number of KL Terms", num_KL);
63 rfParams.
set(
"Mean", mean);
64 rfParams.
set(
"Standard Deviation", std_dev);
67 correlation_length(ndim);
68 for (
int i=0; i<ndim; i++) {
69 domain_upper[i] = xyRight;
70 domain_lower[i] = xyLeft;
71 correlation_length[i] = L;
73 rfParams.
set(
"Domain Upper Bounds", domain_upper);
74 rfParams.
set(
"Domain Lower Bounds", domain_lower);
75 rfParams.
set(
"Correlation Lengths", correlation_length);
81 ~KL_Diffusion_Func() {
84 double operator() (
double x,
double y,
int k)
const {
89 return rf->evaluate_eigenfunction(point, k-1);
95 template <
typename DiffusionFunc>
96 struct Normalized_KL_Diffusion_Func {
97 const DiffusionFunc& func;
101 Normalized_KL_Diffusion_Func(
102 const DiffusionFunc& func_,
105 d(basis.dimension()),
110 for(
int k=0; k<d; k++) {
118 double operator() (
double x,
double y,
int k)
const {
120 double val = func(
x,
y, 0);
121 for (
int i=1; i<=d; i++)
122 val -= psi_0[i]/(psi_1[i]-psi_0[i])*func(
x,
y, i);
127 return 1.0/(psi_1[k]-psi_0[k])*func(
x,
y, k);
132 template <
typename DiffusionFunc>
133 struct LogNormal_Diffusion_Func {
135 const DiffusionFunc& func;
138 LogNormal_Diffusion_Func(
140 const DiffusionFunc& func_,
142 : mean(mean_), func(func_), prodbasis(prodbasis_) {}
144 double operator() (
double x,
double y,
int k)
const {
145 int d = prodbasis->dimension();
148 double sum_g = 0.0, efval;
149 for (
int l=0; l<d; l++) {
152 efval =
std::exp(mean + 0.5*sum_g)/norms[k];
153 for (
int l=0; l<d; l++) {
154 efval *=
std::pow(func(
x,
y,l+1),multiIndex[l]);
172 log_normal(log_normal_),
173 eliminate_bcs(eliminate_bcs_),
174 precParams(precParams_)
176 Kokkos::initialize();
189 h = (xyRight - xyLeft)/((
double)(
n-1));
191 for (
int j=0;
j<
n;
j++) {
192 double y = xyLeft +
j*
h;
193 for (
int i=0; i<
n; i++) {
194 double x = xyLeft + i*
h;
198 if (i == 0 || i ==
n-1 ||
j == 0 ||
j ==
n-1)
199 mesh[idx].boundary =
true;
201 mesh[idx].left = idx-1;
203 mesh[idx].right = idx+1;
205 mesh[idx].down = idx-
n;
207 mesh[idx].up = idx+
n;
214 int n_global_dof = global_dof_indices.
size();
216 int proc_id = comm->
MyPID();
217 int n_my_dof = n_global_dof / n_proc;
218 if (proc_id == n_proc-1)
219 n_my_dof += n_global_dof % n_proc;
220 int *my_dof = global_dof_indices.
getRawPtr() + proc_id*(n_global_dof / n_proc);
240 for (
int i=0;i<d;i++) {
241 std::stringstream ss;
242 ss <<
"KL Random Variable " << i+1;
243 (*p_names)[i] = ss.str();
250 for (
int i=0; i<NumMyElements; ++i ) {
253 int global_idx = MyGlobalElements[i];
256 if (!
mesh[global_idx].boundary) {
277 KL_Diffusion_Func klFunc(xyLeft, xyRight, mu, s, 1.0, d);
280 if (
basis == Teuchos::null) {
284 Normalized_KL_Diffusion_Func<KL_Diffusion_Func> nklFunc(klFunc, *
basis);
294 LogNormal_Diffusion_Func<KL_Diffusion_Func> lnFunc(mu, klFunc, prodbasis);
303 for(
int i=0 ; i<NumMyElements; ++i ) {
304 int global_idx = MyGlobalElements[i];
305 if (
mesh[global_idx].boundary)
311 if (
basis != Teuchos::null) {
317 std::string name =
precParams->
get(
"Preconditioner Type",
"Ifpack");
354 "Error! twoD_diffusion_ME::get_p_map(): " <<
355 "Invalid parameter index l = " << l << std::endl);
367 "Error! twoD_diffusion_ME::get_g_map(): " <<
368 "Invalid parameter index l = " << l << std::endl);
380 "Error! twoD_diffusion_ME::get_p_names(): " <<
381 "Invalid parameter index l = " << l << std::endl);
400 "Error! twoD_diffusion_ME::get_p_init(): " <<
401 "Invalid parameter index l = " << l << std::endl);
424 return Teuchos::rcp(
new EpetraExt::ModelEvaluator::Preconditioner(precOp,
427 return Teuchos::null;
430 EpetraExt::ModelEvaluator::InArgs
435 inArgs.setModelEvalDescription(
"TwoD Diffusion Model Evaluator");
438 inArgs.setSupports(IN_ARG_x,
true);
442 inArgs.setSupports(IN_ARG_x_sg,
true);
443 inArgs.setSupports(IN_ARG_p_sg, 0,
true);
444 inArgs.setSupports(IN_ARG_sg_basis,
true);
445 inArgs.setSupports(IN_ARG_sg_quadrature,
true);
446 inArgs.setSupports(IN_ARG_sg_expansion,
true);
449 inArgs.setSupports(IN_ARG_x_mp,
true);
450 inArgs.setSupports(IN_ARG_p_mp, 0,
true);
455 EpetraExt::ModelEvaluator::OutArgs
459 OutArgsSetup outArgs;
460 outArgs.setModelEvalDescription(
"TwoD Diffusion Model Evaluator");
463 outArgs.set_Np_Ng(1, 1);
464 outArgs.setSupports(OUT_ARG_f,
true);
465 outArgs.setSupports(OUT_ARG_W,
true);
467 outArgs.setSupports(OUT_ARG_WPrec,
true);
470 outArgs.setSupports(OUT_ARG_f_sg,
true);
471 outArgs.setSupports(OUT_ARG_W_sg,
true);
472 outArgs.setSupports(OUT_ARG_g_sg, 0,
true);
475 outArgs.setSupports(OUT_ARG_f_mp,
true);
476 outArgs.setSupports(OUT_ARG_W_mp,
true);
477 outArgs.setSupports(OUT_ARG_g_mp, 0,
true);
484 evalModel(
const InArgs& inArgs,
const OutArgs& outArgs)
const
496 if (p == Teuchos::null)
502 if (
f != Teuchos::null || W != Teuchos::null || WPrec != Teuchos::null) {
503 if (
basis != Teuchos::null) {
509 EpetraExt::MatrixMatrix::Add((*
A_k[k]),
false,
basis_vals[k], *
A, 1.0);
514 EpetraExt::MatrixMatrix::Add((*
A_k[k]),
false, (*p)[k-1], *
A, 1.0);
517 A->OptimizeStorage();
521 if (
f != Teuchos::null) {
523 A->Apply(*det_x,*kx);
524 f->Update(1.0,*kx,-1.0, *
b, 0.0);
528 if (W != Teuchos::null) {
530 Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(W,
true);
537 if (WPrec != Teuchos::null)
542 if (
g != Teuchos::null) {
552 InArgs::sg_const_vector_t x_sg = inArgs.get_x_sg();
555 InArgs::sg_const_vector_t p_sg = inArgs.get_p_sg(0);
558 OutArgs::sg_vector_t f_sg = outArgs.get_f_sg();
559 if (f_sg != Teuchos::null) {
563 inArgs.get_sg_expansion();
576 Cijk_type::k_iterator k_begin = Cijk->k_begin();
577 Cijk_type::k_iterator k_end = Cijk->k_end();
578 for (Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
579 int k = Stokhos::index(k_it);
580 for (Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
581 j_it != Cijk->j_end(k_it); ++j_it) {
582 int j = Stokhos::index(j_it);
585 for (Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
586 j_it != Cijk->j_end(k_it); ++j_it) {
587 int j = Stokhos::index(j_it);
588 for (Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
589 i_it != Cijk->i_end(j_it); ++i_it) {
590 int i = Stokhos::index(i_it);
591 double c = Stokhos::value(i_it);
596 (*f_sg)[0].Update(-1.0,*
b,1.0);
600 OutArgs::sg_operator_t W_sg = outArgs.get_W_sg();
601 if (W_sg != Teuchos::null) {
602 for (
int i=0; i<W_sg->size(); i++) {
604 Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(W_sg->getCoeffPtr(i),
true);
607 jac->OptimizeStorage();
614 if (g_sg != Teuchos::null) {
615 int sz = x_sg->
size();
616 for (
int i=0; i<sz; i++) {
617 (*x_sg)[i].MeanValue(&(*g_sg)[i][0]);
627 mp_const_vector_t x_mp = inArgs.get_x_mp();
630 mp_const_vector_t p_mp = inArgs.get_p_mp(0);
633 mp_vector_t f_mp = outArgs.get_f_mp();
634 mp_operator_t W_mp = outArgs.get_W_mp();
635 if (f_mp != Teuchos::null || W_mp != Teuchos::null) {
636 int num_mp = x_mp->size();
637 for (
int i=0; i<num_mp; i++) {
639 if (
basis != Teuchos::null) {
641 point[k] = (*p_mp)[i][k];
651 EpetraExt::MatrixMatrix::Add((*
A_k[k]),
false, (*p_mp)[i][k-1], *
A,
656 A->OptimizeStorage();
659 if (f_mp != Teuchos::null) {
660 A->Apply((*x_mp)[i], (*f_mp)[i]);
661 (*f_mp)[i].Update(-1.0, *
b, 1.0);
665 if (W_mp != Teuchos::null) {
667 Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(W_mp->getCoeffPtr(i),
671 jac->OptimizeStorage();
677 mp_vector_t g_mp = outArgs.get_g_mp(0);
678 if (g_mp != Teuchos::null) {
679 int sz = x_mp->size();
680 for (
int i=0; i<sz; i++) {
681 (*x_mp)[i].MeanValue(&(*g_mp)[i][0]);