53 #include "Teuchos_GlobalMPISession.hpp"
55 #include "Shards_CellTopology.hpp"
57 #include "Teuchos_RCP.hpp"
60 using namespace Intrepid;
61 using namespace shards;
63 int main(
int argc,
char *argv[]) {
65 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
68 typedef shards::CellTopology CellTopology;
71 <<
"===============================================================================\n" \
73 <<
"| Example use of the CellTools class |\n" \
75 <<
"| 1) Reference edge parametrizations |\n" \
76 <<
"| 2) Reference face parametrizations |\n" \
78 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) |\n" \
79 <<
"| Denis Ridzal (dridzal@sandia.gov), or |\n" \
80 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
82 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
83 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
85 <<
"===============================================================================\n"\
87 <<
"| Reference edge parametrizations map [-1,1] to the edges of reference cells. |\n"\
88 <<
"| They are used to define, e.g., integration points on the edges of 2D and 3D |\n"\
89 <<
"| reference cells. Edge parametrizations for special 2D cells such as Beam |\n"\
90 <<
"| and ShellLine, are also supported. |\n"\
91 <<
"===============================================================================\n";
120 CellTopology edgeParam(shards::getCellTopologyData<shards::Line<2> >() );
124 Teuchos::RCP<Cubature<double> > edgeParamCubature = cubFactory.
create(edgeParam, 6);
128 int cubDim = edgeParamCubature -> getDimension();
129 int numCubPoints = edgeParamCubature -> getNumPoints();
133 edgeParamCubature -> getCubature(edgeParamCubPts, edgeParamCubWts);
138 <<
"===============================================================================\n"\
139 <<
"| EXAMPLE 1.1 |\n"
140 <<
"| Edge parametrizations for standard 2D cells: Triangle |\n"\
141 <<
"===============================================================================\n";
144 CellTopology triangle_3(getCellTopologyData<Triangle<3> >() );
152 for(
int edgeOrd = 0; edgeOrd < (int)triangle_3.getEdgeCount(); edgeOrd++){
154 CellTools::mapToReferenceSubcell(triEdgePoints,
161 CellTools::printSubcellVertices(1, edgeOrd, triangle_3);
163 for(
int pt = 0; pt < numCubPoints; pt++){
164 std::cout <<
"\t Parameter point "
165 << std::setw(12) << std::right << edgeParamCubPts(pt, 0) << std::setw(10) <<
" --> " <<
"("
166 << std::setw(10) << std::right << triEdgePoints(pt, 0) <<
" , "
167 << std::setw(10) << std::right << triEdgePoints(pt, 1) <<
")\n";
175 <<
"===============================================================================\n"\
176 <<
"| EXAMPLE 1.2 |\n"
177 <<
"| Edge parametrizations for standard 2D cells: Quadrilateral |\n"\
178 <<
"===============================================================================\n";
181 CellTopology quad_4(getCellTopologyData<Quadrilateral<4> >() );
189 for(
int edgeOrd = 0; edgeOrd < (int)quad_4.getEdgeCount(); edgeOrd++){
191 CellTools::mapToReferenceSubcell(quadEdgePoints,
198 CellTools::printSubcellVertices(1, edgeOrd, quad_4);
200 for(
int pt = 0; pt < numCubPoints; pt++){
201 std::cout <<
"\t Parameter point "
202 << std::setw(12) << std::right << edgeParamCubPts(pt, 0) << std::setw(10) <<
" --> " <<
"("
203 << std::setw(10) << std::right << quadEdgePoints(pt, 0) <<
" , "
204 << std::setw(10) << std::right << quadEdgePoints(pt, 1) <<
")\n";
238 CellTopology triFaceParam(shards::getCellTopologyData<shards::Triangle<3> >() );
239 CellTopology quadFaceParam(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
243 Teuchos::RCP<Cubature<double> > triFaceParamCubature = cubFactory.
create(triFaceParam, 3);
244 Teuchos::RCP<Cubature<double> > quadFaceParamCubature = cubFactory.
create(quadFaceParam, 3);
248 int triFaceCubDim = triFaceParamCubature -> getDimension();
249 int triFaceNumCubPts = triFaceParamCubature -> getNumPoints();
253 triFaceParamCubature -> getCubature(triFaceParamCubPts, triFaceParamCubWts);
257 int quadFaceCubDim = quadFaceParamCubature -> getDimension();
258 int quadFaceNumCubPts = quadFaceParamCubature -> getNumPoints();
262 quadFaceParamCubature -> getCubature(quadFaceParamCubPts, quadFaceParamCubWts);
267 <<
"===============================================================================\n"\
268 <<
"| EXAMPLE 2.1 |\n"
269 <<
"| Face parametrizations for standard 3D cells: Tetrahedron |\n"\
270 <<
"===============================================================================\n";
273 CellTopology tet_4(getCellTopologyData<Tetrahedron<4> >() );
281 for(
int faceOrd = 0; faceOrd < (int)tet_4.getSideCount(); faceOrd++){
283 CellTools::mapToReferenceSubcell(tetFacePoints,
290 CellTools::printSubcellVertices(2, faceOrd, tet_4);
292 for(
int pt = 0; pt < triFaceNumCubPts; pt++){
293 std::cout <<
"\t Parameter point ("
294 << std::setw(10) << std::right << triFaceParamCubPts(pt, 0) <<
" , "
295 << std::setw(10) << std::right << triFaceParamCubPts(pt, 1) <<
") "
296 << std::setw(10) <<
" --> " <<
"("
297 << std::setw(10) << std::right << tetFacePoints(pt, 0) <<
" , "
298 << std::setw(10) << std::right << tetFacePoints(pt, 1) <<
" , "
299 << std::setw(10) << std::right << tetFacePoints(pt, 2) <<
")\n";
307 <<
"===============================================================================\n"\
308 <<
"| EXAMPLE 2.2 |\n"
309 <<
"| Face parametrizations for standard 3D cells: Wedge |\n"\
310 <<
"| Example of a reference cell that has two different kinds of faces |\n"\
311 <<
"===============================================================================\n";
315 CellTopology wedge_6(getCellTopologyData<Wedge<6> >() );
328 for(
int faceOrd = 0; faceOrd < (int)wedge_6.getSideCount(); faceOrd++){
331 CellTools::printSubcellVertices(2, faceOrd, wedge_6);
333 if( wedge_6.getKey(2, faceOrd) == shards::Triangle<3>::key ){
334 CellTools::mapToReferenceSubcell(wedgeTriFacePoints,
340 for(
int pt = 0; pt < triFaceNumCubPts; pt++){
341 std::cout <<
"\t Parameter point ("
342 << std::setw(10) << std::right << triFaceParamCubPts(pt, 0) <<
" , "
343 << std::setw(10) << std::right << triFaceParamCubPts(pt, 1) <<
") "
344 << std::setw(10) <<
" --> " <<
"("
345 << std::setw(10) << std::right << wedgeTriFacePoints(pt, 0) <<
" , "
346 << std::setw(10) << std::right << wedgeTriFacePoints(pt, 1) <<
" , "
347 << std::setw(10) << std::right << wedgeTriFacePoints(pt, 2) <<
")\n";
351 else if(wedge_6.getKey(2, faceOrd) == shards::Quadrilateral<4>::key) {
352 CellTools::mapToReferenceSubcell(wedgeQuadFacePoints,
358 for(
int pt = 0; pt < quadFaceNumCubPts; pt++){
359 std::cout <<
"\t Parameter point ("
360 << std::setw(10) << std::right << quadFaceParamCubPts(pt, 0) <<
" , "
361 << std::setw(10) << std::right << quadFaceParamCubPts(pt, 1) <<
") "
362 << std::setw(10) <<
" --> " <<
"("
363 << std::setw(10) << std::right << wedgeQuadFacePoints(pt, 0) <<
" , "
364 << std::setw(10) << std::right << wedgeQuadFacePoints(pt, 1) <<
" , "
365 << std::setw(10) << std::right << wedgeQuadFacePoints(pt, 2) <<
")\n";
370 std::cout <<
" Invalid face encountered \n";
Header file for utility class to provide multidimensional containers.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
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.