52 #include "Teuchos_oblackholestream.hpp"
53 #include "Teuchos_RCP.hpp"
54 #include "Teuchos_GlobalMPISession.hpp"
55 #include "Teuchos_ScalarTraits.hpp"
58 using namespace Intrepid;
80 const shards::CellTopology& parentCell,
84 const Teuchos::RCP<std::ostream>& outStream);
89 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double> > >
90 getIntrepidBasis(
const shards::CellTopology &cellTopo)
92 Teuchos::RCP< Basis< double, FieldContainer<double> > > HGRAD_Basis;
94 switch( cellTopo.getKey() ){
97 case shards::Line<2>::key:
101 case shards::Triangle<3>::key:
105 case shards::Quadrilateral<4>::key:
109 case shards::Tetrahedron<4>::key:
113 case shards::Hexahedron<8>::key:
117 case shards::Wedge<6>::key:
121 case shards::Pyramid<5>::key:
126 case shards::Triangle<6>::key:
130 case shards::Quadrilateral<9>::key:
134 case shards::Tetrahedron<10>::key:
138 case shards::Tetrahedron<11>::key:
142 case shards::Hexahedron<20>::key:
146 case shards::Hexahedron<27>::key:
150 case shards::Wedge<15>::key:
154 case shards::Wedge<18>::key:
158 case shards::Pyramid<13>::key:
163 case shards::Quadrilateral<8>::key:
164 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
165 ">>> ERROR (getIntrepidBasis): Cell topology not supported. ");
169 case shards::Line<3>::key:
170 case shards::Beam<2>::key:
171 case shards::Beam<3>::key:
172 case shards::ShellLine<2>::key:
173 case shards::ShellLine<3>::key:
174 case shards::ShellTriangle<3>::key:
175 case shards::ShellTriangle<6>::key:
176 case shards::ShellQuadrilateral<4>::key:
177 case shards::ShellQuadrilateral<8>::key:
178 case shards::ShellQuadrilateral<9>::key:
179 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
180 ">>> ERROR (getIntrepidBasis): Cell topology not supported. ");
183 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::invalid_argument,
184 ">>> ERROR (getIntrepidBasis): Cell topology not supported.");
189 int main(
int argc,
char *argv[]) {
191 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
194 typedef shards::CellTopology CellTopology;
197 int iprint = argc - 1;
199 Teuchos::RCP<std::ostream> outStream;
200 Teuchos::oblackholestream bhs;
203 outStream = Teuchos::rcp(&std::cout,
false);
205 outStream = Teuchos::rcp(&bhs,
false);
208 Teuchos::oblackholestream oldFormatState;
209 oldFormatState.copyfmt(std::cout);
212 <<
"===============================================================================\n" \
214 <<
"| Unit Test CellTools |\n" \
216 <<
"| 1) Edge parametrizations |\n" \
217 <<
"| 2) Face parametrizations |\n" \
218 <<
"| 3) Edge tangents |\n" \
219 <<
"| 4) Face tangents and normals |\n" \
221 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) |\n" \
222 <<
"| Denis Ridzal (dridzal@sandia.gov), or |\n" \
223 <<
"| Kara Peterson (kjpeter@sandia.gov) |\n" \
225 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
226 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
228 <<
"===============================================================================\n";
241 simplex_2(0, 0) = 0.0; simplex_2(0, 1) = 0.0;
242 simplex_2(1, 0) = 1.0; simplex_2(1, 1) = 0.0;
243 simplex_2(2, 0) = 0.0; simplex_2(2, 1) = 1.0;
248 cube_2(0, 0) = -1.0; cube_2(0, 1) = -1.0;
249 cube_2(1, 0) = 1.0; cube_2(1, 1) = -1.0;
250 cube_2(2, 0) = 1.0; cube_2(2, 1) = 1.0;
251 cube_2(3, 0) = -1.0; cube_2(3, 1) = 1.0;
255 std::vector<shards::CellTopology> allTopologies;
256 shards::getTopologies(allTopologies);
266 <<
"===============================================================================\n"\
267 <<
"| Test 1: edge parametrizations: |\n"\
268 <<
"===============================================================================\n\n";
273 for(
int topoOrd = 0; topoOrd < (int)allTopologies.size(); topoOrd++){
276 if(allTopologies[topoOrd].getDimension() > 1 && CellTools::hasReferenceCell(allTopologies[topoOrd]) ){
277 *outStream <<
" Testing edge parametrization for " << allTopologies[topoOrd].getName() <<
"\n";
279 allTopologies[topoOrd],
290 <<
"===============================================================================\n"\
291 <<
"| Test 2: face parametrizations: |\n"\
292 <<
"===============================================================================\n\n";
297 for(
int topoOrd = 0; topoOrd < (int)allTopologies.size(); topoOrd++){
300 if(allTopologies[topoOrd].getDimension() > 2 && CellTools::hasReferenceCell(allTopologies[topoOrd]) ){
301 *outStream <<
" Testing face parametrization for cell topology " << allTopologies[topoOrd].getName() <<
"\n";
303 allTopologies[topoOrd],
319 std::vector<shards::CellTopology> standardBaseTopologies;
320 shards::getTopologies(standardBaseTopologies, 4, shards::STANDARD_CELL, shards::BASE_TOPOLOGY);
323 CellTopology paramEdge (shards::getCellTopologyData<shards::Line<2> >() );
324 CellTopology paramTriFace (shards::getCellTopologyData<shards::Triangle<3> >() );
325 CellTopology paramQuadFace(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
332 <<
"===============================================================================\n"\
333 <<
"| Test 3: edge tangents/normals for stand. cells with base topologies: |\n"\
334 <<
"===============================================================================\n\n";
336 std::vector<shards::CellTopology>::iterator cti;
339 Teuchos::RCP<Cubature<double> > edgeCubature = cubFactory.
create(paramEdge, 6);
340 int cubDim = edgeCubature -> getDimension();
341 int numCubPoints = edgeCubature -> getNumPoints();
346 edgeCubature -> getCubature(paramEdgePoints, paramEdgeWeights);
349 for(cti = standardBaseTopologies.begin(); cti !=standardBaseTopologies.end(); ++cti){
352 if( ( (*cti).getDimension() >= 2) && ( (*cti).getKey() != shards::Pyramid<5>::key) ){
354 int cellDim = (*cti).getDimension();
355 int vCount = (*cti).getVertexCount();
357 CellTools::getReferenceSubcellVertices(refCellVertices, cellDim, 0, (*cti) );
359 *outStream <<
" Testing edge tangents";
360 if(cellDim == 2) { *outStream <<
" and normals"; }
361 *outStream <<
" for cell topology " << (*cti).getName() <<
"\n";
369 for(
int v = 0; v < vCount; v++){
370 for(
int d = 0; d < cellDim; d++){
371 double delta = Teuchos::ScalarTraits<double>::random()/8.0;
372 physCellVertices(0, v, d) = refCellVertices(v, d) + delta;
383 for(
int edgeOrd = 0; edgeOrd < (int)(*cti).getEdgeCount(); edgeOrd++){
390 CellTools::mapToReferenceSubcell(refEdgePoints, paramEdgePoints, 1, edgeOrd, (*cti) );
391 CellTools::setJacobian(edgePointsJacobians, refEdgePoints, physCellVertices, getIntrepidBasis(*cti) );
392 CellTools::getPhysicalEdgeTangents(edgePointTangents, edgePointsJacobians, edgeOrd, (*cti));
400 int v0ord = (*cti).getNodeMap(1, edgeOrd, 0);
401 int v1ord = (*cti).getNodeMap(1, edgeOrd, 1);
403 for(
int pt = 0; pt < numCubPoints; pt++){
408 for(
int d = 0; d < cellDim; d++){
409 edgeBenchmarkTangents(d) = (physCellVertices(0, v1ord, d) - physCellVertices(0, v0ord, d))/2.0;
412 if( abs(edgeBenchmarkTangents(d) - edgePointTangents(0, pt, d)) > INTREPID_THRESHOLD ){
415 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
416 <<
" Edge tangent computation by CellTools failed for: \n"
417 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
418 <<
" Edge ordinal = " << edgeOrd <<
"\n"
419 <<
" Edge point number = " << pt <<
"\n"
420 <<
" Tangent coordinate = " << d <<
"\n"
421 <<
" CellTools value = " << edgePointTangents(0, pt, d) <<
"\n"
422 <<
" Benchmark value = " << edgeBenchmarkTangents(d) <<
"\n\n";
428 CellTools::getPhysicalSideNormals(edgePointNormals, edgePointsJacobians, edgeOrd, (*cti));
429 if( abs(edgeBenchmarkTangents(1) - edgePointNormals(0, pt, 0)) > INTREPID_THRESHOLD ){
432 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
433 <<
" Edge Normal computation by CellTools failed for: \n"
434 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
435 <<
" Edge ordinal = " << edgeOrd <<
"\n"
436 <<
" Edge point number = " << pt <<
"\n"
437 <<
" Normal coordinate = " << 0 <<
"\n"
438 <<
" CellTools value = " << edgePointNormals(0, pt, 0) <<
"\n"
439 <<
" Benchmark value = " << edgeBenchmarkTangents(1) <<
"\n\n";
441 if( abs(edgeBenchmarkTangents(0) + edgePointNormals(0, pt, 1)) > INTREPID_THRESHOLD ){
444 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
445 <<
" Edge Normal computation by CellTools failed for: \n"
446 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
447 <<
" Edge ordinal = " << edgeOrd <<
"\n"
448 <<
" Edge point number = " << pt <<
"\n"
449 <<
" Normal coordinate = " << 1 <<
"\n"
450 <<
" CellTools value = " << edgePointNormals(0, pt, 1) <<
"\n"
451 <<
" Benchmark value = " << -edgeBenchmarkTangents(0) <<
"\n\n";
463 <<
"===============================================================================\n"\
464 <<
"| Test 4: face/side normals for stand. 3D cells with base topologies: | |\n"\
465 <<
"===============================================================================\n\n";
469 Teuchos::RCP<Cubature<double> > triFaceCubature = cubFactory.
create(paramTriFace, 6);
470 Teuchos::RCP<Cubature<double> > quadFaceCubature = cubFactory.
create(paramQuadFace, 6);
472 int faceCubDim = triFaceCubature -> getDimension();
473 int numTriFaceCubPoints = triFaceCubature -> getNumPoints();
474 int numQuadFaceCubPoints = quadFaceCubature -> getNumPoints();
482 triFaceCubature -> getCubature(paramTriFacePoints, paramTriFaceWeights);
483 quadFaceCubature -> getCubature(paramQuadFacePoints, paramQuadFaceWeights);
486 for(cti = standardBaseTopologies.begin(); cti !=standardBaseTopologies.end(); ++cti){
489 if( ( (*cti).getDimension() == 3) && ( (*cti).getKey() != shards::Pyramid<5>::key) ){
491 int cellDim = (*cti).getDimension();
492 int vCount = (*cti).getVertexCount();
494 CellTools::getReferenceSubcellVertices(refCellVertices, cellDim, 0, (*cti) );
496 *outStream <<
" Testing face/side normals for cell topology " << (*cti).getName() <<
"\n";
503 for(
int v = 0; v < vCount; v++){
504 for(
int d = 0; d < cellDim; d++){
505 double delta = Teuchos::ScalarTraits<double>::random()/8.0;
506 physCellVertices(0, v, d) = refCellVertices(v, d) + delta;
523 for(
int faceOrd = 0; faceOrd < (int)(*cti).getSideCount(); faceOrd++){
527 switch( (*cti).getCellTopologyData(2, faceOrd) -> key ) {
529 case shards::Triangle<3>::key:
532 CellTools::mapToReferenceSubcell(refTriFacePoints, paramTriFacePoints, 2, faceOrd, (*cti) );
533 CellTools::setJacobian(triFacePointsJacobians, refTriFacePoints, physCellVertices, getIntrepidBasis(*cti) );
534 CellTools::getPhysicalFaceNormals(triFacePointNormals, triFacePointsJacobians, faceOrd, (*cti));
535 CellTools::getPhysicalSideNormals(triSidePointNormals, triFacePointsJacobians, faceOrd, (*cti));
542 int v0ord = (*cti).getNodeMap(2, faceOrd, 0);
543 int v1ord = (*cti).getNodeMap(2, faceOrd, 1);
544 int v2ord = (*cti).getNodeMap(2, faceOrd, 2);
548 for(
int pt = 0; pt < numTriFaceCubPoints; pt++){
550 for(
int d = 0; d < cellDim; d++){
551 tanX(d) = (physCellVertices(0, v1ord, d) - physCellVertices(0, v0ord, d));
552 tanY(d) = (physCellVertices(0, v2ord, d) - physCellVertices(0, v0ord, d));
558 for(
int d = 0; d < cellDim; d++){
561 if( abs(faceNormal(d) - triFacePointNormals(0, pt, d)) > INTREPID_THRESHOLD ){
564 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
565 <<
" Face normal computation by CellTools failed for: \n"
566 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
567 <<
" Face Topology = " << (*cti).getCellTopologyData(2, faceOrd) -> name <<
"\n"
568 <<
" Face ordinal = " << faceOrd <<
"\n"
569 <<
" Face point number = " << pt <<
"\n"
570 <<
" Normal coordinate = " << d <<
"\n"
571 <<
" CellTools value = " << triFacePointNormals(0, pt, d)
572 <<
" Benchmark value = " << faceNormal(d) <<
"\n\n";
575 if( abs(faceNormal(d) - triSidePointNormals(0, pt, d)) > INTREPID_THRESHOLD ){
578 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
579 <<
" Side normal computation by CellTools failed for: \n"
580 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
581 <<
" Side Topology = " << (*cti).getCellTopologyData(2, faceOrd) -> name <<
"\n"
582 <<
" Side ordinal = " << faceOrd <<
"\n"
583 <<
" Side point number = " << pt <<
"\n"
584 <<
" Normal coordinate = " << d <<
"\n"
585 <<
" CellTools value = " << triSidePointNormals(0, pt, d)
586 <<
" Benchmark value = " << faceNormal(d) <<
"\n\n";
593 case shards::Quadrilateral<4>::key:
596 CellTools::mapToReferenceSubcell(refQuadFacePoints, paramQuadFacePoints, 2, faceOrd, (*cti) );
597 CellTools::setJacobian(quadFacePointsJacobians, refQuadFacePoints, physCellVertices, getIntrepidBasis(*cti) );
598 CellTools::getPhysicalFaceNormals(quadFacePointNormals, quadFacePointsJacobians, faceOrd, (*cti));
599 CellTools::getPhysicalSideNormals(quadSidePointNormals, quadFacePointsJacobians, faceOrd, (*cti));
608 int v0ord = (*cti).getNodeMap(2, faceOrd, 0);
609 int v1ord = (*cti).getNodeMap(2, faceOrd, 1);
610 int v2ord = (*cti).getNodeMap(2, faceOrd, 2);
611 int v3ord = (*cti).getNodeMap(2, faceOrd, 3);
614 for(
int pt = 0; pt < numTriFaceCubPoints; pt++){
616 for(
int d = 0; d < cellDim; d++){
617 tanX(d) = (physCellVertices(0, v0ord, d)*(-1.0 + paramQuadFacePoints(pt,1) ) +
618 physCellVertices(0, v1ord, d)*( 1.0 - paramQuadFacePoints(pt,1) ) +
619 physCellVertices(0, v2ord, d)*( 1.0 + paramQuadFacePoints(pt,1) ) +
620 physCellVertices(0, v3ord, d)*(-1.0 - paramQuadFacePoints(pt,1) ) )/4.0;
622 tanY(d) = (physCellVertices(0, v0ord, d)*(-1.0 + paramQuadFacePoints(pt,0) ) +
623 physCellVertices(0, v1ord, d)*(-1.0 - paramQuadFacePoints(pt,0) ) +
624 physCellVertices(0, v2ord, d)*( 1.0 + paramQuadFacePoints(pt,0) ) +
625 physCellVertices(0, v3ord, d)*( 1.0 - paramQuadFacePoints(pt,0) ) )/4.0;
630 for(
int d = 0; d < cellDim; d++){
633 if( abs(faceNormal(d) - quadFacePointNormals(0, pt, d)) > INTREPID_THRESHOLD ){
636 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
637 <<
" Face normal computation by CellTools failed for: \n"
638 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
639 <<
" Face Topology = " << (*cti).getCellTopologyData(2, faceOrd) -> name <<
"\n"
640 <<
" Face ordinal = " << faceOrd <<
"\n"
641 <<
" Face point number = " << pt <<
"\n"
642 <<
" Normal coordinate = " << d <<
"\n"
643 <<
" CellTools value = " << quadFacePointNormals(0, pt, d)
644 <<
" Benchmark value = " << faceNormal(d) <<
"\n\n";
647 if( abs(faceNormal(d) - quadSidePointNormals(0, pt, d)) > INTREPID_THRESHOLD ){
650 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
651 <<
" Side normal computation by CellTools failed for: \n"
652 <<
" Cell Topology = " << (*cti).getName() <<
"\n"
653 <<
" Side Topology = " << (*cti).getCellTopologyData(2, faceOrd) -> name <<
"\n"
654 <<
" Side ordinal = " << faceOrd <<
"\n"
655 <<
" Side point number = " << pt <<
"\n"
656 <<
" Normal coordinate = " << d <<
"\n"
657 <<
" CellTools value = " << quadSidePointNormals(0, pt, d)
658 <<
" Benchmark value = " << faceNormal(d) <<
"\n\n";
666 *outStream <<
" Face normals test failure: face topology not supported \n\n";
676 catch (std::logic_error err) {
677 *outStream << err.what() <<
"\n";
683 std::cout <<
"End Result: TEST FAILED\n";
685 std::cout <<
"End Result: TEST PASSED\n";
688 std::cout.copyfmt(oldFormatState);
696 const shards::CellTopology& parentCell,
700 const Teuchos::RCP<std::ostream>& outStream){
703 int cellDim = parentCell.getDimension();
704 int subcCount = parentCell.getSubcellCount(subcDim);
708 for(
int subcOrd = 0; subcOrd < subcCount; subcOrd++){
709 int subcVertexCount = parentCell.getVertexCount(subcDim, subcOrd);
731 else if(subcDim == 2) {
734 if(subcVertexCount == 3){
742 else if(subcVertexCount == 4){
752 for(
int subcVertOrd = 0; subcVertOrd < subcVertexCount; subcVertOrd++){
753 for(
int dim = 0; dim < cellDim; dim++){
755 if(mappedParamVertices(subcVertOrd, dim) != refSubcellVertices(subcVertOrd, dim) ) {
758 << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n"
759 <<
" Cell Topology = " << parentCell.getName() <<
"\n"
760 <<
" Parametrization of subcell " << subcOrd <<
" which is "
761 << parentCell.getName(subcDim,subcOrd) <<
" failed for vertex " << subcVertOrd <<
":\n"
762 <<
" parametrization map fails to map correctly coordinate " << dim <<
" of that vertex\n\n";
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.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Header file for utility class to provide multidimensional containers.
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.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
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...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
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.
A factory class that generates specific instances of cubatures.
Teuchos::RCP< Cubature< Scalar, ArrayPoint, ArrayWeight > > create(const shards::CellTopology &cellTopology, const std::vector< int > °ree)
Factory method.
void testSubcellParametrizations(int &errorFlag, const shards::CellTopology &parentCell, const FieldContainer< double > &subcParamVert_A, const FieldContainer< double > &subcParamVert_B, const int subcDim, const Teuchos::RCP< std::ostream > &outStream)
Maps the vertices of the subcell parametrization domain to that subcell.
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...