51 template<
class Scalar,
class ArrayScalar>
54 this -> basisCardinality_ = 12;
55 this -> basisDegree_ = 1;
56 this -> basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<8> >() );
57 this -> basisType_ = BASIS_FEM_DEFAULT;
58 this -> basisCoordinates_ = COORDINATES_CARTESIAN;
59 this -> basisTagsAreSet_ =
false;
62 template<
class Scalar,
class ArrayScalar>
87 Intrepid::setOrdinalTagData(
this -> tagToOrdinal_,
88 this -> ordinalToTag_,
90 this -> basisCardinality_,
99 template<
class Scalar,
class ArrayScalar>
101 const ArrayScalar & inputPoints,
102 const EOperator operatorType)
const {
105 #ifdef HAVE_INTREPID_DEBUG
106 Intrepid::getValues_HCURL_Args<Scalar, ArrayScalar>(outputValues,
109 this -> getBaseCellTopology(),
110 this -> getCardinality() );
114 int dim0 = inputPoints.dimension(0);
121 switch (operatorType) {
123 for (
int i0 = 0; i0 < dim0; i0++) {
124 x = inputPoints(i0, 0);
125 y = inputPoints(i0, 1);
126 z = inputPoints(i0, 2);
129 outputValues(0, i0, 0) = (1.0 - y)*(1.0 - z)/8.0;
130 outputValues(0, i0, 1) = 0.0;
131 outputValues(0, i0, 2) = 0.0;
133 outputValues(1, i0, 0) = 0.0;
134 outputValues(1, i0, 1) = (1.0 + x)*(1.0 - z)/8.0;
135 outputValues(1, i0, 2) = 0.0;
137 outputValues(2, i0, 0) = -(1.0 + y)*(1.0 - z)/8.0;
138 outputValues(2, i0, 1) = 0.0;
139 outputValues(2, i0, 2) = 0.0;
141 outputValues(3, i0, 0) = 0.0;
142 outputValues(3, i0, 1) = -(1.0 - x)*(1.0 - z)/8.0;
143 outputValues(3, i0, 2) = 0.0;
145 outputValues(4, i0, 0) = (1.0 - y)*(1.0 + z)/8.0;
146 outputValues(4, i0, 1) = 0.0;
147 outputValues(4, i0, 2) = 0.0;
149 outputValues(5, i0, 0) = 0.0;
150 outputValues(5, i0, 1) = (1.0 + x)*(1.0 + z)/8.0;
151 outputValues(5, i0, 2) = 0.0;
153 outputValues(6, i0, 0) = -(1.0 + y)*(1.0 + z)/8.0;
154 outputValues(6, i0, 1) = 0.0;
155 outputValues(6, i0, 2) = 0.0;
157 outputValues(7, i0, 0) = 0.0;
158 outputValues(7, i0, 1) = -(1.0 - x)*(1.0 + z)/8.0;
159 outputValues(7, i0, 2) = 0.0;
161 outputValues(8, i0, 0) = 0.0;
162 outputValues(8, i0, 1) = 0.0;
163 outputValues(8, i0, 2) = (1.0 - x)*(1.0 - y)/8.0;
165 outputValues(9, i0, 0) = 0.0;
166 outputValues(9, i0, 1) = 0.0;
167 outputValues(9, i0, 2) = (1.0 + x)*(1.0 - y)/8.0;
169 outputValues(10, i0, 0) = 0.0;
170 outputValues(10, i0, 1) = 0.0;
171 outputValues(10, i0, 2) = (1.0 + x)*(1.0 + y)/8.0;
173 outputValues(11, i0, 0) = 0.0;
174 outputValues(11, i0, 1) = 0.0;
175 outputValues(11, i0, 2) = (1.0 - x)*(1.0 + y)/8.0;
180 for (
int i0 = 0; i0 < dim0; i0++) {
181 x = inputPoints(i0, 0);
182 y = inputPoints(i0, 1);
183 z = inputPoints(i0, 2);
186 outputValues(0, i0, 0) = 0.0;
187 outputValues(0, i0, 1) = -(1.0 - y)/8.0;
188 outputValues(0, i0, 2) = (1.0 - z)/8.0;
190 outputValues(1, i0, 0) = (1.0 + x)/8.0;
191 outputValues(1, i0, 1) = 0.0;
192 outputValues(1, i0, 2) = (1.0 - z)/8.0;
194 outputValues(2, i0, 0) = 0.0;
195 outputValues(2, i0, 1) = (1.0 + y)/8.0;
196 outputValues(2, i0, 2) = (1.0 - z)/8.0;
198 outputValues(3, i0, 0) = -(1.0 - x)/8.0;
199 outputValues(3, i0, 1) = 0.0;
200 outputValues(3, i0, 2) = (1.0 - z)/8.0;
202 outputValues(4, i0, 0) = 0.0;
203 outputValues(4, i0, 1) = (1.0 - y)/8.0;
204 outputValues(4, i0, 2) = (1.0 + z)/8.0;
206 outputValues(5, i0, 0) = -(1.0 + x)/8.0;
207 outputValues(5, i0, 1) = 0.0;
208 outputValues(5, i0, 2) = (1.0 + z)/8.0;
210 outputValues(6, i0, 0) = 0.0;
211 outputValues(6, i0, 1) = -(1.0 + y)/8.0;
212 outputValues(6, i0, 2) = (1.0 + z)/8.0;
214 outputValues(7, i0, 0) = (1.0 - x)/8.0;
215 outputValues(7, i0, 1) = 0.0;
216 outputValues(7, i0, 2) = (1.0 + z)/8.0;
218 outputValues(8, i0, 0) = -(1.0 - x)/8.0;
219 outputValues(8, i0, 1) = (1.0 - y)/8.0;
220 outputValues(8, i0, 2) = 0.0;
222 outputValues(9, i0, 0) = -(1.0 + x)/8.0;
223 outputValues(9, i0, 1) = -(1.0 - y)/8.0;
224 outputValues(9, i0, 2) = 0.0;
226 outputValues(10, i0, 0) = (1.0 + x)/8.0;
227 outputValues(10, i0, 1) = -(1.0 + y)/8.0;
228 outputValues(10, i0, 2) = 0.0;
230 outputValues(11, i0, 0) = (1.0 - x)/8.0;
231 outputValues(11, i0, 1) = (1.0 + y)/8.0;
232 outputValues(11, i0, 2) = 0.0;
237 TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
238 ">>> ERROR (Basis_HCURL_HEX_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
242 TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_GRAD), std::invalid_argument,
243 ">>> ERROR (Basis_HCURL_HEX_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
256 TEUCHOS_TEST_FOR_EXCEPTION( ( (operatorType == OPERATOR_D1) ||
257 (operatorType == OPERATOR_D2) ||
258 (operatorType == OPERATOR_D3) ||
259 (operatorType == OPERATOR_D4) ||
260 (operatorType == OPERATOR_D5) ||
261 (operatorType == OPERATOR_D6) ||
262 (operatorType == OPERATOR_D7) ||
263 (operatorType == OPERATOR_D8) ||
264 (operatorType == OPERATOR_D9) ||
265 (operatorType == OPERATOR_D10) ),
266 std::invalid_argument,
267 ">>> ERROR (Basis_HCURL_HEX_I1_FEM): Invalid operator type");
271 TEUCHOS_TEST_FOR_EXCEPTION( ( (operatorType != OPERATOR_VALUE) &&
272 (operatorType != OPERATOR_GRAD) &&
273 (operatorType != OPERATOR_CURL) &&
274 (operatorType != OPERATOR_DIV) &&
275 (operatorType != OPERATOR_D1) &&
276 (operatorType != OPERATOR_D2) &&
277 (operatorType != OPERATOR_D3) &&
278 (operatorType != OPERATOR_D4) &&
279 (operatorType != OPERATOR_D5) &&
280 (operatorType != OPERATOR_D6) &&
281 (operatorType != OPERATOR_D7) &&
282 (operatorType != OPERATOR_D8) &&
283 (operatorType != OPERATOR_D9) &&
284 (operatorType != OPERATOR_D10) ),
285 std::invalid_argument,
286 ">>> ERROR (Basis_HCURL_HEX_I1_FEM): Invalid operator type");
292 template<
class Scalar,
class ArrayScalar>
294 const ArrayScalar & inputPoints,
295 const ArrayScalar & cellVertices,
296 const EOperator operatorType)
const {
297 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
298 ">>> ERROR (Basis_HCURL_HEX_I1_FEM): FEM Basis calling an FVD member function");
301 template<
class Scalar,
class ArrayScalar>
303 #ifdef HAVE_INTREPID_DEBUG
305 TEUCHOS_TEST_FOR_EXCEPTION( !(DofCoords.rank() == 2), std::invalid_argument,
306 ">>> ERROR: (Intrepid::Basis_HCURL_QUAD_I1_FEM::getDofCoords) rank = 2 required for DofCoords array");
308 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(0) ==
this -> basisCardinality_ ), std::invalid_argument,
309 ">>> ERROR: (Intrepid::Basis_HCURL_QUAD_I1_FEM::getDofCoords) mismatch in number of DoF and 0th dimension of DofCoords array");
311 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(1) == (int)(
this -> basisCellTopology_.getDimension()) ), std::invalid_argument,
312 ">>> ERROR: (Intrepid::Basis_HCURL_QUAD_I1_FEM::getDofCoords) incorrect reference cell (1st) dimension in DofCoords array");
315 DofCoords(0,0) = 0.0; DofCoords(0,1) = -1.0; DofCoords(0,2) = -1.0;
316 DofCoords(1,0) = 1.0; DofCoords(1,1) = 0.0; DofCoords(1,2) = -1.0;
317 DofCoords(2,0) = 0.0; DofCoords(2,1) = 1.0; DofCoords(2,2) = -1.0;
318 DofCoords(3,0) = -1.0; DofCoords(3,1) = 0.0; DofCoords(3,2) = -1.0;
319 DofCoords(4,0) = 0.0; DofCoords(4,1) = -1.0; DofCoords(4,2) = 1.0;
320 DofCoords(5,0) = 1.0; DofCoords(5,1) = 0.0; DofCoords(5,2) = 1.0;
321 DofCoords(6,0) = 0.0; DofCoords(6,1) = 1.0; DofCoords(6,2) = 1.0;
322 DofCoords(7,0) = -1.0; DofCoords(7,1) = 0.0; DofCoords(7,2) = 1.0;
323 DofCoords(8,0) = -1.0; DofCoords(8,1) = -1.0; DofCoords(8,2) = 0.0;
324 DofCoords(9,0) = 1.0; DofCoords(9,1) = -1.0; DofCoords(9,2) = 0.0;
325 DofCoords(10,0) = 1.0; DofCoords(10,1) = 1.0; DofCoords(10,2) = 0.0;
326 DofCoords(11,0) = -1.0; DofCoords(11,1) = 1.0; DofCoords(11,2) = 0.0;
void getDofCoords(ArrayScalar &DofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on a reference Quadrilateral.
void initializeTags()
Initializes tagToOrdinal_ and ordinalToTag_ lookup arrays.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Hexahedron cell.
Basis_HCURL_HEX_I1_FEM()
Constructor.