Example solution of a Poisson equation on a hexahedral mesh using nodal (Hgrad) elements. The system is assembled but not solved.
More...
#include "TrilinosCouplings_config.h"
#include "TrilinosCouplings_Pamgen_Utils.hpp"
#include "Intrepid_FunctionSpaceTools.hpp"
#include "Intrepid_CellTools.hpp"
#include "Intrepid_ArrayTools.hpp"
#include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
#include "Intrepid_RealSpaceTools.hpp"
#include "Intrepid_DefaultCubatureFactory.hpp"
#include "Intrepid_Utils.hpp"
#include "Tpetra_Map.hpp"
#include "Tpetra_Core.hpp"
#include "Tpetra_Export.hpp"
#include "Tpetra_Import.hpp"
#include "Tpetra_CrsGraph.hpp"
#include "Tpetra_CrsMatrix.hpp"
#include "Tpetra_DistObject.hpp"
#include "Tpetra_Vector.hpp"
#include "Tpetra_MultiVector.hpp"
#include "MatrixMarket_Tpetra.hpp"
#include "Tpetra_Operator.hpp"
#include "Teuchos_oblackholestream.hpp"
#include "Teuchos_RCP.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_XMLParameterListHelpers.hpp"
#include "Teuchos_ArrayView.hpp"
#include "Teuchos_ArrayRCP.hpp"
#include "Teuchos_TimeMonitor.hpp"
#include "Shards_CellTopology.hpp"
#include "create_inline_mesh.h"
#include "pamgen_im_exodusII_l.h"
#include "pamgen_im_ne_nemesisI_l.h"
#include "pamgen_extras.h"
#include "RTC_FunctionRTC.hh"
#include "BelosLinearProblem.hpp"
#include "BelosPseudoBlockCGSolMgr.hpp"
#include "BelosConfigDefs.hpp"
#include "BelosTpetraAdapter.hpp"
#include "Sacado_No_Kokkos.hpp"
|
typedef Sacado::Fad::SFad
< double, 3 > | Fad3 |
|
typedef
Intrepid::FunctionSpaceTools | IntrepidFSTools |
|
typedef
Intrepid::RealSpaceTools
< double > | IntrepidRSTools |
|
typedef Intrepid::CellTools
< double > | IntrepidCTools |
|
typedef Tpetra::Map::node_type | Node |
|
typedef double | ST |
|
typedef int | Ordinal |
|
typedef
Tpetra::Map::global_ordinal_type | GlobalOrdinal |
|
typedef Tpetra::Map | Map |
|
typedef Tpetra::Export | export_type |
|
typedef Tpetra::Import | import_type |
|
typedef Teuchos::ArrayView
< Ordinal >::size_type | size_type |
|
typedef Tpetra::CrsMatrix< ST > | sparse_matrix_type |
|
typedef Tpetra::CrsGraph | sparse_graph_type |
|
typedef Tpetra::MultiVector< ST > | MV |
|
typedef Tpetra::Vector< ST > | vector_type |
|
typedef Tpetra::Operator< ST > | OP |
|
typedef Belos::MultiVecTraits
< ST, MV > | MVT |
|
typedef Belos::OperatorTraits
< ST, MV, OP > | OPT |
|
|
template<typename Scalar > |
const Scalar | exactSolution (const Scalar &x, const Scalar &y, const Scalar &z) |
| User-defined exact solution. More...
|
|
template<typename Scalar > |
void | materialTensor (Scalar material[][3], const Scalar &x, const Scalar &y, const Scalar &z) |
| User-defined material tensor. More...
|
|
template<typename Scalar > |
void | exactSolutionGrad (Scalar gradExact[3], const Scalar &x, const Scalar &y, const Scalar &z) |
| Computes gradient of the exact solution. Requires user-defined exact solution. More...
|
|
template<typename Scalar > |
const Scalar | sourceTerm (Scalar &x, Scalar &y, Scalar &z) |
| Computes source term: f = -div(A.grad u). Requires user-defined exact solution and material tensor. More...
|
|
template<class ArrayOut , class ArrayIn > |
void | evaluateMaterialTensor (ArrayOut &worksetMaterialValues, const ArrayIn &evaluationPoints) |
| Computation of the material tensor at array of points in physical space. More...
|
|
template<class ArrayOut , class ArrayIn > |
void | evaluateSourceTerm (ArrayOut &sourceTermValues, const ArrayIn &evaluationPoints) |
| Computation of the source term at array of points in physical space. More...
|
|
template<class ArrayOut , class ArrayIn > |
void | evaluateExactSolution (ArrayOut &exactSolutionValues, const ArrayIn &evaluationPoints) |
| Computation of the exact solution at array of points in physical space. More...
|
|
template<class ArrayOut , class ArrayIn > |
void | evaluateExactSolutionGrad (ArrayOut &exactSolutionGradValues, const ArrayIn &evaluationPoints) |
| Computation of the gradient of the exact solution at array of points in physical space. More...
|
|
int | main (int argc, char *argv[]) |
|
Example solution of a Poisson equation on a hexahedral mesh using nodal (Hgrad) elements. The system is assembled but not solved.
This example uses the following Trilinos packages:
- Pamgen to generate a Hexahedral mesh.
- Sacado to form the source term from user-specified manufactured solution.
- Intrepid to build the discretization matrix and right-hand side.
- Tpetra to handle the global matrix and vector.
Poisson system:
div A grad u = f in Omega
u = g on Gamma
where
A is a symmetric, positive definite material tensor
f is a given source term
Corresponding discrete linear system for nodal coefficients(x):
Kx = b
K - HGrad stiffness matrix
b - right hand side vector
- Author
- Created by P. Bochev, D. Ridzal, K. Peterson, D. Hensinger, C. Siefert. Converted to Tpetra by I. Kalashnikova (ikala.nosp@m.sh@s.nosp@m.andia.nosp@m..gov)
template<class ArrayOut , class ArrayIn >
void evaluateExactSolution |
( |
ArrayOut & |
exactSolutionValues, |
|
|
const ArrayIn & |
evaluationPoints |
|
) |
| |
Computation of the exact solution at array of points in physical space.
- Parameters
-
exactSolutionValues | [out] Rank-2 (C,P) array with the values of the exact solution |
evaluationPoints | [in] Rank-3 (C,P,D) array with the evaluation points in physical frame |
template<class ArrayOut , class ArrayIn >
void evaluateExactSolutionGrad |
( |
ArrayOut & |
exactSolutionGradValues, |
|
|
const ArrayIn & |
evaluationPoints |
|
) |
| |
Computation of the gradient of the exact solution at array of points in physical space.
- Parameters
-
exactSolutionGradValues | [out] Rank-3 (C,P,D) array with the values of the gradient of the exact solution |
evaluationPoints | [in] Rank-3 (C,P,D) array with the evaluation points in physical frame |
template<class ArrayOut , class ArrayIn >
void evaluateMaterialTensor |
( |
ArrayOut & |
worksetMaterialValues, |
|
|
const ArrayIn & |
evaluationPoints |
|
) |
| |
Computation of the material tensor at array of points in physical space.
- Parameters
-
worksetMaterialValues | [out] Rank-2, 3 or 4 array with dimensions (C,P), (C,P,D) or (C,P,D,D) with the values of the material tensor |
evaluationPoints | [in] Rank-3 (C,P,D) array with the evaluation points in physical frame |
template<class ArrayOut , class ArrayIn >
void evaluateSourceTerm |
( |
ArrayOut & |
sourceTermValues, |
|
|
const ArrayIn & |
evaluationPoints |
|
) |
| |
Computation of the source term at array of points in physical space.
- Parameters
-
sourceTermValues | [out] Rank-2 (C,P) array with the values of the source term |
evaluationPoints | [in] Rank-3 (C,P,D) array with the evaluation points in physical frame |
template<typename Scalar >
const Scalar exactSolution |
( |
const Scalar & |
x, |
|
|
const Scalar & |
y, |
|
|
const Scalar & |
z |
|
) |
| |
User-defined exact solution.
- Parameters
-
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |
- Returns
- Value of the exact solution at (x,y,z)
template<typename Scalar >
void exactSolutionGrad |
( |
Scalar |
gradExact[3], |
|
|
const Scalar & |
x, |
|
|
const Scalar & |
y, |
|
|
const Scalar & |
z |
|
) |
| |
Computes gradient of the exact solution. Requires user-defined exact solution.
- Parameters
-
gradExact | [out] gradient of the exact solution evaluated at (x,y,z) |
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |
template<typename Scalar >
void materialTensor |
( |
Scalar |
material[][3], |
|
|
const Scalar & |
x, |
|
|
const Scalar & |
y, |
|
|
const Scalar & |
z |
|
) |
| |
User-defined material tensor.
- Parameters
-
material | [out] 3 x 3 material tensor evaluated at (x,y,z) |
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |
- Warning
- Symmetric and positive definite tensor is required for every (x,y,z).
template<typename Scalar >
const Scalar sourceTerm |
( |
Scalar & |
x, |
|
|
Scalar & |
y, |
|
|
Scalar & |
z |
|
) |
| |
Computes source term: f = -div(A.grad u). Requires user-defined exact solution and material tensor.
- Parameters
-
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |
- Returns
- Source term corresponding to the user-defined exact solution evaluated at (x,y,z)