Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
twoD_diffusion_ME.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TWOD_DIFFUSION_ME_HPP
11 #define TWOD_DIFFUSION_ME_HPP
12 
13 #include "Teuchos_RCP.hpp"
14 #include "Teuchos_Array.hpp"
16 
17 #include "EpetraExt_ModelEvaluator.h"
18 
19 #include "Epetra_Map.h"
20 #include "Epetra_LocalMap.h"
21 #include "Epetra_Import.h"
22 #include "Epetra_CrsGraph.h"
23 #include "Epetra_CrsMatrix.h"
24 
25 #include "Stokhos.hpp"
26 #include "Stokhos_Epetra.hpp"
28 
30 class twoD_diffusion_ME : public EpetraExt::ModelEvaluator {
31 public:
32 
35  const Teuchos::RCP<const Epetra_Comm>& comm, int n, int d,
36  double s = 0.1, double mu = 0.2,
39  bool log_normal = false,
40  bool eliminate_bcs = false,
42 
45 
48 
51 
54 
57 
60 
63  get_p_names(int l) const;
64 
67 
70 
73 
76 
78  InArgs createInArgs() const;
79 
81  OutArgs createOutArgs() const;
82 
84  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
85 
87 
90 
91 protected:
92 
94  template <typename FuncT>
95  void fillMatrices(const FuncT& func, int sz);
96 
97 protected:
98 
99  double h;
100  struct MeshPoint {
101  MeshPoint() : up(-1), down(-1), left(-1), right(-1), boundary(false) {}
102  double x, y;
103  int up, down, left, right;
104  bool boundary;
105  };
108 
112 
115 
118 
121 
124 
127 
130 
133 
136 
139 
142 
145 
148 
151 
154 
157 
158 };
159 
160 template <typename FuncT>
161 void
163 fillMatrices(const FuncT& func, int sz)
164 {
165  int NumMyElements = x_map->NumMyElements();
166  int *MyGlobalElements = x_map->MyGlobalElements();
167  double h2 = h*h;
168  double val;
169 
170  A_k.resize(sz);
171  for (int k=0; k<sz; k++) {
173  for( int i=0 ; i<NumMyElements; ++i ) {
174 
175  // Center
176  int global_idx = MyGlobalElements[i];
177  if (mesh[global_idx].boundary) {
178  if (k == 0)
179  val = 1.0; //Enforce BC in mean matrix
180  else
181  val = 0.0;
182  A_k[k]->ReplaceGlobalValues(global_idx, 1, &val, &global_idx);
183  }
184  else {
185  double a_down =
186  -func(mesh[global_idx].x, mesh[global_idx].y-h/2.0, k)/h2;
187  double a_left =
188  -func(mesh[global_idx].x-h/2.0, mesh[global_idx].y, k)/h2;
189  double a_right =
190  -func(mesh[global_idx].x+h/2.0, mesh[global_idx].y, k)/h2;
191  double a_up =
192  -func(mesh[global_idx].x, mesh[global_idx].y+h/2.0, k)/h2;
193 
194  // Center
195  val = -(a_down + a_left + a_right + a_up);
196  A_k[k]->ReplaceGlobalValues(global_idx, 1, &val, &global_idx);
197 
198  // Down
199  if (!(eliminate_bcs && mesh[mesh[global_idx].down].boundary))
200  A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_down,
201  &mesh[global_idx].down);
202 
203  // Left
204  if (!(eliminate_bcs && mesh[mesh[global_idx].left].boundary))
205  A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_left,
206  &mesh[global_idx].left);
207 
208  // Right
209  if (!(eliminate_bcs && mesh[mesh[global_idx].right].boundary))
210  A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_right,
211  &mesh[global_idx].right);
212 
213  // Up
214  if (!(eliminate_bcs && mesh[mesh[global_idx].up].boundary))
215  A_k[k]->ReplaceGlobalValues(global_idx, 1, &a_up,
216  &mesh[global_idx].up);
217  }
218  }
219  A_k[k]->FillComplete();
220  A_k[k]->OptimizeStorage();
221  }
222 }
223 
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.
expr val()
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.
int n
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.