1 #ifndef INTREPID_HGRAD_TRI_C2_FEMDEF_HPP
2 #define INTREPID_HGRAD_TRI_C2_FEMDEF_HPP
54 template<
class Scalar,
class ArrayScalar>
57 this -> basisCardinality_ = 6;
58 this -> basisDegree_ = 2;
59 this -> basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
60 this -> basisType_ = BASIS_FEM_DEFAULT;
61 this -> basisCoordinates_ = COORDINATES_CARTESIAN;
62 this -> basisTagsAreSet_ =
false;
67 template<
class Scalar,
class ArrayScalar>
77 int tags[] = { 0, 0, 0, 1,
85 Intrepid::setOrdinalTagData(
this -> tagToOrdinal_,
86 this -> ordinalToTag_,
88 this -> basisCardinality_,
97 template<
class Scalar,
class ArrayScalar>
99 const ArrayScalar & inputPoints,
100 const EOperator operatorType)
const {
103 #ifdef HAVE_INTREPID_DEBUG
104 Intrepid::getValues_HGRAD_Args<Scalar, ArrayScalar>(outputValues,
107 this -> getBaseCellTopology(),
108 this -> getCardinality() );
112 int dim0 = inputPoints.dimension(0);
118 switch (operatorType) {
121 for (
int i0 = 0; i0 < dim0; i0++) {
122 x = inputPoints(i0, 0);
123 y = inputPoints(i0, 1);
126 outputValues(0, i0) = (x + y - 1.0)*(2.0*x + 2.0*y - 1.0);
127 outputValues(1, i0) = x*(2.0*x - 1.0);
128 outputValues(2, i0) = y*(2.0*y - 1.0);
129 outputValues(3, i0) = -4.0*x*(x + y - 1.0);
130 outputValues(4, i0) = 4.0*x*y;
131 outputValues(5, i0) = -4.0*y*(x + y - 1.0);
138 for (
int i0 = 0; i0 < dim0; i0++) {
139 x = inputPoints(i0, 0);
140 y = inputPoints(i0, 1);
143 outputValues(0, i0, 0) = 4.0*x + 4.0*y - 3.0;
144 outputValues(0, i0, 1) = 4.0*x + 4.0*y - 3.0;
146 outputValues(1, i0, 0) = 4.0*x - 1.0;
147 outputValues(1, i0, 1) = 0.0;
149 outputValues(2, i0, 0) = 0.0;
150 outputValues(2, i0, 1) = 4.0*y - 1.0;
152 outputValues(3, i0, 0) = -4.0*(2.0*x + y - 1.0);
153 outputValues(3, i0, 1) = -4.0*x;
155 outputValues(4, i0, 0) = 4.0*y;
156 outputValues(4, i0, 1) = 4.0*x;
158 outputValues(5, i0, 0) = -4.0*y;
159 outputValues(5, i0, 1) = -4.0*(x + 2.0*y - 1.0);
164 for (
int i0 = 0; i0 < dim0; i0++) {
165 x = inputPoints(i0, 0);
166 y = inputPoints(i0, 1);
169 outputValues(0, i0, 1) =-(4.0*x + 4.0*y - 3.0);
170 outputValues(0, i0, 0) = 4.0*x + 4.0*y - 3.0;
172 outputValues(1, i0, 1) =-(4.0*x - 1.0);
173 outputValues(1, i0, 0) = 0.0;
175 outputValues(2, i0, 1) = 0.0;
176 outputValues(2, i0, 0) = 4.0*y - 1.0;
178 outputValues(3, i0, 1) = 4.0*(2.0*x + y - 1.0);
179 outputValues(3, i0, 0) = -4.0*x;
181 outputValues(4, i0, 1) = -4.0*y;
182 outputValues(4, i0, 0) = 4.0*x;
184 outputValues(5, i0, 1) = 4.0*y;
185 outputValues(5, i0, 0) = -4.0*(x + 2.0*y - 1.0);
190 TEUCHOS_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
191 ">>> ERROR (Basis_HGRAD_TRI_C2_FEM): DIV is invalid operator for rank-0 (scalar) fields in 2D.");
195 for (
int i0 = 0; i0 < dim0; i0++) {
198 outputValues(0, i0, 0) = 4.0;
199 outputValues(1, i0, 0) = 4.0;
200 outputValues(2, i0, 0) = 0.0;
201 outputValues(3, i0, 0) =-8.0;
202 outputValues(4, i0, 0) = 0.0;
203 outputValues(5, i0, 0) = 0.0;
206 outputValues(0, i0, 1) = 4.0;
207 outputValues(1, i0, 1) = 0.0;
208 outputValues(2, i0, 1) = 0.0;
209 outputValues(3, i0, 1) =-4.0;
210 outputValues(4, i0, 1) = 4.0;
211 outputValues(5, i0, 1) =-4.0;
214 outputValues(0, i0, 2) = 4.0;
215 outputValues(1, i0, 2) = 0.0;
216 outputValues(2, i0, 2) = 4.0;
217 outputValues(3, i0, 2) = 0.0;
218 outputValues(4, i0, 2) = 0.0;
219 outputValues(5, i0, 2) =-8.0;
233 int DkCardinality = Intrepid::getDkCardinality(operatorType,
234 this -> basisCellTopology_.getDimension() );
235 for(
int dofOrd = 0; dofOrd <
this -> basisCardinality_; dofOrd++) {
236 for (
int i0 = 0; i0 < dim0; i0++) {
237 for(
int dkOrd = 0; dkOrd < DkCardinality; dkOrd++){
238 outputValues(dofOrd, i0, dkOrd) = 0.0;
246 TEUCHOS_TEST_FOR_EXCEPTION( !( Intrepid::isValidOperator(operatorType) ), std::invalid_argument,
247 ">>> ERROR (Basis_HGRAD_TRI_C2_FEM): Invalid operator type");
253 template<
class Scalar,
class ArrayScalar>
255 const ArrayScalar & inputPoints,
256 const ArrayScalar & cellVertices,
257 const EOperator operatorType)
const {
258 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
259 ">>> ERROR (Basis_HGRAD_TRI_C2_FEM): FEM Basis calling an FVD member function");
262 template<
class Scalar,
class ArrayScalar>
264 #ifdef HAVE_INTREPID_DEBUG
266 TEUCHOS_TEST_FOR_EXCEPTION( !(DofCoords.rank() == 2), std::invalid_argument,
267 ">>> ERROR: (Intrepid::Basis_HGRAD_TRI_C2_FEM::getDofCoords) rank = 2 required for DofCoords array");
269 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(0) ==
this -> basisCardinality_ ), std::invalid_argument,
270 ">>> ERROR: (Intrepid::Basis_HGRAD_TRI_C2_FEM::getDofCoords) mismatch in number of DoF and 0th dimension of DofCoords array");
272 TEUCHOS_TEST_FOR_EXCEPTION( !( DofCoords.dimension(1) == (int)(
this -> basisCellTopology_.getDimension()) ), std::invalid_argument,
273 ">>> ERROR: (Intrepid::Basis_HGRAD_TRI_C2_FEM::getDofCoords) incorrect reference cell (1st) dimension in DofCoords array");
276 DofCoords(0,0) = 0.0; DofCoords(0,1) = 0.0;
277 DofCoords(1,0) = 1.0; DofCoords(1,1) = 0.0;
278 DofCoords(2,0) = 0.0; DofCoords(2,1) = 1.0;
279 DofCoords(3,0) = 0.5; DofCoords(3,1) = 0.0;
280 DofCoords(4,0) = 0.5; DofCoords(4,1) = 0.5;
281 DofCoords(5,0) = 0.0; DofCoords(5,1) = 0.5;
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 Triangle cell.
Basis_HGRAD_TRI_C2_FEM()
Constructor.
void getDofCoords(ArrayScalar &DofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on a reference Quadrilateral.