Example solution of a Poisson equation on a hexahedral or tetrahedral mesh using nodal (Hgrad) elements. More...
#include <unistd.h>
#include <vector>
#include <map>
#include "TrilinosCouplings_config.h"
#include <Teuchos_CommandLineProcessor.hpp>
#include "Intrepid_FunctionSpaceTools.hpp"
#include "Intrepid_CellTools.hpp"
#include "Intrepid_ArrayTools.hpp"
#include "Intrepid_Basis.hpp"
#include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
#include "Intrepid_HGRAD_TET_C1_FEM.hpp"
#include "Intrepid_HGRAD_QUAD_C1_FEM.hpp"
#include "Intrepid_HGRAD_TRI_C1_FEM.hpp"
#include "Intrepid_RealSpaceTools.hpp"
#include "Intrepid_DefaultCubatureFactory.hpp"
#include "Intrepid_Utils.hpp"
#include "Epetra_Time.h"
#include "Epetra_Map.h"
#include "Epetra_MpiComm.h"
#include "Epetra_FECrsMatrix.h"
#include "Epetra_FEVector.h"
#include "Epetra_Import.h"
#include "Teuchos_oblackholestream.hpp"
#include "Teuchos_RCP.hpp"
#include "Teuchos_BLAS.hpp"
#include "Teuchos_GlobalMPISession.hpp"
#include "Teuchos_XMLParameterListHelpers.hpp"
#include "Shards_CellTopology.hpp"
#include "EpetraExt_RowMatrixOut.h"
#include "EpetraExt_MultiVectorOut.h"
#include "AztecOO.h"
#include "ml_MultiLevelPreconditioner.h"
#include "ml_epetra_utils.h"
#include "Sacado_No_Kokkos.hpp"
#include "Ionit_Initializer.h"
#include "Ioss_SubSystem.h"
#include "stk_io/IossBridge.hpp"
#include "stk_io/StkMeshIoBroker.hpp"
#include "stk_util/parallel/Parallel.hpp"
#include "stk_mesh/base/MetaData.hpp"
#include "stk_mesh/base/CoordinateSystems.hpp"
#include "stk_mesh/base/BulkData.hpp"
#include "stk_mesh/base/Comm.hpp"
#include "stk_mesh/base/Selector.hpp"
#include "stk_mesh/base/GetEntities.hpp"
#include "stk_mesh/base/GetBuckets.hpp"
#include "stk_mesh/base/CreateAdjacentEntities.hpp"
#include <stk_mesh/base/FEMHelpers.hpp>
#include <stk_mesh/base/Field.hpp>
#include "stk_mesh/base/Bucket.hpp"
#include "stk_mesh/base/Entity.hpp"
#include "stk_mesh/base/FieldBase.hpp"
#include "stk_mesh/base/Types.hpp"
Functions | |
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 | TestMultiLevelPreconditioner (char ProblemType[], Teuchos::ParameterList &MLList, Epetra_CrsMatrix &A, const Epetra_MultiVector &xexact, Epetra_MultiVector &b, Epetra_MultiVector &uh, Epetra_MultiVector &coords, double &TotalErrorResidual, double &TotalErrorExactSol) |
void | getBasis (Teuchos::RCP< Intrepid::Basis< double, IntrepidFieldContainer > > &basis, const ShardsCellTopology &cellTopology, int order) |
Simple factory that chooses basis function based on cell topology. More... | |
int | getDimension (const ShardsCellTopology &cellTopology) |
void | mesh_read_write (const std::string &type, const std::string &working_directory, const std::string &filename, stk::io::StkMeshIoBroker &broker, int db_integer_size, stk::io::HeartbeatType hb_type) |
int | main (int argc, char *argv[]) |
Variables | |
int | spaceDim |
Example solution of a Poisson equation on a hexahedral or tetrahedral mesh using nodal (Hgrad) elements.
This example requires a hexahedral or tetrahedral mesh in Exodus format with a nodeset containing boundary nodes. STK is used to read the mesh and populate a mesh database, Intrepid is used to build the stiffness matrix and right-hand side, and ML is used to solve the resulting linear system.
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
void evaluateExactSolution | ( | ArrayOut & | exactSolutionValues, |
const ArrayIn & | evaluationPoints | ||
) |
Computation of the exact solution at array of points in physical space.
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 |
void evaluateExactSolutionGrad | ( | ArrayOut & | exactSolutionGradValues, |
const ArrayIn & | evaluationPoints | ||
) |
Computation of the gradient of the exact solution at array of points in physical space.
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 |
void evaluateMaterialTensor | ( | ArrayOut & | worksetMaterialValues, |
const ArrayIn & | evaluationPoints | ||
) |
Computation of the material tensor at array of points in physical space.
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 |
void evaluateSourceTerm | ( | ArrayOut & | sourceTermValues, |
const ArrayIn & | evaluationPoints | ||
) |
Computation of the source term at array of points in physical space.
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 |
const Scalar exactSolution | ( | const Scalar & | x, |
const Scalar & | y, | ||
const Scalar & | z | ||
) |
User-defined exact solution.
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |
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.
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 |
void getBasis | ( | Teuchos::RCP< Intrepid::Basis< double, IntrepidFieldContainer > > & | basis, |
const ShardsCellTopology & | cellTopology, | ||
int | order | ||
) |
Simple factory that chooses basis function based on cell topology.
cellTopology | [in] Shards cell topology |
order | [in] basis function order, currently unused |
basis | [out] pointer to Intrepid basis |
void materialTensor | ( | Scalar | material[][3], |
const Scalar & | x, | ||
const Scalar & | y, | ||
const Scalar & | z | ||
) |
User-defined material tensor.
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 |
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.
x | [in] x-coordinate of the evaluation point |
y | [in] y-coordinate of the evaluation point |
z | [in] z-coordinate of the evaluation point |