49 #ifndef INTREPID_CELLTOOLSDEF_HPP
50 #define INTREPID_CELLTOOLSDEF_HPP
53 #if defined (__clang__) && !defined (__INTEL_COMPILER)
54 #pragma clang system_header
59 template<
class Scalar>
61 const shards::CellTopology& parentCell){
63 #ifdef HAVE_INTREPID_DEBUG
64 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(parentCell) ), std::invalid_argument,
65 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): the specified cell topology does not have a reference cell.");
67 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <= subcellDim) && (subcellDim <= 2 ) ), std::invalid_argument,
68 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): parametrization defined only for 1 and 2-dimensional subcells.");
115 switch(parentCell.getKey() ) {
118 case shards::Tetrahedron<4>::key:
119 case shards::Tetrahedron<8>::key:
120 case shards::Tetrahedron<10>::key:
121 case shards::Tetrahedron<11>::key:
122 if(subcellDim == 2) {
124 setSubcellParametrization(tetFaces, subcellDim, parentCell);
129 else if(subcellDim == 1) {
131 setSubcellParametrization(tetEdges, subcellDim, parentCell);
137 TEUCHOS_TEST_FOR_EXCEPTION( (subcellDim != 1 || subcellDim != 2), std::invalid_argument,
138 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Tet parametrizations defined for 1 and 2-subcells only");
143 case shards::Hexahedron<8>::key:
144 case shards::Hexahedron<20>::key:
145 case shards::Hexahedron<27>::key:
146 if(subcellDim == 2) {
148 setSubcellParametrization(hexFaces, subcellDim, parentCell);
153 else if(subcellDim == 1) {
155 setSubcellParametrization(hexEdges, subcellDim, parentCell);
161 TEUCHOS_TEST_FOR_EXCEPTION( (subcellDim != 1 || subcellDim != 2), std::invalid_argument,
162 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Hex parametrizations defined for 1 and 2-subcells only");
167 case shards::Pyramid<5>::key:
168 case shards::Pyramid<13>::key:
169 case shards::Pyramid<14>::key:
170 if(subcellDim == 2) {
172 setSubcellParametrization(pyrFaces, subcellDim, parentCell);
177 else if(subcellDim == 1) {
179 setSubcellParametrization(pyrEdges, subcellDim, parentCell);
185 TEUCHOS_TEST_FOR_EXCEPTION( (subcellDim != 1 || subcellDim != 2), std::invalid_argument,
186 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Pyramid parametrizations defined for 1 and 2-subcells only");
191 case shards::Wedge<6>::key:
192 case shards::Wedge<15>::key:
193 case shards::Wedge<18>::key:
194 if(subcellDim == 2) {
196 setSubcellParametrization(wedgeFaces, subcellDim, parentCell);
201 else if(subcellDim == 1) {
203 setSubcellParametrization(wedgeEdges, subcellDim, parentCell);
209 TEUCHOS_TEST_FOR_EXCEPTION( (subcellDim != 1 || subcellDim != 2), std::invalid_argument,
210 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Wedge parametrization defined for 1 and 2-subcells only");
216 case shards::Triangle<3>::key:
217 case shards::Triangle<4>::key:
218 case shards::Triangle<6>::key:
219 if(subcellDim == 1) {
221 setSubcellParametrization(triEdges, subcellDim, parentCell);
227 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
228 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Triangle parametrizations defined for 1-subcells only");
232 case shards::Quadrilateral<4>::key:
233 case shards::Quadrilateral<8>::key:
234 case shards::Quadrilateral<9>::key:
235 if(subcellDim == 1) {
237 setSubcellParametrization(quadEdges, subcellDim, parentCell);
243 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
244 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Quad parametrizations defined for 1-subcells only");
251 case shards::ShellTriangle<3>::key:
252 case shards::ShellTriangle<6>::key:
253 if(subcellDim == 2) {
254 if(!shellTriFacesSet){
255 setSubcellParametrization(shellTriFaces, subcellDim, parentCell);
256 shellTriFacesSet = 1;
258 return shellTriFaces;
260 else if(subcellDim == 1) {
261 if(!shellTriEdgesSet){
262 setSubcellParametrization(shellTriEdges, subcellDim, parentCell);
263 shellTriEdgesSet = 1;
265 return shellTriEdges;
267 else if( subcellDim != 1 || subcellDim != 2){
268 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
269 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Shell Triangle parametrizations defined for 1 and 2-subcells only");
273 case shards::ShellQuadrilateral<4>::key:
274 case shards::ShellQuadrilateral<8>::key:
275 case shards::ShellQuadrilateral<9>::key:
276 if(subcellDim == 2) {
277 if(!shellQuadFacesSet){
278 setSubcellParametrization(shellQuadFaces, subcellDim, parentCell);
279 shellQuadFacesSet = 1;
281 return shellQuadFaces;
283 else if(subcellDim == 1) {
284 if(!shellQuadEdgesSet){
285 setSubcellParametrization(shellQuadEdges, subcellDim, parentCell);
286 shellQuadEdgesSet = 1;
288 return shellQuadEdges;
290 else if( subcellDim != 1 || subcellDim != 2){
291 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
292 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): Shell Quad parametrizations defined for 1 and 2-subcells only");
299 case shards::ShellLine<2>::key:
300 case shards::ShellLine<3>::key:
301 case shards::Beam<2>::key:
302 case shards::Beam<3>::key:
303 if(subcellDim == 1) {
305 setSubcellParametrization(lineEdges, subcellDim, parentCell);
311 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
312 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): shell line/beam parametrizations defined for 1-subcells only");
317 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
318 ">>> ERROR (Intrepid::CellTools::getSubcellParametrization): invalid cell topology.");
326 template<
class Scalar>
328 const int subcellDim,
329 const shards::CellTopology& parentCell)
331 #ifdef HAVE_INTREPID_DEBUG
332 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(parentCell) ), std::invalid_argument,
333 ">>> ERROR (Intrepid::CellTools::setSubcellParametrization): the specified cell topology does not have a reference cell.");
335 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <= subcellDim) && (subcellDim <= 2 ) ), std::invalid_argument,
336 ">>> ERROR (Intrepid::CellTools::setSubcellParametrization): parametrization defined only for 1 and 2-dimensional subcells.");
356 unsigned sc = parentCell.getSubcellCount(subcellDim);
357 unsigned pcd = parentCell.getDimension();
358 unsigned coeff = (subcellDim == 1) ? 2 : 3;
362 subcellParametrization.
resize(sc, pcd, coeff);
366 for(
unsigned subcellOrd = 0; subcellOrd < sc; subcellOrd++){
368 int v0ord = parentCell.getNodeMap(subcellDim, subcellOrd, 0);
369 int v1ord = parentCell.getNodeMap(subcellDim, subcellOrd, 1);
373 const double* v0 = getReferenceVertex(parentCell, v0ord);
374 const double* v1 = getReferenceVertex(parentCell, v1ord);
377 subcellParametrization(subcellOrd, 0, 0) = (v0[0] + v1[0])/2.0;
378 subcellParametrization(subcellOrd, 0, 1) = (v1[0] - v0[0])/2.0;
381 subcellParametrization(subcellOrd, 1, 0) = (v0[1] + v1[1])/2.0;
382 subcellParametrization(subcellOrd, 1, 1) = (v1[1] - v0[1])/2.0;
386 subcellParametrization(subcellOrd, 2, 0) = (v0[2] + v1[2])/2.0;
387 subcellParametrization(subcellOrd, 2, 1) = (v1[2] - v0[2])/2.0;
395 else if(subcellDim == 2) {
396 for(
unsigned subcellOrd = 0; subcellOrd < sc; subcellOrd++){
398 switch(parentCell.getKey(subcellDim,subcellOrd)){
401 case shards::Triangle<3>::key:
402 case shards::Triangle<4>::key:
403 case shards::Triangle<6>::key:
405 int v0ord = parentCell.getNodeMap(subcellDim, subcellOrd, 0);
406 int v1ord = parentCell.getNodeMap(subcellDim, subcellOrd, 1);
407 int v2ord = parentCell.getNodeMap(subcellDim, subcellOrd, 2);
408 const double* v0 = getReferenceVertex(parentCell, v0ord);
409 const double* v1 = getReferenceVertex(parentCell, v1ord);
410 const double* v2 = getReferenceVertex(parentCell, v2ord);
413 subcellParametrization(subcellOrd, 0, 0) = v0[0];
414 subcellParametrization(subcellOrd, 0, 1) = v1[0] - v0[0];
415 subcellParametrization(subcellOrd, 0, 2) = v2[0] - v0[0];
418 subcellParametrization(subcellOrd, 1, 0) = v0[1];
419 subcellParametrization(subcellOrd, 1, 1) = v1[1] - v0[1];
420 subcellParametrization(subcellOrd, 1, 2) = v2[1] - v0[1];
423 subcellParametrization(subcellOrd, 2, 0) = v0[2];
424 subcellParametrization(subcellOrd, 2, 1) = v1[2] - v0[2];
425 subcellParametrization(subcellOrd, 2, 2) = v2[2] - v0[2];
430 case shards::Quadrilateral<4>::key:
431 case shards::Quadrilateral<8>::key:
432 case shards::Quadrilateral<9>::key:
434 int v0ord = parentCell.getNodeMap(subcellDim, subcellOrd, 0);
435 int v1ord = parentCell.getNodeMap(subcellDim, subcellOrd, 1);
436 int v2ord = parentCell.getNodeMap(subcellDim, subcellOrd, 2);
437 int v3ord = parentCell.getNodeMap(subcellDim, subcellOrd, 3);
438 const double* v0 = getReferenceVertex(parentCell, v0ord);
439 const double* v1 = getReferenceVertex(parentCell, v1ord);
440 const double* v2 = getReferenceVertex(parentCell, v2ord);
441 const double* v3 = getReferenceVertex(parentCell, v3ord);
444 subcellParametrization(subcellOrd, 0, 0) = ( v0[0] + v1[0] + v2[0] + v3[0])/4.0;
445 subcellParametrization(subcellOrd, 0, 1) = (-v0[0] + v1[0] + v2[0] - v3[0])/4.0;
446 subcellParametrization(subcellOrd, 0, 2) = (-v0[0] - v1[0] + v2[0] + v3[0])/4.0;
449 subcellParametrization(subcellOrd, 1, 0) = ( v0[1] + v1[1] + v2[1] + v3[1])/4.0;
450 subcellParametrization(subcellOrd, 1, 1) = (-v0[1] + v1[1] + v2[1] - v3[1])/4.0;
451 subcellParametrization(subcellOrd, 1, 2) = (-v0[1] - v1[1] + v2[1] + v3[1])/4.0;
454 subcellParametrization(subcellOrd, 2, 0) = ( v0[2] + v1[2] + v2[2] + v3[2])/4.0;
455 subcellParametrization(subcellOrd, 2, 1) = (-v0[2] + v1[2] + v2[2] - v3[2])/4.0;
456 subcellParametrization(subcellOrd, 2, 2) = (-v0[2] - v1[2] + v2[2] + v3[2])/4.0;
460 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
461 ">>> ERROR (Intrepid::CellTools::setSubcellParametrization): parametrization not defined for the specified face topology.");
470 template<
class Scalar>
472 const int vertexOrd){
474 #ifdef HAVE_INTREPID_DEBUG
475 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(cell) ), std::invalid_argument,
476 ">>> ERROR (Intrepid::CellTools::getReferenceVertex): the specified cell topology does not have a reference cell.");
478 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= vertexOrd) && (vertexOrd < (
int)cell.getVertexCount() ) ), std::invalid_argument,
479 ">>> ERROR (Intrepid::CellTools::getReferenceVertex): invalid node ordinal for the specified cell topology. ");
483 return getReferenceNode(cell.getBaseCellTopologyData(), vertexOrd);
488 template<
class Scalar>
489 template<
class ArraySubcellVert>
491 const int subcellDim,
492 const int subcellOrd,
493 const shards::CellTopology& parentCell){
494 #ifdef HAVE_INTREPID_DEBUG
495 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(parentCell) ), std::invalid_argument,
496 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellVertices): the specified cell topology does not have a reference cell.");
500 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= subcellDim) && (subcellDim <= (
int)parentCell.getDimension()) ), std::invalid_argument,
501 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellVertices): subcell dimension out of range.");
503 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= subcellOrd) && (subcellOrd < (
int)parentCell.getSubcellCount(subcellDim) ) ), std::invalid_argument,
504 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellVertices): subcell ordinal out of range.");
508 std::string errmsg =
">>> ERROR (Intrepid::CellTools::getReferenceSubcellVertices):";
509 TEUCHOS_TEST_FOR_EXCEPTION( !( requireRankRange(errmsg, subcellVertices, 2, 2) ), std::invalid_argument, errmsg);
511 int subcVertexCount = parentCell.getVertexCount(subcellDim, subcellOrd);
512 int spaceDim = parentCell.getDimension();
514 TEUCHOS_TEST_FOR_EXCEPTION( !( requireDimensionRange(errmsg, subcellVertices, 0, subcVertexCount, subcVertexCount) ),
515 std::invalid_argument, errmsg);
517 TEUCHOS_TEST_FOR_EXCEPTION( !( requireDimensionRange(errmsg, subcellVertices, 1, spaceDim, spaceDim) ),
518 std::invalid_argument, errmsg);
523 getReferenceSubcellNodes(subcellVertices, subcellDim, subcellOrd, parentCell.getBaseCellTopologyData() );
528 template<
class Scalar>
532 #ifdef HAVE_INTREPID_DEBUG
533 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(cell) ), std::invalid_argument,
534 ">>> ERROR (Intrepid::CellTools::getReferenceNode): the specified cell topology does not have a reference cell.");
536 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= nodeOrd) && (nodeOrd < (
int)cell.getNodeCount() ) ), std::invalid_argument,
537 ">>> ERROR (Intrepid::CellTools::getReferenceNode): invalid node ordinal for the specified cell topology. ");
542 static const double line[2][3] ={
543 {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}
545 static const double line_3[3][3] = {
546 {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0},
553 static const double triangle[3][3] = {
554 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}
556 static const double triangle_4[4][3] = {
557 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
561 static const double triangle_6[6][3] = {
562 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
564 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}
569 static const double quadrilateral[4][3] = {
570 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}
572 static const double quadrilateral_8[8][3] = {
573 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
575 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0}
577 static const double quadrilateral_9[9][3] = {
578 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
580 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0}, { 0.0, 0.0, 0.0}
585 static const double tetrahedron[4][3] = {
586 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0}
588 static const double tetrahedron_8[8][3] = {
589 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
591 { 1/3, 0.0, 1/3}, { 1/3, 1/3, 1/3}, { 1/3, 1/3, 0.0}, { 0.0, 1/3, 1/3}
593 static const double tetrahedron_10[10][3] = {
594 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
596 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}, { 0.0, 0.0, 0.5}, { 0.5, 0.0, 0.5}, { 0.0, 0.5, 0.5}
599 static const double tetrahedron_11[10][3] = {
600 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
602 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}, { 0.0, 0.0, 0.5}, { 0.5, 0.0, 0.5}, { 0.0, 0.5, 0.5}
607 static const double hexahedron[8][3] = {
608 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
609 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0}
611 static const double hexahedron_20[20][3] = {
612 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
613 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0},
615 { 0.0,-1.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, {-1.0, 0.0,-1.0},
616 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
617 { 0.0,-1.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}, {-1.0, 0.0, 1.0}
619 static const double hexahedron_27[27][3] = {
620 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
621 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0},
623 { 0.0,-1.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, {-1.0, 0.0,-1.0},
624 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
625 { 0.0,-1.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}, {-1.0, 0.0, 1.0},
627 { 0.0, 0.0,-1.0}, { 0.0, 0.0, 1.0}, {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, {0.0,-1.0, 0.0}, {0.0, 1.0, 0.0}
632 static const double pyramid[5][3] = {
633 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0}
635 static const double pyramid_13[13][3] = {
636 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
638 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0},
639 {-0.5,-0.5, 0.5}, { 0.5,-0.5, 0.5}, { 0.5, 0.5, 0.5}, {-0.5, 0.5, 0.5}
641 static const double pyramid_14[14][3] = {
642 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
644 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0},
645 {-0.5,-0.5, 0.5}, { 0.5,-0.5, 0.5}, { 0.5, 0.5, 0.5}, {-0.5, 0.5, 0.5}, { 0.0, 0.0, 0.0}
650 static const double wedge[6][3] = {
651 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}
653 static const double wedge_15[15][3] = {
654 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0},
656 { 0.5, 0.0,-1.0}, { 0.5, 0.5,-1.0}, { 0.0, 0.5,-1.0}, { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
657 { 0.5, 0.0, 1.0}, { 0.5, 0.5, 1.0}, { 0.0, 0.5, 1.0}
659 static const double wedge_18[18][3] = {
660 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0},
662 { 0.5, 0.0,-1.0}, { 0.5, 0.5,-1.0}, { 0.0, 0.5,-1.0}, { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
663 { 0.5, 0.0, 1.0}, { 0.5, 0.5, 1.0}, { 0.0, 0.5, 1.0},
664 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}
668 switch(cell.getKey() ) {
671 case shards::Line<2>::key:
672 case shards::ShellLine<2>::key:
673 case shards::Beam<2>::key:
674 return line[nodeOrd];
678 case shards::Line<3>::key:
679 case shards::ShellLine<3>::key:
680 case shards::Beam<3>::key:
681 return line_3[nodeOrd];
686 case shards::Triangle<3>::key:
687 case shards::ShellTriangle<3>::key:
688 return triangle[nodeOrd];
692 case shards::Triangle<4>::key:
693 return triangle_4[nodeOrd];
695 case shards::Triangle<6>::key:
696 case shards::ShellTriangle<6>::key:
697 return triangle_6[nodeOrd];
702 case shards::Quadrilateral<4>::key:
703 case shards::ShellQuadrilateral<4>::key:
704 return quadrilateral[nodeOrd];
708 case shards::Quadrilateral<8>::key:
709 case shards::ShellQuadrilateral<8>::key:
710 return quadrilateral_8[nodeOrd];
712 case shards::Quadrilateral<9>::key:
713 case shards::ShellQuadrilateral<9>::key:
714 return quadrilateral_9[nodeOrd];
719 case shards::Tetrahedron<4>::key:
720 return tetrahedron[nodeOrd];
724 case shards::Tetrahedron<8>::key:
725 return tetrahedron_8[nodeOrd];
727 case shards::Tetrahedron<10>::key:
728 return tetrahedron_10[nodeOrd];
730 case shards::Tetrahedron<11>::key:
731 return tetrahedron_11[nodeOrd];
736 case shards::Hexahedron<8>::key:
737 return hexahedron[nodeOrd];
741 case shards::Hexahedron<20>::key:
742 return hexahedron_20[nodeOrd];
744 case shards::Hexahedron<27>::key:
745 return hexahedron_27[nodeOrd];
750 case shards::Pyramid<5>::key:
751 return pyramid[nodeOrd];
755 case shards::Pyramid<13>::key:
756 return pyramid_13[nodeOrd];
758 case shards::Pyramid<14>::key:
759 return pyramid_14[nodeOrd];
764 case shards::Wedge<6>::key:
765 return wedge[nodeOrd];
769 case shards::Wedge<15>::key:
770 return wedge_15[nodeOrd];
772 case shards::Wedge<18>::key:
773 return wedge_18[nodeOrd];
777 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
778 ">>> ERROR (Intrepid::CellTools::getReferenceNode): invalid cell topology.");
786 template<
class Scalar>
787 template<
class ArraySubcellNode>
789 const int subcellDim,
790 const int subcellOrd,
791 const shards::CellTopology& parentCell){
792 #ifdef HAVE_INTREPID_DEBUG
793 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(parentCell) ), std::invalid_argument,
794 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellNodes): the specified cell topology does not have a reference cell.");
798 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= subcellDim) && (subcellDim <= (
int)parentCell.getDimension()) ), std::invalid_argument,
799 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellNodes): subcell dimension out of range.");
801 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= subcellOrd) && (subcellOrd < (
int)parentCell.getSubcellCount(subcellDim) ) ), std::invalid_argument,
802 ">>> ERROR (Intrepid::CellTools::getReferenceSubcellNodes): subcell ordinal out of range.");
806 std::string errmsg =
">>> ERROR (Intrepid::CellTools::getReferenceSubcellNodes):";
807 TEUCHOS_TEST_FOR_EXCEPTION( !( requireRankRange(errmsg, subcellNodes, 2, 2) ), std::invalid_argument, errmsg);
809 int subcNodeCount = parentCell.getNodeCount(subcellDim, subcellOrd);
810 int spaceDim = parentCell.getDimension();
812 TEUCHOS_TEST_FOR_EXCEPTION( !( requireDimensionRange(errmsg, subcellNodes, 0, subcNodeCount, subcNodeCount) ),
813 std::invalid_argument, errmsg);
815 TEUCHOS_TEST_FOR_EXCEPTION( !( requireDimensionRange(errmsg, subcellNodes, 1, spaceDim, spaceDim) ),
816 std::invalid_argument, errmsg);
821 int subcNodeCount = parentCell.getNodeCount(subcellDim, subcellOrd);
824 for(
int subcNodeOrd = 0; subcNodeOrd < subcNodeCount; subcNodeOrd++){
827 int cellNodeOrd = parentCell.getNodeMap(subcellDim, subcellOrd, subcNodeOrd);
830 for(
int dim = 0; dim < (int)parentCell.getDimension(); dim++){
838 template<
class Scalar>
841 switch(cell.getKey() ) {
842 case shards::Line<2>::key:
843 case shards::Line<3>::key:
844 case shards::ShellLine<2>::key:
845 case shards::ShellLine<3>::key:
846 case shards::Beam<2>::key:
847 case shards::Beam<3>::key:
849 case shards::Triangle<3>::key:
850 case shards::Triangle<4>::key:
851 case shards::Triangle<6>::key:
852 case shards::ShellTriangle<3>::key:
853 case shards::ShellTriangle<6>::key:
855 case shards::Quadrilateral<4>::key:
856 case shards::Quadrilateral<8>::key:
857 case shards::Quadrilateral<9>::key:
858 case shards::ShellQuadrilateral<4>::key:
859 case shards::ShellQuadrilateral<8>::key:
860 case shards::ShellQuadrilateral<9>::key:
862 case shards::Tetrahedron<4>::key:
863 case shards::Tetrahedron<8>::key:
864 case shards::Tetrahedron<10>::key:
865 case shards::Tetrahedron<11>::key:
867 case shards::Hexahedron<8>::key:
868 case shards::Hexahedron<20>::key:
869 case shards::Hexahedron<27>::key:
871 case shards::Pyramid<5>::key:
872 case shards::Pyramid<13>::key:
873 case shards::Pyramid<14>::key:
875 case shards::Wedge<6>::key:
876 case shards::Wedge<15>::key:
877 case shards::Wedge<18>::key:
896 template<
class Scalar>
897 template<
class ArrayJac,
class ArrayPo
int,
class ArrayCell>
899 const ArrayPoint & points,
900 const ArrayCell & cellWorkset,
901 const shards::CellTopology & cellTopo,
902 const int & whichCell)
904 INTREPID_VALIDATE( validateArguments_setJacobian(jacobian, points, cellWorkset, whichCell, cellTopo) );
909 int spaceDim = (size_t)cellTopo.getDimension();
910 size_t numCells =
static_cast<size_t>(cellWorkset.dimension(0));
912 size_t numPoints = (getrank(points) == 2) ? static_cast<size_t>(points.dimension(0)) : static_cast<size_t>(points.dimension(1));
915 Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis;
918 switch( cellTopo.getKey() ){
921 case shards::Line<2>::key:
925 case shards::Triangle<3>::key:
929 case shards::Quadrilateral<4>::key:
933 case shards::Tetrahedron<4>::key:
937 case shards::Hexahedron<8>::key:
941 case shards::Wedge<6>::key:
945 case shards::Pyramid<5>::key:
950 case shards::Triangle<6>::key:
953 case shards::Quadrilateral<9>::key:
957 case shards::Tetrahedron<10>::key:
961 case shards::Tetrahedron<11>::key:
965 case shards::Hexahedron<20>::key:
969 case shards::Hexahedron<27>::key:
973 case shards::Wedge<15>::key:
977 case shards::Wedge<18>::key:
981 case shards::Pyramid<13>::key:
986 case shards::Quadrilateral<8>::key:
987 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
988 ">>> ERROR (Intrepid::CellTools::setJacobian): Cell topology not supported. ");
992 case shards::Line<3>::key:
993 case shards::Beam<2>::key:
994 case shards::Beam<3>::key:
995 case shards::ShellLine<2>::key:
996 case shards::ShellLine<3>::key:
997 case shards::ShellTriangle<3>::key:
998 case shards::ShellTriangle<6>::key:
999 case shards::ShellQuadrilateral<4>::key:
1000 case shards::ShellQuadrilateral<8>::key:
1001 case shards::ShellQuadrilateral<9>::key:
1002 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1003 ">>> ERROR (Intrepid::CellTools::setJacobian): Cell topology not supported. ");
1006 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1007 ">>> ERROR (Intrepid::CellTools::setJacobian): Cell topology not supported.");
1011 int basisCardinality = HGRAD_Basis -> getCardinality();
1012 FieldContainer<Scalar> basisGrads(basisCardinality, numPoints, spaceDim);
1015 if(getrank(jacobian)==4){
1016 for (
size_t i=0; i< static_cast<size_t>(jacobian.dimension(0)); i++){
1017 for (
size_t j=0; j< static_cast<size_t>(jacobian.dimension(1)); j++){
1018 for (
size_t k=0; k< static_cast<size_t>(jacobian.dimension(2)); k++){
1019 for (
size_t l=0; l< static_cast<size_t>(jacobian.dimension(3)); l++){
1020 jacobianWrap(i,j,k,l)=0.0;
1027 if(getrank(jacobian)==3){
1028 for (
size_t i=0; i< static_cast<size_t>(jacobian.dimension(0)); i++){
1029 for (
size_t j=0; j< static_cast<size_t>(jacobian.dimension(1)); j++){
1030 for (
size_t k=0; k< static_cast<size_t>(jacobian.dimension(2)); k++){
1031 jacobianWrap(i,j,k)=0.0;
1037 switch(getrank(points)) {
1043 FieldContainer<Scalar> tempPoints( static_cast<size_t>(points.dimension(0)), static_cast<size_t>(points.dimension(1)) );
1045 for(
size_t pt = 0; pt < static_cast<size_t>(points.dimension(0)); pt++){
1046 for(
size_t dm = 0; dm < static_cast<size_t>(points.dimension(1)) ; dm++){
1047 tempPoints(pt, dm) = pointsWrap(pt, dm);
1051 HGRAD_Basis -> getValues(basisGrads, tempPoints, OPERATOR_GRAD);
1055 size_t cellLoop = (whichCell == -1) ? numCells : 1 ;
1057 if(whichCell == -1) {
1058 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1059 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1060 for(
int row = 0; row < spaceDim; row++){
1061 for(
int col = 0; col < spaceDim; col++){
1064 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1065 jacobianWrap(cellOrd, pointOrd, row, col) += cellWorksetWrap(cellOrd, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1076 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1077 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1078 for(
int row = 0; row < spaceDim; row++){
1079 for(
int col = 0; col < spaceDim; col++){
1082 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1083 jacobianWrap(pointOrd, row, col) += cellWorksetWrap(whichCell, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1098 FieldContainer<Scalar> tempPoints( static_cast<size_t>(points.dimension(1)), static_cast<size_t>(points.dimension(2)) );
1099 for(
size_t cellOrd = 0; cellOrd < numCells; cellOrd++) {
1102 for(
size_t pt = 0; pt < static_cast<size_t>(points.dimension(1)); pt++){
1103 for(
size_t dm = 0; dm < static_cast<size_t>(points.dimension(2)) ; dm++){
1104 tempPoints(pt, dm) = pointsWrap(cellOrd, pt, dm);
1109 HGRAD_Basis -> getValues(basisGrads, tempPoints, OPERATOR_GRAD);
1112 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1113 for(
int row = 0; row < spaceDim; row++){
1114 for(
int col = 0; col < spaceDim; col++){
1117 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1118 jacobianWrap(cellOrd, pointOrd, row, col) += cellWorksetWrap(cellOrd, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1130 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(points) == 2) && (getrank(points) == 3) ), std::invalid_argument,
1131 ">>> ERROR (Intrepid::CellTools::setJacobian): rank 2 or 3 required for points array. ");
1136 template<
class Scalar>
1137 template<
class ArrayJac,
class ArrayPo
int,
class ArrayCell>
1138 void CellTools<Scalar>::setJacobian(ArrayJac & jacobian,
1139 const ArrayPoint & points,
1140 const ArrayCell & cellWorkset,
1141 const Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis,
1142 const int & whichCell)
1150 int spaceDim = (size_t)HGRAD_Basis->getBaseCellTopology().getDimension();
1151 size_t numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1153 size_t numPoints = (getrank(points) == 2) ? static_cast<size_t>(points.dimension(0)) : static_cast<size_t>(points.dimension(1));
1156 int basisCardinality = HGRAD_Basis -> getCardinality();
1157 FieldContainer<Scalar> basisGrads(basisCardinality, numPoints, spaceDim);
1160 if(getrank(jacobian)==4){
1161 for (
size_t i=0; i< static_cast<size_t>(jacobian.dimension(0)); i++){
1162 for (
size_t j=0; j< static_cast<size_t>(jacobian.dimension(1)); j++){
1163 for (
size_t k=0; k< static_cast<size_t>(jacobian.dimension(2)); k++){
1164 for (
size_t l=0; l< static_cast<size_t>(jacobian.dimension(3)); l++){
1165 jacobianWrap(i,j,k,l)=0.0;
1172 if(getrank(jacobian)==3){
1173 for (
size_t i=0; i< static_cast<size_t>(jacobian.dimension(0)); i++){
1174 for (
size_t j=0; j< static_cast<size_t>(jacobian.dimension(1)); j++){
1175 for (
size_t k=0; k< static_cast<size_t>(jacobian.dimension(2)); k++){
1176 jacobianWrap(i,j,k)=0.0;
1182 switch(getrank(points)) {
1188 FieldContainer<Scalar> tempPoints( static_cast<size_t>(points.dimension(0)), static_cast<size_t>(points.dimension(1)) );
1190 for(
size_t pt = 0; pt < static_cast<size_t>(points.dimension(0)); pt++){
1191 for(
size_t dm = 0; dm < static_cast<size_t>(points.dimension(1)) ; dm++){
1192 tempPoints(pt, dm) = pointsWrap(pt, dm);
1196 HGRAD_Basis -> getValues(basisGrads, tempPoints, OPERATOR_GRAD);
1200 size_t cellLoop = (whichCell == -1) ? numCells : 1 ;
1202 if(whichCell == -1) {
1203 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1204 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1205 for(
int row = 0; row < spaceDim; row++){
1206 for(
int col = 0; col < spaceDim; col++){
1209 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1210 jacobianWrap(cellOrd, pointOrd, row, col) += cellWorksetWrap(cellOrd, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1221 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1222 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1223 for(
int row = 0; row < spaceDim; row++){
1224 for(
int col = 0; col < spaceDim; col++){
1227 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1228 jacobianWrap(pointOrd, row, col) += cellWorksetWrap(whichCell, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1243 FieldContainer<Scalar> tempPoints( static_cast<size_t>(points.dimension(1)), static_cast<size_t>(points.dimension(2)) );
1244 for(
size_t cellOrd = 0; cellOrd < numCells; cellOrd++) {
1247 for(
size_t pt = 0; pt < static_cast<size_t>(points.dimension(1)); pt++){
1248 for(
size_t dm = 0; dm < static_cast<size_t>(points.dimension(2)) ; dm++){
1249 tempPoints(pt, dm) = pointsWrap(cellOrd, pt, dm);
1254 HGRAD_Basis -> getValues(basisGrads, tempPoints, OPERATOR_GRAD);
1257 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1258 for(
int row = 0; row < spaceDim; row++){
1259 for(
int col = 0; col < spaceDim; col++){
1262 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1263 jacobianWrap(cellOrd, pointOrd, row, col) += cellWorksetWrap(cellOrd, bfOrd, row)*basisGrads(bfOrd, pointOrd, col);
1275 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(points) == 2) && (getrank(points) == 3) ), std::invalid_argument,
1276 ">>> ERROR (Intrepid::CellTools::setJacobian): rank 2 or 3 required for points array. ");
1281 template<
class Scalar>
1282 template<
class ArrayJacInv,
class ArrayJac>
1284 const ArrayJac & jacobian)
1286 INTREPID_VALIDATE( validateArguments_setJacobianInv(jacobianInv, jacobian) );
1312 template<
class Scalar>
1313 template<
class ArrayJacDet,
class ArrayJac>
1315 const ArrayJac & jacobian)
1317 INTREPID_VALIDATE( validateArguments_setJacobianDetArgs(jacobianDet, jacobian) );
1327 template<
class Scalar>
1328 template<
class ArrayPhysPo
int,
class ArrayRefPo
int,
class ArrayCell>
1330 const ArrayRefPoint & refPoints,
1331 const ArrayCell & cellWorkset,
1333 const int & whichCell)
1341 size_t spaceDim = (size_t)HGRAD_Basis->getBaseCellTopology().getDimension();
1343 size_t numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1345 size_t numPoints = (getrank(refPoints) == 2) ? static_cast<size_t>(refPoints.dimension(0)) : static_cast<size_t>(refPoints.dimension(1));
1348 int basisCardinality = HGRAD_Basis -> getCardinality();
1352 if(getrank(physPoints)==3){
1353 for(
size_t i = 0; i < static_cast<size_t>(physPoints.dimension(0)); i++) {
1354 for(
size_t j = 0; j < static_cast<size_t>(physPoints.dimension(1)); j++){
1355 for(
size_t k = 0; k < static_cast<size_t>(physPoints.dimension(2)); k++){
1356 physPointsWrap(i,j,k) = 0.0;
1360 }
else if(getrank(physPoints)==2){
1361 for(
size_t i = 0; i < static_cast<size_t>(physPoints.dimension(0)); i++){
1362 for(
size_t j = 0; j < static_cast<size_t>(physPoints.dimension(1)); j++){
1363 physPointsWrap(i,j) = 0.0;
1370 switch(getrank(refPoints)) {
1377 FieldContainer<Scalar> tempPoints( static_cast<size_t>(refPoints.dimension(0)), static_cast<size_t>(refPoints.dimension(1)) );
1379 for(
size_t pt = 0; pt < static_cast<size_t>(refPoints.dimension(0)); pt++){
1380 for(
size_t dm = 0; dm < static_cast<size_t>(refPoints.dimension(1)) ; dm++){
1381 tempPoints(pt, dm) = refPointsWrap(pt, dm);
1384 HGRAD_Basis -> getValues(basisVals, tempPoints, OPERATOR_VALUE);
1387 size_t cellLoop = (whichCell == -1) ? numCells : 1 ;
1390 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1391 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1392 for(
size_t dim = 0; dim < spaceDim; dim++){
1393 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1395 if(whichCell == -1){
1396 physPointsWrap(cellOrd, pointOrd, dim) += cellWorksetWrap(cellOrd, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1399 physPointsWrap(pointOrd, dim) += cellWorksetWrap(whichCell, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1414 FieldContainer<Scalar> tempPoints( static_cast<size_t>(refPoints.dimension(1)), static_cast<size_t>(refPoints.dimension(2)) );
1417 for(
size_t cellOrd = 0; cellOrd < numCells; cellOrd++) {
1420 for(
size_t pt = 0; pt < static_cast<size_t>(refPoints.dimension(1)); pt++){
1421 for(
size_t dm = 0; dm < static_cast<size_t>(refPoints.dimension(2)) ; dm++){
1422 tempPoints(pt, dm) = refPointsWrap(cellOrd, pt, dm);
1427 HGRAD_Basis -> getValues(basisVals, tempPoints, OPERATOR_VALUE);
1429 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1430 for(
size_t dim = 0; dim < spaceDim; dim++){
1431 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1433 physPointsWrap(cellOrd, pointOrd, dim) += cellWorksetWrap(cellOrd, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1445 template<
class Scalar>
1446 template<
class ArrayPhysPo
int,
class ArrayRefPo
int,
class ArrayCell>
1448 const ArrayRefPoint & refPoints,
1449 const ArrayCell & cellWorkset,
1450 const shards::CellTopology & cellTopo,
1451 const int & whichCell)
1453 INTREPID_VALIDATE(validateArguments_mapToPhysicalFrame( physPoints, refPoints, cellWorkset, cellTopo, whichCell) );
1459 size_t spaceDim = (size_t)cellTopo.getDimension();
1460 size_t numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1462 size_t numPoints = (getrank(refPoints) == 2) ? static_cast<size_t>(refPoints.dimension(0)) : static_cast<size_t>(refPoints.dimension(1));
1465 Teuchos::RCP<Basis<Scalar, FieldContainer<Scalar> > > HGRAD_Basis;
1468 switch( cellTopo.getKey() ){
1471 case shards::Line<2>::key:
1475 case shards::Triangle<3>::key:
1479 case shards::Quadrilateral<4>::key:
1483 case shards::Tetrahedron<4>::key:
1487 case shards::Hexahedron<8>::key:
1491 case shards::Wedge<6>::key:
1495 case shards::Pyramid<5>::key:
1500 case shards::Triangle<6>::key:
1504 case shards::Quadrilateral<9>::key:
1508 case shards::Tetrahedron<10>::key:
1512 case shards::Tetrahedron<11>::key:
1516 case shards::Hexahedron<20>::key:
1520 case shards::Hexahedron<27>::key:
1524 case shards::Wedge<15>::key:
1528 case shards::Wedge<18>::key:
1532 case shards::Pyramid<13>::key:
1537 case shards::Quadrilateral<8>::key:
1538 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1539 ">>> ERROR (Intrepid::CellTools::mapToPhysicalFrame): Cell topology not supported. ");
1543 case shards::Line<3>::key:
1544 case shards::Beam<2>::key:
1545 case shards::Beam<3>::key:
1546 case shards::ShellLine<2>::key:
1547 case shards::ShellLine<3>::key:
1548 case shards::ShellTriangle<3>::key:
1549 case shards::ShellTriangle<6>::key:
1550 case shards::ShellQuadrilateral<4>::key:
1551 case shards::ShellQuadrilateral<8>::key:
1552 case shards::ShellQuadrilateral<9>::key:
1553 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1554 ">>> ERROR (Intrepid::CellTools::mapToPhysicalFrame): Cell topology not supported. ");
1557 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1558 ">>> ERROR (Intrepid::CellTools::mapToPhysicalFrame): Cell topology not supported.");
1562 int basisCardinality = HGRAD_Basis -> getCardinality();
1566 if(getrank(physPoints)==3){
1567 for(
size_t i = 0; i < static_cast<size_t>(physPoints.dimension(0)); i++) {
1568 for(
size_t j = 0; j < static_cast<size_t>(physPoints.dimension(1)); j++){
1569 for(
size_t k = 0; k < static_cast<size_t>(physPoints.dimension(2)); k++){
1570 physPointsWrap(i,j,k) = 0.0;
1574 }
else if(getrank(physPoints)==2){
1575 for(
size_t i = 0; i < static_cast<size_t>(physPoints.dimension(0)); i++){
1576 for(
size_t j = 0; j < static_cast<size_t>(physPoints.dimension(1)); j++){
1577 physPointsWrap(i,j) = 0.0;
1584 switch(getrank(refPoints)) {
1591 FieldContainer<Scalar> tempPoints( static_cast<size_t>(refPoints.dimension(0)), static_cast<size_t>(refPoints.dimension(1)) );
1593 for(
size_t pt = 0; pt < static_cast<size_t>(refPoints.dimension(0)); pt++){
1594 for(
size_t dm = 0; dm < static_cast<size_t>(refPoints.dimension(1)) ; dm++){
1595 tempPoints(pt, dm) = refPointsWrap(pt, dm);
1598 HGRAD_Basis -> getValues(basisVals, tempPoints, OPERATOR_VALUE);
1601 size_t cellLoop = (whichCell == -1) ? numCells : 1 ;
1604 for(
size_t cellOrd = 0; cellOrd < cellLoop; cellOrd++) {
1605 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1606 for(
size_t dim = 0; dim < spaceDim; dim++){
1607 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1609 if(whichCell == -1){
1610 physPointsWrap(cellOrd, pointOrd, dim) += cellWorksetWrap(cellOrd, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1613 physPointsWrap(pointOrd, dim) += cellWorksetWrap(whichCell, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1628 FieldContainer<Scalar> tempPoints( static_cast<size_t>(refPoints.dimension(1)), static_cast<size_t>(refPoints.dimension(2)) );
1631 for(
size_t cellOrd = 0; cellOrd < numCells; cellOrd++) {
1634 for(
size_t pt = 0; pt < static_cast<size_t>(refPoints.dimension(1)); pt++){
1635 for(
size_t dm = 0; dm < static_cast<size_t>(refPoints.dimension(2)) ; dm++){
1636 tempPoints(pt, dm) = refPointsWrap(cellOrd, pt, dm);
1641 HGRAD_Basis -> getValues(basisVals, tempPoints, OPERATOR_VALUE);
1643 for(
size_t pointOrd = 0; pointOrd < numPoints; pointOrd++) {
1644 for(
size_t dim = 0; dim < spaceDim; dim++){
1645 for(
int bfOrd = 0; bfOrd < basisCardinality; bfOrd++){
1647 physPointsWrap(cellOrd, pointOrd, dim) += cellWorksetWrap(cellOrd, bfOrd, dim)*basisVals(bfOrd, pointOrd);
1660 template<
class Scalar>
1661 template<
class ArrayRefPo
int,
class ArrayPhysPo
int,
class ArrayCell>
1663 const ArrayPhysPoint & physPoints,
1664 const ArrayCell & cellWorkset,
1665 const shards::CellTopology & cellTopo,
1666 const int & whichCell)
1668 INTREPID_VALIDATE( validateArguments_mapToReferenceFrame(refPoints, physPoints, cellWorkset, cellTopo, whichCell) );
1670 size_t spaceDim = (size_t)cellTopo.getDimension();
1676 switch( cellTopo.getKey() ){
1678 case shards::Line<2>::key:
1679 cellCenter(0) = 0.0;
break;
1681 case shards::Triangle<3>::key:
1682 case shards::Triangle<6>::key:
1683 cellCenter(0) = 1./3.; cellCenter(1) = 1./3.;
break;
1685 case shards::Quadrilateral<4>::key:
1686 case shards::Quadrilateral<9>::key:
1687 cellCenter(0) = 0.0; cellCenter(1) = 0.0;
break;
1689 case shards::Tetrahedron<4>::key:
1690 case shards::Tetrahedron<10>::key:
1691 case shards::Tetrahedron<11>::key:
1692 cellCenter(0) = 1./6.; cellCenter(1) = 1./6.; cellCenter(2) = 1./6.;
break;
1694 case shards::Hexahedron<8>::key:
1695 case shards::Hexahedron<20>::key:
1696 case shards::Hexahedron<27>::key:
1697 cellCenter(0) = 0.0; cellCenter(1) = 0.0; cellCenter(2) = 0.0;
break;
1699 case shards::Wedge<6>::key:
1700 case shards::Wedge<15>::key:
1701 case shards::Wedge<18>::key:
1702 cellCenter(0) = 1./3.; cellCenter(1) = 1./3.; cellCenter(2) = 0.0;
break;
1704 case shards::Pyramid<5>::key:
1705 case shards::Pyramid<13>::key:
1706 cellCenter(0) = 0.; cellCenter(1) = 0.; cellCenter(2) = 0.25;
break;
1709 case shards::Quadrilateral<8>::key:
1710 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1711 ">>> ERROR (Intrepid::CellTools::mapToReferenceFrame): Cell topology not supported. ");
1715 case shards::Line<3>::key:
1716 case shards::Beam<2>::key:
1717 case shards::Beam<3>::key:
1718 case shards::ShellLine<2>::key:
1719 case shards::ShellLine<3>::key:
1720 case shards::ShellTriangle<3>::key:
1721 case shards::ShellTriangle<6>::key:
1722 case shards::ShellQuadrilateral<4>::key:
1723 case shards::ShellQuadrilateral<8>::key:
1724 case shards::ShellQuadrilateral<9>::key:
1725 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1726 ">>> ERROR (Intrepid::CellTools::mapToReferenceFrame): Cell topology not supported. ");
1729 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
1730 ">>> ERROR (Intrepid::CellTools::mapToReferenceFrame): Cell topology not supported.");
1737 if(whichCell == -1){
1738 numPoints =
static_cast<size_t>(physPoints.dimension(1));
1739 numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1740 initGuess.
resize(numCells, numPoints, spaceDim);
1742 for(
size_t c = 0; c < numCells; c++){
1743 for(
size_t p = 0; p < numPoints; p++){
1744 for(
size_t d = 0; d < spaceDim; d++){
1745 initGuess(c, p, d) = cellCenter(d);
1752 numPoints =
static_cast<size_t>(physPoints.dimension(0));
1753 initGuess.
resize(numPoints, spaceDim);
1755 for(
size_t p = 0; p < numPoints; p++){
1756 for(
size_t d = 0; d < spaceDim; d++){
1757 initGuess(p, d) = cellCenter(d);
1762 mapToReferenceFrameInitGuess(refPoints, initGuess, physPoints, cellWorkset, cellTopo, whichCell);
1767 template<
class Scalar>
1768 template<
class ArrayRefPo
int,
class ArrayInitGuess,
class ArrayPhysPo
int,
class ArrayCell>
1770 const ArrayInitGuess & initGuess,
1771 const ArrayPhysPoint & physPoints,
1772 const ArrayCell & cellWorkset,
1774 const int & whichCell)
1779 size_t spaceDim = (size_t)HGRAD_Basis->getBaseCellTopology().getDimension();
1792 if(whichCell == -1){
1793 numPoints =
static_cast<size_t>(physPoints.dimension(1));
1794 numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1795 xOld.
resize(numCells, numPoints, spaceDim);
1796 xTem.
resize(numCells, numPoints, spaceDim);
1797 jacobian.
resize(numCells,numPoints, spaceDim, spaceDim);
1798 jacobInv.
resize(numCells,numPoints, spaceDim, spaceDim);
1799 error.
resize(numCells,numPoints);
1801 for(
size_t c = 0; c < numCells; c++){
1802 for(
size_t p = 0; p < numPoints; p++){
1803 for(
size_t d = 0; d < spaceDim; d++){
1804 xOld(c, p, d) = initGuessWrap(c, p, d);
1811 numPoints =
static_cast<size_t>(physPoints.dimension(0));
1812 xOld.
resize(numPoints, spaceDim);
1813 xTem.
resize(numPoints, spaceDim);
1814 jacobian.
resize(numPoints, spaceDim, spaceDim);
1815 jacobInv.
resize(numPoints, spaceDim, spaceDim);
1818 for(
size_t c = 0; c < numCells; c++){
1819 for(
size_t p = 0; p < numPoints; p++){
1820 for(
size_t d = 0; d < spaceDim; d++){
1821 xOld(c, p, d) = initGuessWrap(c, p, d);
1832 setJacobian(jacobian, xOld, cellWorkset, HGRAD_Basis, whichCell);
1833 setJacobianInv(jacobInv, jacobian);
1835 mapToPhysicalFrame( xTem, xOld, cellWorkset, HGRAD_Basis->getBaseCellTopology(), whichCell );
1846 if(whichCell == -1) {
1847 FieldContainer<Scalar> cellWiseError(numCells);
1857 totalError = totalError;
1861 if (totalError < INTREPID_TOL) {
1864 else if ( iter > INTREPID_MAX_NEWTON) {
1865 INTREPID_VALIDATE(std::cout <<
" Intrepid::CellTools::mapToReferenceFrameInitGuess failed to converge to desired tolerance within "
1866 << INTREPID_MAX_NEWTON <<
" iterations\n" );
1872 int refPointsRank=getrank(refPoints);
1873 if (refPointsRank==3){
1874 for(
size_t i=0;i<static_cast<size_t>(refPoints.dimension(0));i++){
1875 for(
size_t j=0;j<static_cast<size_t>(refPoints.dimension(1));j++){
1876 for(
size_t k=0;k<static_cast<size_t>(refPoints.dimension(2));k++){
1877 xOld(i,j,k) = refPointsWrap(i,j,k);
1881 }
else if(refPointsRank==2){
1882 for(
size_t i=0;i<static_cast<size_t>(refPoints.dimension(0));i++){
1883 for(
size_t j=0;j<static_cast<size_t>(refPoints.dimension(1));j++){
1884 xOld(i,j) = refPointsWrap(i,j);
1897 template<
class Scalar>
1898 template<
class ArrayRefPo
int,
class ArrayInitGuess,
class ArrayPhysPo
int,
class ArrayCell>
1900 const ArrayInitGuess & initGuess,
1901 const ArrayPhysPoint & physPoints,
1902 const ArrayCell & cellWorkset,
1903 const shards::CellTopology & cellTopo,
1904 const int & whichCell)
1908 INTREPID_VALIDATE( validateArguments_mapToReferenceFrame(refPoints, initGuess, physPoints, cellWorkset, cellTopo, whichCell) );
1909 size_t spaceDim = (size_t)cellTopo.getDimension();
1922 if(whichCell == -1){
1923 numPoints =
static_cast<size_t>(physPoints.dimension(1));
1924 numCells =
static_cast<size_t>(cellWorkset.dimension(0));
1925 xOld.
resize(numCells, numPoints, spaceDim);
1926 xTem.
resize(numCells, numPoints, spaceDim);
1927 jacobian.
resize(numCells,numPoints, spaceDim, spaceDim);
1928 jacobInv.
resize(numCells,numPoints, spaceDim, spaceDim);
1929 error.
resize(numCells,numPoints);
1931 for(
size_t c = 0; c < numCells; c++){
1932 for(
size_t p = 0; p < numPoints; p++){
1933 for(
size_t d = 0; d < spaceDim; d++){
1934 xOld(c, p, d) = initGuessWrap(c, p, d);
1941 numPoints =
static_cast<size_t>(physPoints.dimension(0));
1942 xOld.
resize(numPoints, spaceDim);
1943 xTem.
resize(numPoints, spaceDim);
1944 jacobian.
resize(numPoints, spaceDim, spaceDim);
1945 jacobInv.
resize(numPoints, spaceDim, spaceDim);
1948 for(
size_t p = 0; p < numPoints; p++){
1949 for(
size_t d = 0; d < spaceDim; d++){
1950 xOld(p, d) = initGuessWrap(p, d);
1960 setJacobian(jacobian, xOld, cellWorkset, cellTopo, whichCell);
1961 setJacobianInv(jacobInv, jacobian);
1963 mapToPhysicalFrame( xTem, xOld, cellWorkset, cellTopo, whichCell );
1974 if(whichCell == -1) {
1985 totalError = totalError;
1989 if (totalError < INTREPID_TOL) {
1992 else if ( iter > INTREPID_MAX_NEWTON) {
1993 INTREPID_VALIDATE(std::cout <<
" Intrepid::CellTools::mapToReferenceFrameInitGuess failed to converge to desired tolerance within "
1994 << INTREPID_MAX_NEWTON <<
" iterations\n" );
2000 int refPointsRank=getrank(refPoints);
2001 if (refPointsRank==3){
2002 for(
size_t i=0;i<static_cast<size_t>(refPoints.dimension(0));i++){
2003 for(
size_t j=0;j<static_cast<size_t>(refPoints.dimension(1));j++){
2004 for(
size_t k=0;k<static_cast<size_t>(refPoints.dimension(2));k++){
2005 xOld(i,j,k) = refPointsWrap(i,j,k);
2009 }
else if(refPointsRank==2){
2010 for(
size_t i=0;i<static_cast<size_t>(refPoints.dimension(0));i++){
2011 for(
size_t j=0;j<static_cast<size_t>(refPoints.dimension(1));j++){
2012 xOld(i,j) = refPointsWrap(i,j);
2025 template<
class Scalar>
2026 template<
class ArraySubcellPo
int,
class ArrayParamPo
int>
2028 const ArrayParamPoint & paramPoints,
2029 const int subcellDim,
2030 const int subcellOrd,
2031 const shards::CellTopology & parentCell){
2033 int cellDim = parentCell.getDimension();
2034 size_t numPts =
static_cast<size_t>(paramPoints.dimension(0));
2036 #ifdef HAVE_INTREPID_DEBUG
2037 TEUCHOS_TEST_FOR_EXCEPTION( !(hasReferenceCell(parentCell) ), std::invalid_argument,
2038 ">>> ERROR (Intrepid::CellTools::mapToReferenceSubcell): the specified cell topology does not have a reference cell.");
2040 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <= subcellDim) && (subcellDim <= 2 ) ), std::invalid_argument,
2041 ">>> ERROR (Intrepid::CellTools::mapToReferenceSubcell): method defined only for 1 and 2-dimensional subcells.");
2043 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= subcellOrd) && (subcellOrd < (
int)parentCell.getSubcellCount(subcellDim) ) ), std::invalid_argument,
2044 ">>> ERROR (Intrepid::CellTools::mapToReferenceSubcell): subcell ordinal out of range.");
2047 std::string errmsg =
">>> ERROR (Intrepid::mapToReferenceSubcell):";
2048 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, refSubcellPoints, 2,2), std::invalid_argument, errmsg);
2049 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, refSubcellPoints, 1, cellDim, cellDim), std::invalid_argument, errmsg);
2052 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, paramPoints, 2,2), std::invalid_argument, errmsg);
2053 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, paramPoints, 1, subcellDim, subcellDim), std::invalid_argument, errmsg);
2056 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, refSubcellPoints, 0, paramPoints, 0), std::invalid_argument, errmsg);
2064 if(subcellDim == 2) {
2065 for(
size_t pt = 0; pt < numPts; pt++){
2066 double u = paramPoints(pt,0);
2067 double v = paramPoints(pt,1);
2070 for(
int dim = 0; dim < cellDim; dim++){
2071 refSubcellPoints(pt, dim) = subcellMap(subcellOrd, dim, 0) + \
2072 subcellMap(subcellOrd, dim, 1)*u + \
2073 subcellMap(subcellOrd, dim, 2)*v;
2077 else if(subcellDim == 1) {
2078 for(
size_t pt = 0; pt < numPts; pt++){
2079 for(
int dim = 0; dim < cellDim; dim++) {
2080 refSubcellPoints(pt, dim) = subcellMap(subcellOrd, dim, 0) + subcellMap(subcellOrd, dim, 1)*paramPoints(pt,0);
2085 TEUCHOS_TEST_FOR_EXCEPTION( !( (subcellDim == 1) || (subcellDim == 2) ), std::invalid_argument,
2086 ">>> ERROR (Intrepid::CellTools::mapToReferenceSubcell): method defined only for 1 and 2-subcells");
2092 template<
class Scalar>
2093 template<
class ArrayEdgeTangent>
2095 const int & edgeOrd,
2096 const shards::CellTopology & parentCell){
2098 int spaceDim = parentCell.getDimension();
2100 #ifdef HAVE_INTREPID_DEBUG
2102 TEUCHOS_TEST_FOR_EXCEPTION( !( (spaceDim == 2) || (spaceDim == 3) ), std::invalid_argument,
2103 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): two or three-dimensional parent cell required");
2105 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= edgeOrd) && (edgeOrd < (
int)parentCell.getSubcellCount(1) ) ), std::invalid_argument,
2106 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): edge ordinal out of bounds");
2108 TEUCHOS_TEST_FOR_EXCEPTION( !( refEdgeTangent.size() == spaceDim ), std::invalid_argument,
2109 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): output array size is required to match space dimension");
2117 refEdgeTangent(0) = edgeMap(edgeOrd, 0, 1);
2118 refEdgeTangent(1) = edgeMap(edgeOrd, 1, 1);
2122 refEdgeTangent(2) = edgeMap(edgeOrd, 2, 1);
2128 template<
class Scalar>
2129 template<
class ArrayFaceTangentU,
class ArrayFaceTangentV>
2131 ArrayFaceTangentV & vTan,
2132 const int & faceOrd,
2133 const shards::CellTopology & parentCell){
2135 #ifdef HAVE_INTREPID_DEBUG
2136 int spaceDim = parentCell.getDimension();
2137 TEUCHOS_TEST_FOR_EXCEPTION( !(spaceDim == 3), std::invalid_argument,
2138 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): three-dimensional parent cell required");
2140 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= faceOrd) && (faceOrd < (
int)parentCell.getSubcellCount(2) ) ), std::invalid_argument,
2141 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): face ordinal out of bounds");
2143 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(uTan) == 1) && (getrank(vTan) == 1) ), std::invalid_argument,
2144 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): rank = 1 required for output arrays");
2146 TEUCHOS_TEST_FOR_EXCEPTION( !( uTan.dimension(0) == spaceDim ), std::invalid_argument,
2147 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): dim0 (spatial dim) must match that of parent cell");
2149 TEUCHOS_TEST_FOR_EXCEPTION( !( vTan.dimension(0) == spaceDim ), std::invalid_argument,
2150 ">>> ERROR (Intrepid::CellTools::getReferenceFaceTangents): dim0 (spatial dim) must match that of parent cell");
2161 uTan(0) = faceMap(faceOrd, 0, 1);
2162 uTan(1) = faceMap(faceOrd, 1, 1);
2163 uTan(2) = faceMap(faceOrd, 2, 1);
2166 vTan(0) = faceMap(faceOrd, 0, 2);
2167 vTan(1) = faceMap(faceOrd, 1, 2);
2168 vTan(2) = faceMap(faceOrd, 2, 2);
2173 template<
class Scalar>
2174 template<
class ArrayS
ideNormal>
2176 const int & sideOrd,
2177 const shards::CellTopology & parentCell){
2178 int spaceDim = parentCell.getDimension();
2180 #ifdef HAVE_INTREPID_DEBUG
2182 TEUCHOS_TEST_FOR_EXCEPTION( !( (spaceDim == 2) || (spaceDim == 3) ), std::invalid_argument,
2183 ">>> ERROR (Intrepid::CellTools::getReferenceSideNormal): two or three-dimensional parent cell required");
2186 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= sideOrd) && (sideOrd < (
int)parentCell.getSubcellCount(spaceDim - 1) ) ), std::invalid_argument,
2187 ">>> ERROR (Intrepid::CellTools::getReferenceSideNormal): side ordinal out of bounds");
2192 getReferenceEdgeTangent(refSideNormal, sideOrd, parentCell);
2195 Scalar temp = refSideNormal(0);
2196 refSideNormal(0) = refSideNormal(1);
2197 refSideNormal(1) = -temp;
2201 getReferenceFaceNormal(refSideNormal, sideOrd, parentCell);
2207 template<
class Scalar>
2208 template<
class ArrayFaceNormal>
2210 const int & faceOrd,
2211 const shards::CellTopology & parentCell){
2212 int spaceDim = parentCell.getDimension();
2213 #ifdef HAVE_INTREPID_DEBUG
2215 TEUCHOS_TEST_FOR_EXCEPTION( !(spaceDim == 3), std::invalid_argument,
2216 ">>> ERROR (Intrepid::CellTools::getReferenceFaceNormal): three-dimensional parent cell required");
2218 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= faceOrd) && (faceOrd < (
int)parentCell.getSubcellCount(2) ) ), std::invalid_argument,
2219 ">>> ERROR (Intrepid::CellTools::getReferenceFaceNormal): face ordinal out of bounds");
2221 TEUCHOS_TEST_FOR_EXCEPTION( !( getrank(refFaceNormal) == 1 ), std::invalid_argument,
2222 ">>> ERROR (Intrepid::CellTools::getReferenceFaceNormal): rank = 1 required for output array");
2224 TEUCHOS_TEST_FOR_EXCEPTION( !( static_cast<size_t>(refFaceNormal.dimension(0)) == static_cast<size_t>(spaceDim) ), std::invalid_argument,
2225 ">>> ERROR (Intrepid::CellTools::getReferenceFaceNormal): dim0 (spatial dim) must match that of parent cell");
2231 getReferenceFaceTangents(uTan, vTan, faceOrd, parentCell);
2237 template<
class Scalar>
2238 template<
class ArrayEdgeTangent,
class ArrayJac>
2240 const ArrayJac & worksetJacobians,
2241 const int & worksetEdgeOrd,
2242 const shards::CellTopology & parentCell){
2243 size_t worksetSize =
static_cast<size_t>(worksetJacobians.dimension(0));
2244 size_t edgePtCount =
static_cast<size_t>(worksetJacobians.dimension(1));
2245 int pCellDim = parentCell.getDimension();
2246 #ifdef HAVE_INTREPID_DEBUG
2247 std::string errmsg =
">>> ERROR (Intrepid::CellTools::getPhysicalEdgeTangents):";
2249 TEUCHOS_TEST_FOR_EXCEPTION( !( (pCellDim == 3) || (pCellDim == 2) ), std::invalid_argument,
2250 ">>> ERROR (Intrepid::CellTools::getPhysicalEdgeTangents): 2D or 3D parent cell required");
2253 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, edgeTangents, 3,3), std::invalid_argument, errmsg);
2254 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, edgeTangents, 2, 2,3), std::invalid_argument, errmsg);
2257 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, worksetJacobians, 4,4), std::invalid_argument, errmsg);
2258 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 2, 2,3), std::invalid_argument, errmsg);
2259 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 3, 2,3), std::invalid_argument, errmsg);
2262 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, edgeTangents, 0,1,2,2, worksetJacobians, 0,1,2,3), std::invalid_argument, errmsg);
2268 getReferenceEdgeTangent(refEdgeTan, worksetEdgeOrd, parentCell);
2271 for(
size_t pCell = 0; pCell < worksetSize; pCell++){
2272 for(
size_t pt = 0; pt < edgePtCount; pt++){
2275 for(
int i = 0; i < pCellDim; i++){
2276 edgeTangents(pCell, pt, i) = 0.0;
2277 for(
int j = 0; j < pCellDim; j++){
2278 edgeTangents(pCell, pt, i) += worksetJacobians(pCell, pt, i, j)*refEdgeTan(j);
2284 template<
class Scalar>
2285 template<
class ArrayFaceTangentU,
class ArrayFaceTangentV,
class ArrayJac>
2287 ArrayFaceTangentV & faceTanV,
2288 const ArrayJac & worksetJacobians,
2289 const int & worksetFaceOrd,
2290 const shards::CellTopology & parentCell){
2291 size_t worksetSize =
static_cast<size_t>(worksetJacobians.dimension(0));
2292 size_t facePtCount =
static_cast<size_t>(worksetJacobians.dimension(1));
2293 int pCellDim = parentCell.getDimension();
2294 #ifdef HAVE_INTREPID_DEBUG
2295 std::string errmsg =
">>> ERROR (Intrepid::CellTools::getPhysicalFaceTangents):";
2297 TEUCHOS_TEST_FOR_EXCEPTION( !(pCellDim == 3), std::invalid_argument,
2298 ">>> ERROR (Intrepid::CellTools::getPhysicalFaceTangents): three-dimensional parent cell required");
2301 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, faceTanU, 3,3), std::invalid_argument, errmsg);
2302 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, faceTanU, 2, 3,3), std::invalid_argument, errmsg);
2304 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, faceTanV, 3,3), std::invalid_argument, errmsg);
2305 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, faceTanV, 2, 3,3), std::invalid_argument, errmsg);
2307 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, faceTanU, faceTanV), std::invalid_argument, errmsg);
2310 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, worksetJacobians, 4,4), std::invalid_argument, errmsg);
2311 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 2, 3,3), std::invalid_argument, errmsg);
2312 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 3, 3,3), std::invalid_argument, errmsg);
2315 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, faceTanU, 0,1,2,2, worksetJacobians, 0,1,2,3), std::invalid_argument, errmsg);
2322 getReferenceFaceTangents(refFaceTanU, refFaceTanV, worksetFaceOrd, parentCell);
2325 for(
size_t pCell = 0; pCell < worksetSize; pCell++){
2326 for(
size_t pt = 0; pt < facePtCount; pt++){
2329 for(
int dim = 0; dim < pCellDim; dim++){
2330 faceTanU(pCell, pt, dim) = 0.0;
2331 faceTanV(pCell, pt, dim) = 0.0;
2334 faceTanU(pCell, pt, dim) = \
2335 worksetJacobians(pCell, pt, dim, 0)*refFaceTanU(0) + \
2336 worksetJacobians(pCell, pt, dim, 1)*refFaceTanU(1) + \
2337 worksetJacobians(pCell, pt, dim, 2)*refFaceTanU(2);
2338 faceTanV(pCell, pt, dim) = \
2339 worksetJacobians(pCell, pt, dim, 0)*refFaceTanV(0) + \
2340 worksetJacobians(pCell, pt, dim, 1)*refFaceTanV(1) + \
2341 worksetJacobians(pCell, pt, dim, 2)*refFaceTanV(2);
2347 template<
class Scalar>
2348 template<
class ArrayS
ideNormal,
class ArrayJac>
2350 const ArrayJac & worksetJacobians,
2351 const int & worksetSideOrd,
2352 const shards::CellTopology & parentCell){
2353 size_t worksetSize =
static_cast<size_t>(worksetJacobians.dimension(0));
2354 size_t sidePtCount =
static_cast<size_t>(worksetJacobians.dimension(1));
2355 int spaceDim = parentCell.getDimension();
2356 #ifdef HAVE_INTREPID_DEBUG
2357 TEUCHOS_TEST_FOR_EXCEPTION( !( (spaceDim == 2) || (spaceDim == 3) ), std::invalid_argument,
2358 ">>> ERROR (Intrepid::CellTools::getPhysicalSideNormals): two or three-dimensional parent cell required");
2361 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 <= worksetSideOrd) && (worksetSideOrd < (
int)parentCell.getSubcellCount(spaceDim - 1) ) ), std::invalid_argument,
2362 ">>> ERROR (Intrepid::CellTools::getPhysicalSideNormals): side ordinal out of bounds");
2368 getPhysicalEdgeTangents(sideNormals, worksetJacobians, worksetSideOrd, parentCell);
2371 for(
size_t cell = 0; cell < worksetSize; cell++){
2372 for(
size_t pt = 0; pt < sidePtCount; pt++){
2373 Scalar temp = sideNormals(cell, pt, 0);
2374 sideNormals(cell, pt, 0) = sideNormals(cell, pt, 1);
2375 sideNormals(cell, pt, 1) = -temp;
2381 getPhysicalFaceNormals(sideNormals, worksetJacobians, worksetSideOrd, parentCell);
2386 template<
class Scalar>
2387 template<
class ArrayFaceNormal,
class ArrayJac>
2389 const ArrayJac & worksetJacobians,
2390 const int & worksetFaceOrd,
2391 const shards::CellTopology & parentCell){
2392 size_t worksetSize =
static_cast<size_t>(worksetJacobians.dimension(0));
2393 size_t facePtCount =
static_cast<size_t>(worksetJacobians.dimension(1));
2394 int pCellDim = parentCell.getDimension();
2395 #ifdef HAVE_INTREPID_DEBUG
2396 std::string errmsg =
">>> ERROR (Intrepid::CellTools::getPhysicalFaceNormals):";
2398 TEUCHOS_TEST_FOR_EXCEPTION( !(pCellDim == 3), std::invalid_argument,
2399 ">>> ERROR (Intrepid::CellTools::getPhysicalFaceNormals): three-dimensional parent cell required");
2402 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, faceNormals, 3,3), std::invalid_argument, errmsg);
2403 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, faceNormals, 2, 3,3), std::invalid_argument, errmsg);
2406 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, worksetJacobians, 4,4), std::invalid_argument, errmsg);
2407 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 2, 3,3), std::invalid_argument, errmsg);
2408 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionRange(errmsg, worksetJacobians, 3, 3,3), std::invalid_argument, errmsg);
2411 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, faceNormals, 0,1,2,2, worksetJacobians, 0,1,2,3), std::invalid_argument, errmsg);
2417 getPhysicalFaceTangents(faceTanU, faceTanV, worksetJacobians, worksetFaceOrd, parentCell);
2431 template<
class Scalar>
2434 const shards::CellTopology & cellTopo,
2435 const double & threshold) {
2436 #ifdef HAVE_INTREPID_DEBUG
2437 TEUCHOS_TEST_FOR_EXCEPTION( !(pointDim == (
int)cellTopo.getDimension() ), std::invalid_argument,
2438 ">>> ERROR (Intrepid::CellTools::checkPointInclusion): Point and cell dimensions do not match. ");
2443 double minus_one = -1.0 - threshold;
2444 double plus_one = 1.0 + threshold;
2445 double minus_zero = - threshold;
2450 unsigned key = cellTopo.getBaseCellTopologyData() -> key ;
2453 case shards::Line<>::key :
2454 if( !(minus_one <= point[0] && point[0] <= plus_one)) testResult = 0;
2457 case shards::Triangle<>::key : {
2458 Scalar distance = std::max( std::max( -point[0], -point[1] ), point[0] + point[1] - 1.0 );
2459 if( distance > threshold ) testResult = 0;
2463 case shards::Quadrilateral<>::key :
2464 if(!( (minus_one <= point[0] && point[0] <= plus_one) && \
2465 (minus_one <= point[1] && point[1] <= plus_one) ) ) testResult = 0;
2468 case shards::Tetrahedron<>::key : {
2469 Scalar distance = std::max( std::max(-point[0],-point[1]), \
2470 std::max(-point[2], point[0] + point[1] + point[2] - 1) );
2471 if( distance > threshold ) testResult = 0;
2475 case shards::Hexahedron<>::key :
2476 if(!((minus_one <= point[0] && point[0] <= plus_one) && \
2477 (minus_one <= point[1] && point[1] <= plus_one) && \
2478 (minus_one <= point[2] && point[2] <= plus_one))) \
2484 case shards::Wedge<>::key : {
2485 Scalar distance = std::max( std::max( -point[0], -point[1] ), point[0] + point[1] - 1 );
2486 if( distance > threshold || \
2487 point[2] < minus_one || point[2] > plus_one) \
2495 case shards::Pyramid<>::key : {
2496 Scalar left = minus_one + point[2];
2497 Scalar right = plus_one - point[2];
2498 if(!( (left <= point[0] && point[0] <= right) && \
2499 (left <= point[1] && point[1] <= right) &&
2500 (minus_zero <= point[2] && point[2] <= plus_one) ) ) \
2506 TEUCHOS_TEST_FOR_EXCEPTION( !( (key == shards::Line<>::key ) ||
2507 (key == shards::Triangle<>::key) ||
2508 (key == shards::Quadrilateral<>::key) ||
2509 (key == shards::Tetrahedron<>::key) ||
2510 (key == shards::Hexahedron<>::key) ||
2511 (key == shards::Wedge<>::key) ||
2512 (key == shards::Pyramid<>::key) ),
2513 std::invalid_argument,
2514 ">>> ERROR (Intrepid::CellTools::checkPointInclusion): Invalid cell topology. ");
2521 template<
class Scalar>
2522 template<
class ArrayPo
int>
2524 const shards::CellTopology & cellTopo,
2525 const double & threshold) {
2527 int rank = points.rank();
2529 #ifdef HAVE_INTREPID_DEBUG
2530 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <=getrank(points) ) && (getrank(points) <= 3) ), std::invalid_argument,
2531 ">>> ERROR (Intrepid::CellTools::checkPointsetInclusion): rank-1, 2 or 3 required for input points array. ");
2534 TEUCHOS_TEST_FOR_EXCEPTION( !((
size_t) points.dimension(rank - 1) == (size_t)cellTopo.getDimension() ), std::invalid_argument,
2535 ">>> ERROR (Intrepid::CellTools::checkPointsetInclusion): Point and cell dimensions do not match. ");
2541 case 1: inRefCell.
resize(1);
break;
2542 case 2: inRefCell.
resize( static_cast<size_t>(points.dimension(0)) );
break;
2543 case 3: inRefCell.
resize( static_cast<size_t>(points.dimension(0)), static_cast<size_t>(points.dimension(1)) );
break;
2547 checkPointwiseInclusion(inRefCell, points, cellTopo, threshold);
2551 for(
int i = 0; i < inRefCell.
size(); i++ ){
2552 if (inRefCell[i] == 0) {
2562 template<
class Scalar>
2563 template<
class ArrayIncl,
class ArrayPo
int>
2565 const ArrayPoint & points,
2566 const shards::CellTopology & cellTopo,
2567 const double & threshold) {
2568 int apRank = points.rank();
2570 #ifdef HAVE_INTREPID_DEBUG
2573 std::string errmsg =
">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion):";
2574 if(getrank(points) == 1) {
2575 TEUCHOS_TEST_FOR_EXCEPTION( !(getrank(inRefCell) == 1 ), std::invalid_argument,
2576 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-1 input array requires rank-1 output array.");
2577 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(inRefCell.dimension(0)) == 1), std::invalid_argument,
2578 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-1 input array requires dim0 = 1 for output array.");
2580 else if(getrank(points) == 2){
2581 TEUCHOS_TEST_FOR_EXCEPTION( !(getrank(inRefCell) == 1 ), std::invalid_argument,
2582 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-2 input array requires rank-1 output array.");
2584 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch( errmsg, inRefCell, 0, points, 0), std::invalid_argument, errmsg);
2586 else if (getrank(points) == 3) {
2587 TEUCHOS_TEST_FOR_EXCEPTION( !(getrank(inRefCell) == 2 ), std::invalid_argument,
2588 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-3 input array requires rank-2 output array.");
2590 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch( errmsg, inRefCell, 0,1, points, 0,1), std::invalid_argument, errmsg);
2593 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(points) == 1) || (getrank(points) == 2) || (getrank(points) == 3) ), std::invalid_argument,
2594 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-1, 2 or 3 required for input points array. ");
2598 TEUCHOS_TEST_FOR_EXCEPTION( !((
size_t)points.dimension(apRank - 1) == (size_t)cellTopo.getDimension() ), std::invalid_argument,
2599 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): Point and cell dimensions do not match. ");
2609 pointDim =
static_cast<size_t>(points.dimension(0));
2612 dim1 =
static_cast<size_t>(points.dimension(0));
2613 pointDim =
static_cast<size_t>(points.dimension(1));
2616 dim0 =
static_cast<size_t>(points.dimension(0));
2617 dim1 =
static_cast<size_t>(points.dimension(1));
2618 pointDim =
static_cast<size_t>(points.dimension(2));
2621 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <= getrank(points) ) && (getrank(points) <= 3) ), std::invalid_argument,
2622 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): rank-1, 2 or 3 required for input points array. ");
2634 Scalar point[3] = {0.0, 0.0, 0.0};
2636 for(
int i0 = 0; i0 < dim0; i0++){
2638 inPtr0 = outPtr0*pointDim;
2640 for(
int i1 = 0; i1 < dim1; i1++) {
2641 inPtr1 = inPtr0 + i1*pointDim;
2642 point[0] = points[inPtr1];
2644 point[1] = points[inPtr1 + 1];
2646 point[2] = points[inPtr1 + 2];
2648 TEUCHOS_TEST_FOR_EXCEPTION( !( (1 <= pointDim) && (pointDim <= 3)), std::invalid_argument,
2649 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): Input array specifies invalid point dimension ");
2653 inRefCell[outPtr0 + i1] = checkPointInclusion(point, pointDim, cellTopo, threshold);
2660 template<
class Scalar>
2661 template<
class ArrayIncl,
class ArrayPo
int,
class ArrayCell>
2663 const ArrayPoint & points,
2664 const ArrayCell & cellWorkset,
2665 const shards::CellTopology & cell,
2666 const int & whichCell,
2667 const double & threshold)
2669 INTREPID_VALIDATE( validateArguments_checkPointwiseInclusion(inCell, points, cellWorkset, whichCell, cell) );
2673 unsigned baseKey = cell.getBaseCellTopologyData() -> key;
2677 case shards::Line<>::key :
2678 case shards::Triangle<>::key:
2679 case shards::Quadrilateral<>::key :
2680 case shards::Tetrahedron<>::key :
2681 case shards::Hexahedron<>::key :
2682 case shards::Wedge<>::key :
2683 case shards::Pyramid<>::key :
2687 if(getrank(points) == 2){
2688 refPoints.
resize(static_cast<size_t>(points.dimension(0)), static_cast<size_t>(points.dimension(1)) );
2689 mapToReferenceFrame(refPoints, points, cellWorkset, cell, whichCell);
2690 checkPointwiseInclusion(inCell, refPoints, cell, threshold );
2692 else if(getrank(points) == 3){
2693 refPoints.
resize(static_cast<size_t>(points.dimension(0)), static_cast<size_t>(points.dimension(1)), static_cast<size_t>(points.dimension(2)) );
2694 mapToReferenceFrame(refPoints, points, cellWorkset, cell, whichCell);
2695 checkPointwiseInclusion(inCell, refPoints, cell, threshold );
2700 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
2701 ">>> ERROR (Intrepid::CellTools::checkPointwiseInclusion): cell topology not supported");
2713 template<
class Scalar>
2714 template<
class ArrayJac,
class ArrayPo
int,
class ArrayCell>
2716 const ArrayPoint & points,
2717 const ArrayCell & cellWorkset,
2718 const int & whichCell,
2719 const shards::CellTopology & cellTopo){
2722 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(cellWorkset) != 3), std::invalid_argument,
2723 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): rank = 3 required for cellWorkset array");
2725 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(0)) <= 0), std::invalid_argument,
2726 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of cells) >= 1 required for cellWorkset array");
2728 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(1)) != (
size_t)cellTopo.getSubcellCount(0) ), std::invalid_argument,
2729 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (number of cell nodes) of cellWorkset array does not match cell topology");
2731 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(2)) != (
size_t)cellTopo.getDimension() ), std::invalid_argument,
2732 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 (spatial dimension) of cellWorkset array does not match cell dimension");
2735 TEUCHOS_TEST_FOR_EXCEPTION( !( ( (0 <= whichCell ) && (static_cast<size_t>(whichCell) < static_cast<size_t>(cellWorkset.dimension(0)) ) ) || (whichCell == -1) ), std::invalid_argument,
2736 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): whichCell = -1 or a valid cell ordinal is required.");
2741 if(getrank(points) == 2) {
2742 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(points.dimension(0)) <= 0), std::invalid_argument,
2743 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of points) >= 1 required for points array ");
2745 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(points.dimension(1)) != (
size_t)cellTopo.getDimension() ), std::invalid_argument,
2746 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (spatial dimension) of points array does not match cell dimension");
2749 if(whichCell == -1) {
2750 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(jacobian) != 4), std::invalid_argument,
2751 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): rank = 4 required for jacobian array");
2753 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(0)) != static_cast<size_t>(cellWorkset.dimension(0))), std::invalid_argument,
2754 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of cells) of jacobian array must equal dim 0 of cellWorkset array");
2756 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(1)) != static_cast<size_t>(points.dimension(0))), std::invalid_argument,
2757 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (number of points) of jacobian array must equal dim 0 of points array");
2759 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(2)) != static_cast<size_t>(points.dimension(1))), std::invalid_argument,
2760 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 (spatial dimension) of jacobian array must equal dim 1 of points array");
2762 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobian.dimension(2)) == static_cast<size_t>(jacobian.dimension(3)) ), std::invalid_argument,
2763 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 = dim 3 (same spatial dimensions) required for jacobian array. ");
2765 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 < static_cast<size_t>(jacobian.dimension(3)) ) && (
static_cast<size_t>(jacobian.dimension(3)) < 4) ), std::invalid_argument,
2766 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 and dim 3 (spatial dimensions) must be between 1 and 3. ");
2770 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(jacobian) != 3), std::invalid_argument,
2771 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): rank = 3 required for jacobian array");
2773 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(0)) != static_cast<size_t>(points.dimension(0))), std::invalid_argument,
2774 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of points) of jacobian array must equal dim 0 of points array");
2776 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(1)) != static_cast<size_t>(points.dimension(1))), std::invalid_argument,
2777 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (spatial dimension) of jacobian array must equal dim 1 of points array");
2779 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobian.dimension(1)) == static_cast<size_t>(jacobian.dimension(2)) ), std::invalid_argument,
2780 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 = dim 2 (same spatial dimensions) required for jacobian array. ");
2782 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 < static_cast<size_t>(jacobian.dimension(1)) ) && (
static_cast<size_t>(jacobian.dimension(1)) < 4) ), std::invalid_argument,
2783 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 and dim 2 (spatial dimensions) must be between 1 and 3. ");
2787 else if(getrank(points) ==3){
2788 std::string errmsg =
">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian):";
2789 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(points.dimension(0)) != static_cast<size_t>(cellWorkset.dimension(0)) ), std::invalid_argument,
2790 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of cells) of points array must equal dim 0 of cellWorkset array");
2792 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(points.dimension(1)) <= 0), std::invalid_argument,
2793 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (number of points) >= 1 required for points array ");
2795 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(points.dimension(2)) != (
size_t)cellTopo.getDimension() ), std::invalid_argument,
2796 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 (spatial dimension) of points array does not match cell dimension");
2798 TEUCHOS_TEST_FOR_EXCEPTION( (whichCell != -1), std::invalid_argument,
2799 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): default value whichCell=-1 required for rank-3 input points");
2802 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg, jacobian, 4, 4), std::invalid_argument,errmsg);
2804 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(0)) != static_cast<size_t>(points.dimension(0))), std::invalid_argument,
2805 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 0 (number of cells) of jacobian array must equal dim 0 of points array");
2807 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(1)) != static_cast<size_t>(points.dimension(1))), std::invalid_argument,
2808 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 1 (number of points) of jacobian array must equal dim 1 of points array");
2810 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(jacobian.dimension(2)) != static_cast<size_t>(points.dimension(2))), std::invalid_argument,
2811 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 (spatial dimension) of jacobian array must equal dim 2 of points array");
2813 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobian.dimension(2)) == static_cast<size_t>(jacobian.dimension(3)) ), std::invalid_argument,
2814 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 = dim 3 (same spatial dimensions) required for jacobian array. ");
2816 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 < static_cast<size_t>(jacobian.dimension(3)) ) && (
static_cast<size_t>(jacobian.dimension(3)) < 4) ), std::invalid_argument,
2817 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): dim 2 and dim 3 (spatial dimensions) must be between 1 and 3. ");
2820 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(points) == 2) && (getrank(points) ==3) ), std::invalid_argument,
2821 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobian): rank = 2 or 3 required for points array");
2827 template<
class Scalar>
2828 template<
class ArrayJacInv,
class ArrayJac>
2830 const ArrayJac & jacobian)
2834 int jacobRank = getrank(jacobian);
2835 TEUCHOS_TEST_FOR_EXCEPTION( !( (jacobRank == 4) || (jacobRank == 3) ), std::invalid_argument,
2836 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): rank = 4 or 3 required for jacobian array. ");
2839 TEUCHOS_TEST_FOR_EXCEPTION( !(jacobian.dimension(jacobRank - 1) == jacobian.dimension(jacobRank - 2) ), std::invalid_argument,
2840 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): dim(rank-2) = dim(rank-2) (same spatial dimensions) required for jacobian array. ");
2842 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 < jacobian.dimension(jacobRank - 1) ) && (jacobian.dimension(jacobRank - 1) < 4) ), std::invalid_argument,
2843 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): dim(rank-1) and dim(rank-2) (spatial dimensions) must be between 1 and 3. ");
2846 std::string errmsg =
">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv):";
2848 TEUCHOS_TEST_FOR_EXCEPTION( !(requireRankMatch(errmsg, jacobianInv, jacobian) ), std::invalid_argument, errmsg);
2850 TEUCHOS_TEST_FOR_EXCEPTION( !(requireDimensionMatch(errmsg, jacobianInv, jacobian) ), std::invalid_argument, errmsg);
2855 template<
class Scalar>
2856 template<
class ArrayJacDet,
class ArrayJac>
2858 const ArrayJac & jacobian)
2862 int jacobRank = getrank(jacobian);
2863 TEUCHOS_TEST_FOR_EXCEPTION( !( (jacobRank == 4) || (jacobRank == 3) ), std::invalid_argument,
2864 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): rank = 4 or 3 required for jacobian array. ");
2867 TEUCHOS_TEST_FOR_EXCEPTION( !(jacobian.dimension(jacobRank - 1) == jacobian.dimension(jacobRank - 2) ), std::invalid_argument,
2868 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): dim(rank-2) = dim(rank-2) (same spatial dimensions) required for jacobian array. ");
2870 TEUCHOS_TEST_FOR_EXCEPTION( !( (0 < jacobian.dimension(jacobRank - 1) ) && (jacobian.dimension(jacobRank - 1) < 4) ), std::invalid_argument,
2871 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianInv): dim(rank-1) and dim(rank-2) (spatial dimensions) must be between 1 and 3. ");
2876 TEUCHOS_TEST_FOR_EXCEPTION( !(getrank(jacobianDet) == 2), std::invalid_argument,
2877 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianDetArgs): rank = 2 required for jacobianDet if jacobian is rank-4. ");
2879 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobianDet.dimension(0)) == static_cast<size_t>(jacobian.dimension(0)) ), std::invalid_argument,
2880 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianDetArgs): dim 0 (number of cells) of jacobianDet array must equal dim 0 of jacobian array. ");
2882 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobianDet.dimension(1)) == static_cast<size_t>(jacobian.dimension(1)) ), std::invalid_argument,
2883 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianDetArgs): dim 1 (number of points) of jacobianDet array must equal dim 1 of jacobian array.");
2888 TEUCHOS_TEST_FOR_EXCEPTION( !(getrank(jacobianDet) == 1), std::invalid_argument,
2889 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianDetArgs): rank = 1 required for jacobianDet if jacobian is rank-3. ");
2891 TEUCHOS_TEST_FOR_EXCEPTION( !(static_cast<size_t>(jacobianDet.dimension(0)) == static_cast<size_t>(jacobian.dimension(0)) ), std::invalid_argument,
2892 ">>> ERROR (Intrepid::CellTools::validateArguments_setJacobianDetArgs): dim 0 (number of points) of jacobianDet array must equal dim 0 of jacobian array.");
2898 template<
class Scalar>
2899 template<
class ArrayPhysPo
int,
class ArrayRefPo
int,
class ArrayCell>
2901 const ArrayRefPoint & refPoints,
2902 const ArrayCell & cellWorkset,
2903 const shards::CellTopology & cellTopo,
2904 const int& whichCell)
2906 std::string errmsg =
">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame):";
2909 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(cellWorkset) != 3), std::invalid_argument,
2910 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): rank = 3 required for cellWorkset array");
2912 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(0)) <= 0), std::invalid_argument,
2913 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 0 (number of cells) >= 1 required for cellWorkset array");
2915 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(1)) != (
size_t)cellTopo.getSubcellCount(0) ), std::invalid_argument,
2916 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 1 (number of cell nodes) of cellWorkset array does not match cell topology");
2918 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(2)) != (
size_t)cellTopo.getDimension() ), std::invalid_argument,
2919 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 2 (spatial dimension) of cellWorkset array does not match cell dimension");
2924 TEUCHOS_TEST_FOR_EXCEPTION( !( ( (0 <= whichCell ) && ((
size_t)whichCell < (
size_t)cellWorkset.dimension(0) ) ) || (whichCell == -1) ), std::invalid_argument,
2925 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): whichCell = -1 or a valid cell ordinal is required.");
2929 if(getrank(refPoints) == 2) {
2930 TEUCHOS_TEST_FOR_EXCEPTION( (refPoints.dimension(0) <= 0), std::invalid_argument,
2931 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 0 (number of points) >= 1 required for refPoints array ");
2933 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)refPoints.dimension(1) != (size_t)cellTopo.getDimension() ), std::invalid_argument,
2934 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 1 (spatial dimension) of refPoints array does not match cell dimension");
2937 if(whichCell == -1) {
2938 TEUCHOS_TEST_FOR_EXCEPTION( ( (getrank(physPoints) != 3) && (whichCell == -1) ), std::invalid_argument,
2939 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): rank = 3 required for physPoints array for the default whichCell value");
2941 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)physPoints.dimension(0) != (size_t)cellWorkset.dimension(0)), std::invalid_argument,
2942 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 0 (number of cells) of physPoints array must equal dim 0 of cellWorkset array");
2944 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)physPoints.dimension(1) != (size_t)refPoints.dimension(0)), std::invalid_argument,
2945 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 1 (number of points) of physPoints array must equal dim 0 of refPoints array");
2947 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)physPoints.dimension(2) != (size_t)cellTopo.getDimension()), std::invalid_argument,
2948 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 2 (spatial dimension) does not match cell dimension ");
2952 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(physPoints) != 2), std::invalid_argument,
2953 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): rank = 2 required for physPoints array");
2955 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)physPoints.dimension(0) != (size_t)refPoints.dimension(0)), std::invalid_argument,
2956 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 0 (number of points) of physPoints array must equal dim 0 of refPoints array");
2958 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)physPoints.dimension(1) != (size_t)cellTopo.getDimension()), std::invalid_argument,
2959 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 1 (spatial dimension) does not match cell dimension ");
2963 else if(getrank(refPoints) == 3) {
2966 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)refPoints.dimension(0) !=(size_t) cellWorkset.dimension(0) ), std::invalid_argument,
2967 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 0 (number of cells) of refPoints and cellWorkset arraya are required to match ");
2969 TEUCHOS_TEST_FOR_EXCEPTION( (refPoints.dimension(1) <= 0), std::invalid_argument,
2970 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 1 (number of points) >= 1 required for refPoints array ");
2972 TEUCHOS_TEST_FOR_EXCEPTION( ((
size_t)refPoints.dimension(2) != (size_t)cellTopo.getDimension() ), std::invalid_argument,
2973 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): dim 2 (spatial dimension) of refPoints array does not match cell dimension");
2976 TEUCHOS_TEST_FOR_EXCEPTION( (whichCell != -1), std::invalid_argument,
2977 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): default value is required for rank-3 refPoints array");
2980 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankMatch(errmsg, refPoints, physPoints), std::invalid_argument, errmsg );
2981 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, refPoints, physPoints), std::invalid_argument, errmsg);
2985 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(refPoints) == 2) || (getrank(refPoints) == 3) ), std::invalid_argument,
2986 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToPhysicalFrame): rank = 2 or 3 required for refPoints array");
2989 template<
class Scalar>
2990 template<
class ArrayRefPo
int,
class ArrayPhysPo
int,
class ArrayCell>
2992 const ArrayPhysPoint & physPoints,
2993 const ArrayCell & cellWorkset,
2994 const shards::CellTopology & cellTopo,
2995 const int& whichCell)
2997 std::string errmsg =
">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame):";
2998 std::string errmsg1 =
">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame):";
3001 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(cellWorkset) != 3), std::invalid_argument,
3002 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame): rank = 3 required for cellWorkset array");
3004 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(0)) <= 0), std::invalid_argument,
3005 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame): dim 0 (number of cells) >= 1 required for cellWorkset array");
3007 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(1)) != (
size_t)cellTopo.getSubcellCount(0) ), std::invalid_argument,
3008 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame): dim 1 (number of cell nodes) of cellWorkset array does not match cell topology");
3010 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(2)) != (
size_t)cellTopo.getDimension() ), std::invalid_argument,
3011 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame): dim 2 (spatial dimension) of cellWorkset array does not match cell dimension");
3014 TEUCHOS_TEST_FOR_EXCEPTION( !( ( (0 <= whichCell ) && ((
size_t)whichCell <(
size_t) cellWorkset.dimension(0) ) ) || (whichCell == -1) ), std::invalid_argument,
3015 ">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame): whichCell = -1 or a valid cell ordinal is required.");
3019 if(whichCell == -1) {
3021 errmsg1 +=
" default value of whichCell requires rank-3 arrays:";
3025 errmsg1 +=
" rank-2 arrays required when whichCell is valid cell ordinal";
3028 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankRange(errmsg1, refPoints, validRank,validRank), std::invalid_argument, errmsg1);
3029 TEUCHOS_TEST_FOR_EXCEPTION( !requireRankMatch(errmsg1, physPoints, refPoints), std::invalid_argument, errmsg1);
3030 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg1, refPoints, physPoints), std::invalid_argument, errmsg1);
3035 template<
class Scalar>
3036 template<
class ArrayRefPo
int,
class ArrayInitGuess,
class ArrayPhysPo
int,
class ArrayCell>
3038 const ArrayInitGuess & initGuess,
3039 const ArrayPhysPoint & physPoints,
3040 const ArrayCell & cellWorkset,
3041 const shards::CellTopology & cellTopo,
3042 const int& whichCell)
3045 validateArguments_mapToReferenceFrame(refPoints, physPoints, cellWorkset, cellTopo, whichCell);
3048 std::string errmsg =
">>> ERROR (Intrepid::CellTools::validateArguments_mapToReferenceFrame):";
3049 TEUCHOS_TEST_FOR_EXCEPTION( !requireDimensionMatch(errmsg, initGuess, physPoints), std::invalid_argument, errmsg);
3053 template<
class Scalar>
3054 template<
class ArrayIncl,
class ArrayPo
int,
class ArrayCell>
3056 const ArrayPoint & physPoints,
3057 const ArrayCell & cellWorkset,
3058 const int & whichCell,
3059 const shards::CellTopology & cell)
3062 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(cellWorkset) != 3), std::invalid_argument,
3063 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): rank = 3 required for cellWorkset array");
3065 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(0)) <= 0), std::invalid_argument,
3066 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 0 (number of cells) >= 1 required for cellWorkset array");
3068 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(1)) != (
size_t)cell.getSubcellCount(0) ), std::invalid_argument,
3069 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 1 (number of cell nodes) of cellWorkset array does not match cell topology");
3071 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(cellWorkset.dimension(2)) != (
size_t)cell.getDimension() ), std::invalid_argument,
3072 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 2 (spatial dimension) of cellWorkset array does not match cell dimension");
3076 TEUCHOS_TEST_FOR_EXCEPTION( !( ( (0 <= whichCell ) && (whichCell < cellWorkset.dimension(0) ) ) || (whichCell == -1) ), std::invalid_argument,
3077 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): whichCell = -1 or a valid cell ordinal is required.");
3081 if(getrank(physPoints) == 2) {
3083 TEUCHOS_TEST_FOR_EXCEPTION( (whichCell == -1), std::invalid_argument,
3084 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): whichCell = a valid cell ordinal is required with rank-2 input array.");
3086 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(physPoints.dimension(0)) <= 0), std::invalid_argument,
3087 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 0 (number of points) >= 1 required for physPoints array ");
3089 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(physPoints.dimension(1)) != (
size_t)cell.getDimension() ), std::invalid_argument,
3090 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 1 (spatial dimension) of physPoints array does not match cell dimension");
3093 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(inCell) != 1), std::invalid_argument,
3094 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): rank = 1 required for inCell array");
3096 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(inCell.dimension(0)) != static_cast<size_t>(physPoints.dimension(0))), std::invalid_argument,
3097 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 0 (number of points) of inCell array must equal dim 0 of physPoints array");
3100 else if (getrank(physPoints) == 3){
3102 TEUCHOS_TEST_FOR_EXCEPTION( !(whichCell == -1), std::invalid_argument,
3103 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): whichCell = -1 is required with rank-3 input array.");
3105 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(physPoints.dimension(0)) != static_cast<size_t>(cellWorkset.dimension(0)) ), std::invalid_argument,
3106 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 0 (number of cells) of physPoints array must equal dim 0 of cellWorkset array ");
3108 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(physPoints.dimension(1)) <= 0), std::invalid_argument,
3109 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 1 (number of points) >= 1 required for physPoints array ");
3111 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(physPoints.dimension(2)) != (
size_t)cell.getDimension() ), std::invalid_argument,
3112 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 2 (spatial dimension) of physPoints array does not match cell dimension");
3115 TEUCHOS_TEST_FOR_EXCEPTION( (getrank(inCell) != 2), std::invalid_argument,
3116 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): rank = 2 required for inCell array");
3118 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(inCell.dimension(0)) != static_cast<size_t>(physPoints.dimension(0))), std::invalid_argument,
3119 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 0 (number of cells) of inCell array must equal dim 0 of physPoints array");
3121 TEUCHOS_TEST_FOR_EXCEPTION( (static_cast<size_t>(inCell.dimension(1)) != static_cast<size_t>(physPoints.dimension(1))), std::invalid_argument,
3122 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): dim 1 (number of points) of inCell array must equal dim 1 of physPoints array");
3125 TEUCHOS_TEST_FOR_EXCEPTION( !( (getrank(physPoints) == 2) && (getrank(physPoints) ==3) ), std::invalid_argument,
3126 ">>> ERROR (Intrepid::CellTools::validateArguments_checkPointwiseInclusion): rank = 2 or 3 required for points array");
3139 template<
class Scalar>
3141 const int subcellOrd,
3142 const shards::CellTopology & parentCell){
3145 int subcVertexCount = parentCell.getVertexCount(subcellDim, subcellOrd);
3146 int cellDim = parentCell.getDimension();
3152 getReferenceSubcellVertices(subcellVertices,
3159 <<
" Subcell " << std::setw(2) << subcellOrd
3160 <<
" is " << parentCell.getName(subcellDim, subcellOrd) <<
" with vertices = {";
3163 for(
int subcVertOrd = 0; subcVertOrd < subcVertexCount; subcVertOrd++){
3167 for(
int dim = 0; dim < (int)parentCell.getDimension(); dim++){
3168 std::cout << subcellVertices(subcVertOrd, dim);
3169 if(dim < (
int)parentCell.getDimension()-1 ) { std::cout <<
","; }
3172 if(subcVertOrd < subcVertexCount - 1) { std::cout <<
", "; }
3178 template<
class Scalar>
3179 template<
class ArrayCell>
3181 const shards::CellTopology & parentCell,
3182 const int& pCellOrd,
3183 const int& subcellDim,
3184 const int& subcellOrd,
3185 const int& fieldWidth){
3188 int subcNodeCount = parentCell.getNodeCount(subcellDim, subcellOrd);
3189 int pCellDim = parentCell.getDimension();
3190 std::vector<int> subcNodeOrdinals(subcNodeCount);
3192 for(
int i = 0; i < subcNodeCount; i++){
3193 subcNodeOrdinals[i] = parentCell.getNodeMap(subcellDim, subcellOrd, i);
3199 <<
" Subcell " << subcellOrd <<
" on parent cell " << pCellOrd <<
" is "
3200 << parentCell.getName(subcellDim, subcellOrd) <<
" with node(s) \n ({";
3202 for(
int i = 0; i < subcNodeCount; i++){
3205 for(
int dim = 0; dim < pCellDim; dim++){
3207 << std::setw(fieldWidth) << std::right << cellWorkset(pCellOrd, subcNodeOrdinals[i], dim);
3208 if(dim < pCellDim - 1){ std::cout <<
","; }
3211 if(i < subcNodeCount - 1){ std::cout <<
", {"; }
3213 std::cout <<
")\n\n";
3221 template<
class Scalar>
3222 template<
class ArrayCVCoord,
class ArrayCellCoord>
3224 const ArrayCellCoord & cellCoords,
3225 const shards::CellTopology& primaryCell)
3229 int numCells = cellCoords.dimension(0);
3230 int numNodesPerCell = cellCoords.dimension(1);
3231 int spaceDim = cellCoords.dimension(2);
3234 int numEdgesPerCell = primaryCell.getEdgeCount();
3237 int numFacesPerCell = 0;
3239 numFacesPerCell = primaryCell.getFaceCount();
3244 getBarycenter(barycenter,cellCoords);
3247 for (
int icell = 0; icell < numCells; icell++){
3251 for (
int iedge = 0; iedge < numEdgesPerCell; iedge++){
3252 for (
int idim = 0; idim < spaceDim; idim++){
3254 int node0 = primaryCell.getNodeMap(1,iedge,0);
3255 int node1 = primaryCell.getNodeMap(1,iedge,1);
3256 edgeMidpts(iedge,idim) = (cellCoords(icell,node0,idim) +
3257 cellCoords(icell,node1,idim))/2.0;
3263 int numNodesPerFace;
3266 for (
int iface = 0; iface < numFacesPerCell; iface++){
3267 numNodesPerFace = primaryCell.getNodeCount(2,iface);
3269 for (
int idim = 0; idim < spaceDim; idim++){
3271 for (
int inode0 = 0; inode0 < numNodesPerFace; inode0++) {
3272 int node1 = primaryCell.getNodeMap(2,iface,inode0);
3273 faceMidpts(iface,idim) += cellCoords(icell,node1,idim)/numNodesPerFace;
3281 switch(primaryCell.getKey() ) {
3284 case shards::Triangle<3>::key:
3285 case shards::Quadrilateral<4>::key:
3287 for (
int inode = 0; inode < numNodesPerCell; inode++){
3288 for (
int idim = 0; idim < spaceDim; idim++){
3291 subCVCoords(icell,inode,0,idim) = cellCoords(icell,inode,idim);
3294 subCVCoords(icell,inode,1,idim) = edgeMidpts(inode,idim);
3297 subCVCoords(icell,inode,2,idim) = barycenter(icell,idim);
3300 int jnode = numNodesPerCell-1;
3301 if (inode > 0) jnode = inode - 1;
3302 subCVCoords(icell,inode,3,idim) = edgeMidpts(jnode,idim);
3309 case shards::Hexahedron<8>::key:
3311 for (
int idim = 0; idim < spaceDim; idim++){
3314 for (
int icount = 0; icount < 4; icount++){
3318 subCVCoords(icell,icount,0,idim) = cellCoords(icell,icount,idim);
3319 subCVCoords(icell,icount+4,4,idim) = cellCoords(icell,icount+4,idim);
3323 subCVCoords(icell,icount,1,idim) = edgeMidpts(icount,idim);
3324 subCVCoords(icell,icount+4,5,idim) = edgeMidpts(icount+4,idim);
3328 subCVCoords(icell,icount,2,idim) = faceMidpts(4,idim);
3329 subCVCoords(icell,icount+4,6,idim) = faceMidpts(5,idim);
3334 if (icount > 0) jcount = icount - 1;
3335 subCVCoords(icell,icount,3,idim) = edgeMidpts(jcount,idim);
3336 subCVCoords(icell,icount+4,7,idim) = edgeMidpts(jcount+4,idim);
3340 subCVCoords(icell,icount,4,idim) = edgeMidpts(icount+numNodesPerCell,idim);
3341 subCVCoords(icell,icount+4,0,idim) = edgeMidpts(icount+numNodesPerCell,idim);
3345 subCVCoords(icell,icount,5,idim) = faceMidpts(icount,idim);
3346 subCVCoords(icell,icount+4,1,idim) = faceMidpts(icount,idim);
3350 subCVCoords(icell,icount,6,idim) = barycenter(icell,idim);
3351 subCVCoords(icell,icount+4,2,idim) = barycenter(icell,idim);
3356 if (icount > 0) jcount = icount - 1;
3357 subCVCoords(icell,icount,7,idim) = faceMidpts(jcount,idim);
3358 subCVCoords(icell,icount+4,3,idim) = faceMidpts(jcount,idim);
3366 case shards::Tetrahedron<4>::key:
3368 for (
int idim = 0; idim < spaceDim; idim++){
3371 for (
int icount = 0; icount < 3; icount++){
3374 subCVCoords(icell,icount,0,idim) = cellCoords(icell,icount,idim);
3377 subCVCoords(icell,icount,1,idim) = edgeMidpts(icount,idim);
3380 subCVCoords(icell,icount,2,idim) = faceMidpts(3,idim);
3384 if (icount > 0) jcount = icount - 1;
3385 subCVCoords(icell,icount,3,idim) = edgeMidpts(jcount,idim);
3388 subCVCoords(icell,icount,4,idim) = edgeMidpts(icount+3,idim);
3391 subCVCoords(icell,icount,5,idim) = faceMidpts(icount,idim);
3394 subCVCoords(icell,icount,6,idim) = barycenter(icell,idim);
3398 if (icount > 0) jcount = icount - 1;
3399 subCVCoords(icell,icount,7,idim) = faceMidpts(jcount,idim);
3405 subCVCoords(icell,3,0,idim) = cellCoords(icell,3,idim);
3408 subCVCoords(icell,3,1,idim) = edgeMidpts(3,idim);
3411 subCVCoords(icell,3,2,idim) = faceMidpts(2,idim);
3414 subCVCoords(icell,3,3,idim) = edgeMidpts(5,idim);
3417 subCVCoords(icell,3,4,idim) = edgeMidpts(4,idim);
3420 subCVCoords(icell,3,5,idim) = faceMidpts(0,idim);
3423 subCVCoords(icell,3,6,idim) = barycenter(icell,idim);
3426 subCVCoords(icell,3,7,idim) = faceMidpts(1,idim);
3433 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
3434 ">>> ERROR (getSubCVCoords: invalid cell topology.");
3441 template<
class Scalar>
3442 template<
class ArrayCent,
class ArrayCellCoord>
3446 int numCells = cellCoords.dimension(0);
3447 int numVertsPerCell = cellCoords.dimension(1);
3448 int spaceDim = cellCoords.dimension(2);
3453 for (
int icell = 0; icell < numCells; icell++){
3458 for (
int inode = 0; inode < numVertsPerCell; inode++){
3460 int jnode = inode + 1;
3461 if (jnode >= numVertsPerCell) {
3465 Scalar area_mult = cellCoords(icell,inode,0)*cellCoords(icell,jnode,1)
3466 - cellCoords(icell,jnode,0)*cellCoords(icell,inode,1);
3467 cell_centroid(0) += (cellCoords(icell,inode,0) + cellCoords(icell,jnode,0))*area_mult;
3468 cell_centroid(1) += (cellCoords(icell,inode,1) + cellCoords(icell,jnode,1))*area_mult;
3470 area += 0.5*area_mult;
3473 barycenter(icell,0) = cell_centroid(0)/(6.0*area);
3474 barycenter(icell,1) = cell_centroid(1)/(6.0*area);
3483 for (
int icell = 0; icell < numCells; icell++){
3487 for (
int inode = 0; inode < numVertsPerCell; inode++){
3488 for (
int idim = 0; idim < spaceDim; idim++){
3489 cell_centroid(idim) += cellCoords(icell,inode,idim)/numVertsPerCell;
3492 for (
int idim = 0; idim < spaceDim; idim++){
3493 barycenter(icell,idim) = cell_centroid(idim);
3502 #if defined(Intrepid_SHOW_DEPRECATED_WARNINGS)
3504 #warning "The Intrepid package is deprecated"
int size() const
Returns size of the FieldContainer defined as the product of its dimensions.
Implementation of the serendipity-family H(grad)-compatible FEM basis of degree 2 on a Hexahedron cel...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
#define INTREPID_MAX_NEWTON
Maximum number of Newton iterations used internally in methods such as computing the action of the in...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
Implementation of an H(grad)-compatible FEM basis of degree 2 on Wedge cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Pyramid cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell...
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell...