1 #ifndef INTREPID_HGRAD_TET_CN_FEMDEF_HPP
2 #define INTREPID_HGRAD_TET_CN_FEMDEF_HPP
53 template<
class Scalar,
class ArrayScalar>
55 const EPointType pointType ):
57 V((n+1)*(n+2)*(n+3)/6,(n+1)*(n+2)*(n+3)/6),
58 Vinv((n+1)*(n+2)*(n+3)/6,(n+1)*(n+2)*(n+3)/6),
59 latticePts( (n+1)*(n+2)*(n+3)/6 , 3 )
61 const int N = (n+1)*(n+2)*(n+3)/6;
64 this ->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
71 PointTools::getLattice<Scalar,FieldContainer<Scalar> >(
latticePts ,
81 Phis.getValues(
V , latticePts , OPERATOR_VALUE );
84 Teuchos::SerialDenseMatrix<int,Scalar> Vsdm(N,N);
85 for (
int i=0;i<N;i++) {
86 for (
int j=0;j<N;j++) {
92 Teuchos::SerialDenseSolver<int,Scalar> solver;
93 solver.setMatrix( rcp( &Vsdm ,
false ) );
97 for (
int i=0;i<N;i++) {
98 for (
int j=0;j<N;j++) {
99 Vinv(i,j) = Vsdm(j,i);
106 template<
class Scalar,
class ArrayScalar>
117 int *tags =
new int[ tagSize * this->getCardinality() ];
119 const int degree = this->getDegree();
120 const int numEdgeDof = degree - 1;
121 const int numFaceDof =
PointTools::getLatticeSize( shards::CellTopology( shards::getCellTopologyData<shards::Triangle<3> >() ) ,
127 int edge_dof_cur[] = {0,0,0,0,0,0};
128 int face_dof_cur[] = {0,0,0,0};
129 int cell_dof_cur = 0;
134 tag_cur[0] = 0; tag_cur[1] = 0; tag_cur[2] = 0; tag_cur[3] = 1;
139 for (
int i=1;i<degree;i++) {
140 tag_cur[0] = 1; tag_cur[1] = 0; tag_cur[2] = edge_dof_cur[0]; tag_cur[3] = numEdgeDof;
147 tag_cur[0] = 0; tag_cur[1] = 1; tag_cur[2] = 0; tag_cur[3] = 1;
152 for (
int i=1;i<degree;i++) {
154 tag_cur[0] = 1; tag_cur[1] = 2; tag_cur[2] = edge_dof_cur[2]; tag_cur[3] = numEdgeDof;
160 for (
int j=1;j<degree-i;j++) {
161 tag_cur[0] = 2; tag_cur[1] = 3; tag_cur[2] = face_dof_cur[3]; tag_cur[3] = numFaceDof;
168 tag_cur[0] = 1; tag_cur[1] = 1; tag_cur[2] = edge_dof_cur[1]; tag_cur[3] = numEdgeDof;
176 tag_cur[0] = 0; tag_cur[1] = 2; tag_cur[2] = 0; tag_cur[3] = 1;
183 for (
int i=1;i<degree;i++) {
186 tag_cur[0] = 1; tag_cur[1] = 3; tag_cur[2] = edge_dof_cur[3]; tag_cur[3] = numEdgeDof;
191 for (
int j=1;j<degree-i;j++) {
192 tag_cur[0] = 2; tag_cur[1] = 0; tag_cur[2] = face_dof_cur[0]; tag_cur[3] = numFaceDof;
198 tag_cur[0] = 1; tag_cur[1] = 4; tag_cur[2] = edge_dof_cur[4]; tag_cur[3] = numEdgeDof;
205 for (
int j=1;j<degree-i;j++) {
207 tag_cur[0] = 2; tag_cur[1] = 2; tag_cur[2] = face_dof_cur[2]; tag_cur[3] = numFaceDof;
212 for (
int k=1;k<degree-i-j;k++) {
213 tag_cur[0] = 3; tag_cur[1] = 0; tag_cur[2] = cell_dof_cur; tag_cur[3] = numCellDof;
219 tag_cur[0] = 2; tag_cur[1] = 1; tag_cur[2] = face_dof_cur[1]; tag_cur[3] = numFaceDof;
226 tag_cur[0] = 1; tag_cur[1] = 5; tag_cur[2] = edge_dof_cur[5]; tag_cur[3] = numEdgeDof;
234 tag_cur[0] = 0; tag_cur[1] = 3; tag_cur[2] = 0; tag_cur[3] = 1;
241 Intrepid::setOrdinalTagData(
this -> tagToOrdinal_,
242 this -> ordinalToTag_,
244 this -> basisCardinality_,
256 template<
class Scalar,
class ArrayScalar>
258 const ArrayScalar & inputPoints,
259 const EOperator operatorType)
const {
262 #ifdef HAVE_INTREPID_DEBUG
263 Intrepid::getValues_HGRAD_Args<Scalar, ArrayScalar>(outputValues,
266 this -> getBaseCellTopology(),
267 this -> getCardinality() );
269 const int numPts = inputPoints.dimension(0);
270 const int numBf = this->getCardinality();
273 switch (operatorType) {
277 Phis.getValues( phisCur , inputPoints , operatorType );
278 for (
int i=0;i<outputValues.dimension(0);i++) {
279 for (
int j=0;j<outputValues.dimension(1);j++) {
280 outputValues(i,j) = 0.0;
281 for (
int k=0;k<this->getCardinality();k++) {
282 outputValues(i,j) += this->Vinv(k,i) * phisCur(k,j);
301 (operatorType == OPERATOR_GRAD)? getDkCardinality(OPERATOR_D1,3): getDkCardinality(operatorType,3);
304 Phis.getValues( phisCur , inputPoints , operatorType );
306 for (
int i=0;i<outputValues.dimension(0);i++) {
307 for (
int j=0;j<outputValues.dimension(1);j++) {
308 for (
int k=0;k<outputValues.dimension(2);k++) {
309 outputValues(i,j,k) = 0.0;
310 for (
int l=0;l<this->getCardinality();l++) {
311 outputValues(i,j,k) += this->Vinv(l,i) * phisCur(l,j,k);
320 TEUCHOS_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
321 ">>> ERROR (Basis_HGRAD_TET_Cn_FEM): Operator type not implemented");
325 catch (std::invalid_argument &exception){
326 TEUCHOS_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
327 ">>> ERROR (Basis_HGRAD_TET_Cn_FEM): Operator type not implemented");
334 template<
class Scalar,
class ArrayScalar>
336 const ArrayScalar & inputPoints,
337 const ArrayScalar & cellVertices,
338 const EOperator operatorType)
const {
339 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
340 ">>> ERROR (Basis_HGRAD_TET_Cn_FEM): FEM Basis calling an FVD member function");
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
Basis_HGRAD_TET_Cn_FEM_ORTH< Scalar, FieldContainer< Scalar > > Phis
The orthogonal basis on triangles, out of which the nodal basis is constructed.
EBasis basisType_
Type of the basis.
Basis_HGRAD_TET_Cn_FEM(const int n, const EPointType pointType)
Constructor.
bool basisTagsAreSet_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
FieldContainer< Scalar > Vinv
The inverse of V. The columns of Vinv express the Lagrange basis in terms of the orthogonal basis...
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package http://trilinos...
FieldContainer< Scalar > latticePts
stores the points at which degrees of freedom are located.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
int basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
FieldContainer< Scalar > V
The Vandermonde matrix with V_{ij} = phi_i(x_j), where x_j is the j_th point in the lattice...
int basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...