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_problem_tpetra.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_PROBLEM_TPETRA_HPP
11 #define TWOD_DIFFUSION_PROBLEM_TPETRA_HPP
12 
13 #include "Teuchos_RCP.hpp"
14 #include "Teuchos_Array.hpp"
15 #include "Teuchos_Comm.hpp"
16 
17 #include "Tpetra_Map.hpp"
18 #include "Tpetra_Vector.hpp"
19 #include "Tpetra_Import.hpp"
20 #include "Tpetra_CrsGraph.hpp"
21 #include "Tpetra_CrsMatrix.hpp"
22 
23 #include "Stokhos.hpp"
25 //#include "Stokhos_Epetra.hpp"
26 
28 template <typename Scalar,
29  typename MeshScalar,
30  typename BasisScalar,
31  typename LocalOrdinal,
32  typename GlobalOrdinal,
33  typename Node>
35 public:
36 
37  typedef Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> Tpetra_Map;
38  typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> Tpetra_Vector;
39  typedef Tpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> Tpetra_Operator;
40  typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Tpetra_CrsMatrix;
41  typedef Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> Tpetra_CrsGraph;
42  typedef Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> Tpetra_Import;
43 
46  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
47  LocalOrdinal n, LocalOrdinal d,
48  BasisScalar s = 0.1, BasisScalar mu = 0.2,
49  bool log_normal = false,
50  bool eliminate_bcs = false);
51 
54 
57 
60 
62  Teuchos::RCP<const Tpetra_Map> get_p_map(LocalOrdinal l) const;
63 
65  Teuchos::RCP<const Tpetra_Map> get_g_map(LocalOrdinal j) const;
66 
69  get_p_names(LocalOrdinal l) const;
70 
73 
75  Teuchos::RCP<const Tpetra_Vector> get_p_init(LocalOrdinal l) const;
76 
79 
81  void computeResidual(const Tpetra_Vector& x,
82  const Tpetra_Vector& p,
83  Tpetra_Vector& f);
84 
86  void computeJacobian(const Tpetra_Vector& x,
87  const Tpetra_Vector& p,
88  Tpetra_CrsMatrix& J);
89 
91  void computeResponse(const Tpetra_Vector& x,
92  const Tpetra_Vector& p,
93  Tpetra_Vector& g);
94 
96 
97 protected:
98 
100  template <typename FuncT>
101  void computeA(const FuncT& func, const Tpetra_Vector& p,
102  Tpetra_CrsMatrix& jac);
103 
104 protected:
105 
106  // Mesh
107  MeshScalar h;
108  struct MeshPoint {
109  MeshPoint() : up(-1), down(-1), left(-1), right(-1), boundary(false) {}
110  MeshScalar x, y;
111  LocalOrdinal up, down, left, right;
112  bool boundary;
113  };
116  bool log_normal;
117  bool eliminate_bcs;
118 
121 
124 
127 
130 
133 
136 
139 
142 
145 
148 
149  // Functor representing a diffusion function given by a KL expansion
153  KL_Diffusion_Func(MeshScalar xyLeft, MeshScalar xyRight,
154  BasisScalar mean, BasisScalar std_dev,
155  MeshScalar L, LocalOrdinal num_KL);
156  Scalar operator() (MeshScalar x, MeshScalar y,
157  const Teuchos::Array<Scalar>& rv) const;
158  };
159 
160  // Functor representing a diffusion function given by a log-normal
161  // PC expansion
162  template <typename DiffusionFunc>
164  const DiffusionFunc& func;
165  LogNormal_Diffusion_Func(const DiffusionFunc& func_) : func(func_) {}
166  Scalar operator() (MeshScalar x, MeshScalar y,
167  const Teuchos::Array<Scalar>& rv) const {
168  return std::exp(func(x,y,rv));
169  }
170  };
171 
174 
175 };
176 
178 
179 #endif // TWOD_DIFFUSION_ME_HPP
Scalar operator()(MeshScalar x, MeshScalar y, const Teuchos::Array< Scalar > &rv) const
Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Tpetra_Map
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
KL_Diffusion_Func(MeshScalar xyLeft, MeshScalar xyRight, BasisScalar mean, BasisScalar std_dev, MeshScalar L, LocalOrdinal num_KL)
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
Return response function map.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra_Operator
Teuchos::RCP< const Tpetra_Map > p_map
Parameter vector map.
void computeResidual(const Epetra_Vector &x, const Epetra_Vector &p, Epetra_Vector &f)
Compute residual.
void computeResponse(const Epetra_Vector &x, const Epetra_Vector &p, Epetra_Vector &g)
Compute response.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< LogNormal_Diffusion_Func< KL_Diffusion_Func > > lnFunc
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
void computeJacobian(const Epetra_Vector &x, const Epetra_Vector &p, Epetra_Operator &J)
Compute Jacobian.
Teuchos::RCP< Tpetra_Vector > x_init
Initial guess.
Scalar operator()(MeshScalar x, MeshScalar y, const Teuchos::Array< Scalar > &rv) const
Teuchos::RCP< const Tpetra_Map > g_map
Response vector map.
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > Tpetra_Import
Teuchos::Array< MeshPoint > mesh
Teuchos::RCP< Stokhos::KL::ExponentialRandomField< MeshScalar > > rf
Teuchos::RCP< KL_Diffusion_Func > klFunc
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > Tpetra_CrsGraph
Teuchos::RCP< const Tpetra_Map > x_map
Solution vector map.
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
twoD_diffusion_problem(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)
Constructor.
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
Teuchos::Array< GlobalOrdinal > bcIndices
Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra_Vector
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< Tpetra_CrsMatrix > A
Matrix to store operator.
Teuchos::RCP< Tpetra_Import > importer
Importer to overlapped distribution.
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
Teuchos::RCP< Tpetra_Vector > b
RHS.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType Node
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra_CrsMatrix
Teuchos::RCP< Teuchos::Array< std::string > > p_names
Parameter names.
Teuchos::RCP< Tpetra_Vector > p_init
Initial parameters.
void computeA(const FuncT &func, const Tpetra_Vector &p, Tpetra_CrsMatrix &jac)
Fill coefficient matrix given function to evaluate diffusion coefficient.
Teuchos::RCP< Tpetra_CrsGraph > graph
Jacobian graph.
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)
A linear 2-D diffusion problem.