1 #ifndef INTREPID_HGRAD_TET_C2_FEMDEF_HPP
2 #define INTREPID_HGRAD_TET_C2_FEMDEF_HPP
53 template<
class Scalar,
class ArrayScalar>
56 this -> basisCardinality_ = 10;
57 this -> basisDegree_ = 2;
58 this -> basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<4> >() );
59 this -> basisType_ = BASIS_FEM_DEFAULT;
60 this -> basisCoordinates_ = COORDINATES_CARTESIAN;
61 this -> basisTagsAreSet_ =
false;
65 template<
class Scalar,
class ArrayScalar>
75 int tags[] = { 0, 0, 0, 1,
88 Intrepid::setOrdinalTagData(
this -> tagToOrdinal_,
89 this -> ordinalToTag_,
91 this -> basisCardinality_,
100 template<
class Scalar,
class ArrayScalar>
102 const ArrayScalar & inputPoints,
103 const EOperator operatorType)
const {
106 #ifdef HAVE_INTREPID_DEBUG
107 Intrepid::getValues_HGRAD_Args<Scalar, ArrayScalar>(outputValues,
110 this -> getBaseCellTopology(),
111 this -> getCardinality() );
115 int dim0 = inputPoints.dimension(0);
122 switch (operatorType) {
125 for (
int i0 = 0; i0 < dim0; i0++) {
126 x = inputPoints(i0, 0);
127 y = inputPoints(i0, 1);
128 z = inputPoints(i0, 2);
131 outputValues(0, i0) = (-1. + x + y + z)*(-1. + 2.*x + 2.*y + 2.*z);
132 outputValues(1, i0) = x*(-1. + 2.*x);
133 outputValues(2, i0) = y*(-1. + 2.*y);
134 outputValues(3, i0) = z*(-1. + 2.*z);
136 outputValues(4, i0) = -4.*x*(-1. + x + y + z);
137 outputValues(5, i0) = 4.*x*y;
138 outputValues(6, i0) = -4.*y*(-1. + x + y + z);
139 outputValues(7, i0) = -4.*z*(-1. + x + y + z);
140 outputValues(8, i0) = 4.*x*z;
141 outputValues(9, i0) = 4.*y*z;
148 for (
int i0 = 0; i0 < dim0; i0++) {
149 x = inputPoints(i0, 0);
150 y = inputPoints(i0, 1);
151 z = inputPoints(i0, 2);
154 outputValues(0, i0, 0) = -3.+ 4.*x + 4.*y + 4.*z;
155 outputValues(0, i0, 1) = -3.+ 4.*x + 4.*y + 4.*z;
156 outputValues(0, i0, 2) = -3.+ 4.*x + 4.*y + 4.*z;
158 outputValues(1, i0, 0) = -1.+ 4.*x;
159 outputValues(1, i0, 1) = 0.;
160 outputValues(1, i0, 2) = 0.;
162 outputValues(2, i0, 0) = 0.;
163 outputValues(2, i0, 1) = -1.+ 4.*y;
164 outputValues(2, i0, 2) = 0.;
166 outputValues(3, i0, 0) = 0.;
167 outputValues(3, i0, 1) = 0.;
168 outputValues(3, i0, 2) = -1.+ 4.*z;
171 outputValues(4, i0, 0) = -4.*(-1.+ 2*x + y + z);
172 outputValues(4, i0, 1) = -4.*x;
173 outputValues(4, i0, 2) = -4.*x;
175 outputValues(5, i0, 0) = 4.*y;
176 outputValues(5, i0, 1) = 4.*x;
177 outputValues(5, i0, 2) = 0.;
179 outputValues(6, i0, 0) = -4.*y;
180 outputValues(6, i0, 1) = -4.*(-1.+ x + 2*y + z);
181 outputValues(6, i0, 2) = -4.*y;
183 outputValues(7, i0, 0) = -4.*z;
184 outputValues(7, i0, 1) = -4.*z;
185 outputValues(7, i0, 2) = -4.*(-1.+ x + y + 2*z);
187 outputValues(8, i0, 0) = 4.*z;
188 outputValues(8, i0, 1) = 0.;
189 outputValues(8, i0, 2) = 4.*x;
191 outputValues(9, i0, 0) = 0.;
192 outputValues(9, i0, 1) = 4.*z;
193 outputValues(9, i0, 2) = 4.*y;
199 TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_CURL), std::invalid_argument,
200 ">>> ERROR (Basis_HGRAD_TET_C2_FEM): CURL is invalid operator for rank-0 (scalar) functions in 3D");
204 TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
205 ">>> ERROR (Basis_HGRAD_TET_C2_FEM): DIV is invalid operator for rank-0 (scalar) functions in 3D");
209 for (
int i0 = 0; i0 < dim0; i0++) {
211 outputValues(0, i0, 0) = 4.;
212 outputValues(0, i0, 1) = 4.;
213 outputValues(0, i0, 2) = 4.;
214 outputValues(0, i0, 3) = 4.;
215 outputValues(0, i0, 4) = 4.;
216 outputValues(0, i0, 5) = 4.;
218 outputValues(1, i0, 0) = 4.;
219 outputValues(1, i0, 1) = 0.;
220 outputValues(1, i0, 2) = 0.;
221 outputValues(1, i0, 3) = 0.;
222 outputValues(1, i0, 4) = 0.;
223 outputValues(1, i0, 5) = 0.;
225 outputValues(2, i0, 0) = 0.;
226 outputValues(2, i0, 1) = 0.;
227 outputValues(2, i0, 2) = 0.;
228 outputValues(2, i0, 3) = 4.;
229 outputValues(2, i0, 4) = 0.;
230 outputValues(2, i0, 5) = 0.;
232 outputValues(3, i0, 0) = 0.;
233 outputValues(3, i0, 1) = 0.;
234 outputValues(3, i0, 2) = 0.;
235 outputValues(3, i0, 3) = 0.;
236 outputValues(3, i0, 4) = 0.;
237 outputValues(3, i0, 5) = 4.;
239 outputValues(4, i0, 0) = -8.;
240 outputValues(4, i0, 1) = -4.;
241 outputValues(4, i0, 2) = -4.;
242 outputValues(4, i0, 3) = 0.;
243 outputValues(4, i0, 4) = 0.;
244 outputValues(4, i0, 5) = 0.;
246 outputValues(5, i0, 0) = 0.;
247 outputValues(5, i0, 1) = 4.;
248 outputValues(5, i0, 2) = 0.;
249 outputValues(5, i0, 3) = 0.;
250 outputValues(5, i0, 4) = 0.;
251 outputValues(5, i0, 5) = 0.;
253 outputValues(6, i0, 0) = 0.;
254 outputValues(6, i0, 1) = -4.;
255 outputValues(6, i0, 2) = 0.;
256 outputValues(6, i0, 3) = -8.;
257 outputValues(6, i0, 4) = -4.;
258 outputValues(6, i0, 5) = 0;
260 outputValues(7, i0, 0) = 0.;
261 outputValues(7, i0, 1) = 0.;
262 outputValues(7, i0, 2) = -4.;
263 outputValues(7, i0, 3) = 0.;
264 outputValues(7, i0, 4) = -4.;
265 outputValues(7, i0, 5) = -8.;
267 outputValues(8, i0, 0) = 0.;
268 outputValues(8, i0, 1) = 0.;
269 outputValues(8, i0, 2) = 4.;
270 outputValues(8, i0, 3) = 0.;
271 outputValues(8, i0, 4) = 0.;
272 outputValues(8, i0, 5) = 0.;
274 outputValues(9, i0, 0) = 0.;
275 outputValues(9, i0, 1) = 0.;
276 outputValues(9, i0, 2) = 0.;
277 outputValues(9, i0, 3) = 0.;
278 outputValues(9, i0, 4) = 4.;
279 outputValues(9, i0, 5) = 0.;
293 int DkCardinality = Intrepid::getDkCardinality(operatorType,
294 this -> basisCellTopology_.getDimension() );
295 for(
int dofOrd = 0; dofOrd <
this -> basisCardinality_; dofOrd++) {
296 for (
int i0 = 0; i0 < dim0; i0++) {
297 for(
int dkOrd = 0; dkOrd < DkCardinality; dkOrd++){
298 outputValues(dofOrd, i0, dkOrd) = 0.0;
305 TEUCHOS_TEST_FOR_EXCEPTION( !( Intrepid::isValidOperator(operatorType) ), std::invalid_argument,
306 ">>> ERROR (Basis_HGRAD_TET_C2_FEM): Invalid operator type");
312 template<
class Scalar,
class ArrayScalar>
314 const ArrayScalar & inputPoints,
315 const ArrayScalar & cellVertices,
316 const EOperator operatorType)
const {
317 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
318 ">>> ERROR (Basis_HGRAD_TET_C2_FEM): FEM Basis calling an FVD member function");
323 template<
class Scalar,
class ArrayScalar>
325 #ifdef HAVE_INTREPID_DEBUG
327 TEUCHOS_TEST_FOR_EXCEPTION( !(DofCoords.rank() == 2), std::invalid_argument,
328 ">>> ERROR: (Intrepid::Basis_HGRAD_TET_C2_FEM::getDofCoords) rank = 2 required for DofCoords array");
330 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(0) ==
this -> basisCardinality_ ), std::invalid_argument,
331 ">>> ERROR: (Intrepid::Basis_HGRAD_TET_C2_FEM::getDofCoords) mismatch in number of DoF and 0th dimension of DofCoords array");
333 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(1) == (int)(
this -> basisCellTopology_.getDimension()) ), std::invalid_argument,
334 ">>> ERROR: (Intrepid::Basis_HGRAD_TET_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in DofCoords array");
337 DofCoords(0,0) = 0.0; DofCoords(0,1) = 0.0; DofCoords(0,2) = 0.0;
338 DofCoords(1,0) = 1.0; DofCoords(1,1) = 0.0; DofCoords(1,2) = 0.0;
339 DofCoords(2,0) = 0.0; DofCoords(2,1) = 1.0; DofCoords(2,2) = 0.0;
340 DofCoords(3,0) = 0.0; DofCoords(3,1) = 0.0; DofCoords(3,2) = 1.0;
342 DofCoords(4,0) = 0.5; DofCoords(4,1) = 0.0; DofCoords(4,2) = 0.0;
343 DofCoords(5,0) = 0.5; DofCoords(5,1) = 0.5; DofCoords(5,2) = 0.0;
344 DofCoords(6,0) = 0.0; DofCoords(6,1) = 0.5; DofCoords(6,2) = 0.0;
345 DofCoords(7,0) = 0.0; DofCoords(7,1) = 0.0; DofCoords(7,2) = 0.5;
346 DofCoords(8,0) = 0.5; DofCoords(8,1) = 0.0; DofCoords(8,2) = 0.5;
347 DofCoords(9,0) = 0.0; DofCoords(9,1) = 0.5; DofCoords(9,2) = 0.5;
353 #if defined(Intrepid_SHOW_DEPRECATED_WARNINGS)
355 #warning "The Intrepid package is deprecated"
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
void getDofCoords(ArrayScalar &DofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on a reference Tetrahedron.
Basis_HGRAD_TET_C2_FEM()
Constructor.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
FEM basis evaluation on a reference Tetrahedron cell.