12 #include "../../gridfunc.hpp"
28 return CEED_TOPOLOGY_LINE;
30 return CEED_TOPOLOGY_TRIANGLE;
32 return CEED_TOPOLOGY_QUAD;
34 return CEED_TOPOLOGY_TET;
36 return CEED_TOPOLOGY_HEX;
38 return CEED_TOPOLOGY_PRISM;
40 return CEED_TOPOLOGY_PYRAMID;
42 MFEM_ABORT(
"This type of element is not supported");
43 return CEED_TOPOLOGY_PRISM;
50 Ceed ceed, CeedBasis *basis)
55 const int ndofs = maps.
ndof;
56 const int nqpts = maps.
nqpt;
59 for (
int i = 0; i < nqpts; i++)
63 if (dim>1) { qX(1,i) = ip.
y; }
64 if (dim>2) { qX(2,i) = ip.
z; }
67 CeedBasisCreateH1(ceed, GetCeedTopology(fe.
GetGeomType()),
70 qX.GetData(), qW.GetData(), basis);
76 Ceed ceed, CeedBasis *basis)
80 const int ndofs = maps.
ndof;
81 const int nqpts = maps.
nqpt;
87 for (
int i = 0; i < nqpts; i++)
101 static void InitBasisImpl(
const FiniteElementSpace &fes,
102 const FiniteElement &fe,
103 const IntegrationRule &ir,
104 Ceed ceed, CeedBasis *basis)
107 const int P = fe.GetDof();
108 const int Q = ir.GetNPoints();
109 const int ncomp = fes.GetVDim();
110 BasisKey basis_key(&fes, &ir, ncomp, P, Q);
111 auto basis_itr = mfem::internal::ceed_basis_map.find(basis_key);
116 if (basis_itr == mfem::internal::ceed_basis_map.end())
120 InitTensorBasis(fes, fe, ir, ceed, basis);
124 InitNonTensorBasis(fes, fe, ir, ceed, basis);
126 mfem::internal::ceed_basis_map[basis_key] = *basis;
130 *basis = basis_itr->second;
136 Ceed ceed, CeedBasis *basis)
139 InitBasisImpl(fes, fe, ir, ceed, basis);
146 Ceed ceed, CeedBasis *basis)
149 InitBasisImpl(fes, fe, ir, ceed, basis);
void InitBasisWithIndices(const FiniteElementSpace &fes, const IntegrationRule &ir, int nelem, const int *indices, Ceed ceed, CeedBasis *basis)
Initialize a CeedBasis for mixed meshes.
Abstract class for all finite elements.
void InitBasis(const FiniteElementSpace &fes, const IntegrationRule &ir, Ceed ceed, CeedBasis *basis)
Initialize a CeedBasis for non-mixed meshes.
Class for an integration rule - an Array of IntegrationPoint.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
T * GetData()
Returns the data.
Data type dense matrix using column-major storage.
int nqpt
Number of quadrature points. When mode is TENSOR, this is the 1D number.
Array< double > Gt
Transpose of G.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
std::tuple< const mfem::FiniteElementSpace *, const mfem::IntegrationRule *, int, int, int > BasisKey
Mesh * GetMesh() const
Returns the mesh.
int GetVDim() const
Returns vector dimension.
Array< double > Bt
Transpose of B.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Structure representing the matrices/tensors needed to evaluate (in reference space) the values...
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Class for integration point with weight.
Full multidimensional representation which does not use tensor product structure. The ordering of the...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...