10 #ifndef TWOD_DIFFUSION_ME_HPP
11 #define TWOD_DIFFUSION_ME_HPP
17 #include "EpetraExt_ModelEvaluator.h"
36 double s = 0.1,
double mu = 0.2,
84 void evalModel(
const InArgs& inArgs,
const OutArgs& outArgs)
const;
94 template <
typename FuncT>
160 template <
typename FuncT>
171 for (
int k=0; k<sz; k++) {
173 for(
int i=0 ; i<NumMyElements; ++i ) {
176 int global_idx = MyGlobalElements[i];
177 if (
mesh[global_idx].boundary) {
182 A_k[k]->ReplaceGlobalValues(global_idx, 1, &val, &global_idx);
186 -func(
mesh[global_idx].x,
mesh[global_idx].y-h/2.0, k)/h2;
188 -func(
mesh[global_idx].x-h/2.0,
mesh[global_idx].y, k)/h2;
190 -func(
mesh[global_idx].x+h/2.0,
mesh[global_idx].y, k)/h2;
192 -func(
mesh[global_idx].x,
mesh[global_idx].y+h/2.0, k)/h2;
195 val = -(a_down + a_left + a_right + a_up);
196 A_k[k]->ReplaceGlobalValues(global_idx, 1, &val, &global_idx);
200 A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_down,
201 &
mesh[global_idx].down);
205 A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_left,
206 &
mesh[global_idx].left);
210 A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_right,
211 &
mesh[global_idx].right);
215 A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_up,
216 &
mesh[global_idx].up);
219 A_k[k]->FillComplete();
220 A_k[k]->OptimizeStorage();
224 #endif // TWOD_DIFFUSION_ME_HPP
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::Array< int > bcIndices
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
int MyGlobalElements(int *MyGlobalElementList) const
Teuchos::RCP< Epetra_CrsMatrix > get_mean() const
Get mean matrix.
Teuchos::Array< Teuchos::RCP< Epetra_CrsMatrix > > A_k
KL coefficients of operator.
Teuchos::RCP< Epetra_Vector > x_init
Initial guess.
Teuchos::RCP< Epetra_Map > p_map
Parameter vector map.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
int NumMyElements() const
Teuchos::RCP< Stokhos::AbstractPreconditionerFactory > precFactory
Teuchos::RCP< Epetra_CrsGraph > graph
Jacobian graph.
ModelEvaluator for a linear 2-D diffusion problem.
Teuchos::Array< MeshPoint > mesh
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< Teuchos::ParameterList > precParams
Teuchos::Array< double > basis_vals
Array to store values of basis at a point.
void resize(size_type new_size, const value_type &x=value_type())
~twoD_diffusion_ME()
Destructor.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
Return response function map.
Teuchos::RCP< Epetra_Map > x_map
Solution vector map.
Teuchos::RCP< Epetra_CrsMatrix > A
Matrix to store deterministic operator.
twoD_diffusion_ME(const Teuchos::RCP< const Epetra_Comm > &comm, int n, int d, double s=0.1, double mu=0.2, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &basis=Teuchos::null, bool log_normal=false, bool eliminate_bcs=false, const Teuchos::RCP< Teuchos::ParameterList > &precParams=Teuchos::null)
Constructor.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
InArgs createInArgs() const
Create InArgs.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::RCP< Epetra_Vector > p_init
Initial parameters.
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.
Teuchos::RCP< Epetra_Vector > b
Deterministic RHS.
Teuchos::RCP< Teuchos::Array< std::string > > p_names
Parameter names.
Teuchos::Array< double > point
Array to store a point for basis evaluation.
void fillMatrices(const FuncT &func, int sz)
Fill coefficient matrix given function to evaluate diffusion coefficient.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< Epetra_Map > g_map
Response vector map.
Teuchos::RCP< Epetra_Import > importer
Importer to overlapped distribution.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > basis
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > sg_kx_vec_all
Vectors to store matrix-vector products in SG residual calculation.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.