00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00035 #include "Intrepid_CellTools.hpp"
00036 #include "Intrepid_FieldContainer.hpp"
00037 #include "Intrepid_DefaultCubatureFactory.hpp"
00038 #include "Shards_CellTopology.hpp"
00039 #include "Teuchos_GlobalMPISession.hpp"
00040
00041 using namespace std;
00042 using namespace Intrepid;
00043
00044 int main(int argc, char *argv[]) {
00045
00046 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00047
00048 typedef CellTools<double> CellTools;
00049 typedef shards::CellTopology CellTopology;
00050
00051 std::cout \
00052 << "===============================================================================\n" \
00053 << "| |\n" \
00054 << "| Example use of the CellTools class |\n" \
00055 << "| |\n" \
00056 << "| 1) Definition of integration points on edge and face worksets |\n" \
00057 << "| 2) Computation of face normals and edge tangents on face and edge worksets |\n" \
00058 << "| 2) Computation of side normals on face and edge worksets |\n" \
00059 << "| |\n" \
00060 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) |\n" \
00061 << "| Denis Ridzal (dridzal@sandia.gov), or |\n" \
00062 << "| Kara Peterson (kjpeter@sandia.gov) |\n" \
00063 << "| |\n" \
00064 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00065 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00066 << "| |\n" \
00067 << "===============================================================================\n"\
00068 << "| EXAMPLE 1: Definition of integration points on a Triangle edge workset |\n"\
00069 << "===============================================================================\n";
00070
00071
00072
00073
00074
00075
00076 DefaultCubatureFactory<double> cubFactory;
00077
00078
00079 CellTopology paramEdge(shards::getCellTopologyData<shards::Line<2> >() );
00080
00081
00082 FieldContainer<double> paramEdgeWeights;
00083 FieldContainer<double> paramEdgePoints;
00084
00085
00086 FieldContainer<double> refEdgePoints;
00087
00088
00089 FieldContainer<double> worksetEdgePoints;
00090
00091
00092
00093
00094
00095
00096
00097
00098 CellTopology triangle_3( shards::getCellTopologyData<shards::Triangle<3> >() );
00099
00100
00101 int worksetSize = 4;
00102 int pCellNodeCount = triangle_3.getVertexCount();
00103 int pCellDim = triangle_3.getDimension();
00104
00105 FieldContainer<double> triNodes(worksetSize, pCellNodeCount, pCellDim);
00106
00107 triNodes(0, 0, 0) = 0.5; triNodes(0, 0, 1) = 2.0;
00108 triNodes(0, 1, 0) = 0.0; triNodes(0, 1, 1) = 1.0;
00109 triNodes(0, 2, 0) = 1.0; triNodes(0, 2, 1) = 1.0;
00110
00111 triNodes(1, 0, 0) = 1.0; triNodes(1, 0, 1) = 1.0;
00112 triNodes(1, 1, 0) = 1.0; triNodes(1, 1, 1) = 0.0;
00113 triNodes(1, 2, 0) = 0.0; triNodes(1, 2, 1) = 0.0;
00114
00115 triNodes(2, 0, 0) = 0.0; triNodes(2, 0, 1) = 0.0;
00116 triNodes(2, 1, 0) = 1.0; triNodes(2, 1, 1) = 1.0;
00117 triNodes(2, 2, 0) = 0.0; triNodes(2, 2, 1) = 1.0;
00118
00119 triNodes(3, 0, 0) = 1.0; triNodes(3, 0, 1) = 1.0;
00120 triNodes(3, 1, 0) = 2.5; triNodes(3, 1, 1) = 1.5;
00121 triNodes(3, 2, 0) = 0.5; triNodes(3, 2, 1) = 2.0;
00122
00123
00124 int subcellDim = 1;
00125 int subcellOrd = 2;
00126
00127
00128
00129
00130
00131
00132
00133
00134 Teuchos::RCP<Cubature<double> > triEdgeCubature = cubFactory.create(paramEdge, 6);
00135
00136
00137 int cubDim = triEdgeCubature -> getDimension();
00138 int numCubPoints = triEdgeCubature -> getNumPoints();
00139
00140
00141 paramEdgePoints.resize(numCubPoints, cubDim);
00142 paramEdgeWeights.resize(numCubPoints);
00143 triEdgeCubature -> getCubature(paramEdgePoints, paramEdgeWeights);
00144
00145
00146
00147
00148
00149
00150
00151
00152 refEdgePoints.resize(numCubPoints, pCellDim);
00153
00154
00155 worksetEdgePoints.resize(worksetSize, numCubPoints, pCellDim);
00156
00157
00158 CellTools::mapToReferenceSubcell(refEdgePoints,
00159 paramEdgePoints,
00160 subcellDim,
00161 subcellOrd,
00162 triangle_3);
00163
00164
00165 CellTools::mapToPhysicalFrame(worksetEdgePoints,
00166 refEdgePoints,
00167 triNodes,
00168 triangle_3);
00169
00170
00171
00172
00173
00174
00175 for(int pCell = 0; pCell < worksetSize; pCell++){
00176
00177 CellTools::printWorksetSubcell(triNodes, triangle_3, pCell, subcellDim, subcellOrd);
00178
00179 for(int pt = 0; pt < numCubPoints; pt++){
00180 std::cout << "\t 1D Gauss point "
00181 << std::setw(12) << std::right << paramEdgePoints(pt, 0) << std::setw(10) << " --> " << "("
00182 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00183 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 1) << ")\n";
00184 }
00185 std::cout << "\n";
00186 }
00187
00188
00189
00190 std::cout << "\n" \
00191 << "===============================================================================\n"\
00192 << "| EXAMPLE 2: Definition of integration points on a Quadrilateral edge workset|\n"\
00193 << "===============================================================================\n";
00194
00195
00196
00197
00198
00199
00200
00201 CellTopology quadrilateral_4( shards::getCellTopologyData<shards::Quadrilateral<4> >() );
00202
00203
00204 worksetSize = 2;
00205 pCellNodeCount = quadrilateral_4.getVertexCount();
00206 pCellDim = quadrilateral_4.getDimension();
00207
00208 FieldContainer<double> quadNodes(worksetSize, pCellNodeCount, pCellDim);
00209
00210 quadNodes(0, 0, 0) = 1.00; quadNodes(0, 0, 1) = 1.00;
00211 quadNodes(0, 1, 0) = 2.00; quadNodes(0, 1, 1) = 0.75;
00212 quadNodes(0, 2, 0) = 1.75; quadNodes(0, 2, 1) = 2.00;
00213 quadNodes(0, 3, 0) = 1.25; quadNodes(0, 3, 1) = 2.00;
00214
00215 quadNodes(1, 0, 0) = 2.00; quadNodes(1, 0, 1) = 0.75;
00216 quadNodes(1, 1, 0) = 3.00; quadNodes(1, 1, 1) = 1.25;
00217 quadNodes(1, 2, 0) = 2.75; quadNodes(1, 2, 1) = 2.25;
00218 quadNodes(1, 3, 0) = 1.75; quadNodes(1, 3, 1) = 2.00;
00219
00220
00221 subcellDim = 1;
00222 subcellOrd = 2;
00223
00224
00225
00226
00227
00228
00229 Teuchos::RCP<Cubature<double> > quadEdgeCubature = cubFactory.create(paramEdge, 6);
00230
00231
00232 cubDim = quadEdgeCubature -> getDimension();
00233 numCubPoints = quadEdgeCubature -> getNumPoints();
00234
00235
00236 paramEdgePoints.resize(numCubPoints, cubDim);
00237 paramEdgeWeights.resize(numCubPoints);
00238 quadEdgeCubature -> getCubature(paramEdgePoints, paramEdgeWeights);
00239
00240
00241
00242
00243
00244
00245 refEdgePoints.resize(numCubPoints, pCellDim);
00246
00247
00248 worksetEdgePoints.resize(worksetSize, numCubPoints, pCellDim);
00249
00250
00251 CellTools::mapToReferenceSubcell(refEdgePoints,
00252 paramEdgePoints,
00253 subcellDim,
00254 subcellOrd,
00255 quadrilateral_4);
00256
00257
00258 CellTools::mapToPhysicalFrame(worksetEdgePoints,
00259 refEdgePoints,
00260 quadNodes,
00261 quadrilateral_4);
00262
00263
00264
00265
00266 for(int pCell = 0; pCell < worksetSize; pCell++){
00267
00268 CellTools::printWorksetSubcell(quadNodes, quadrilateral_4, pCell, subcellDim, subcellOrd);
00269
00270 for(int pt = 0; pt < numCubPoints; pt++){
00271 std::cout << "\t 1D Gauss point "
00272 << std::setw(12) << std::right << paramEdgePoints(pt, 0) << std::setw(10) << " --> " << "("
00273 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00274 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 1) << ")\n";
00275 }
00276 std::cout << "\n";
00277 }
00278
00279
00280
00281 std::cout << "\n" \
00282 << "===============================================================================\n"\
00283 << "| EXAMPLE 3: Edge tangents at Gauss points on a Quadrilateral edge workset |\n"\
00284 << "===============================================================================\n";
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303 FieldContainer<double> worksetJacobians(worksetSize, numCubPoints, pCellDim, pCellDim);
00304
00305
00306 CellTools::setJacobian(worksetJacobians,
00307 refEdgePoints,
00308 quadNodes,
00309 quadrilateral_4);
00310
00311
00312
00313
00314 FieldContainer<double> edgeTangents(worksetSize, numCubPoints, pCellDim);
00315
00316
00317 CellTools::getPhysicalEdgeTangents(edgeTangents,
00318 worksetJacobians,
00319 subcellOrd,
00320 quadrilateral_4);
00321
00322
00323 std::cout
00324 << "Edge tangents computed by CellTools::getPhysicalEdgeTangents.\n"
00325 << "Local edge ordinal = " << subcellOrd <<"\n";
00326
00327 for(int pCell = 0; pCell < worksetSize; pCell++){
00328
00329 CellTools::printWorksetSubcell(quadNodes, quadrilateral_4, pCell, subcellDim, subcellOrd);
00330
00331 for(int pt = 0; pt < numCubPoints; pt++){
00332 std::cout << "\t 2D Gauss point: ("
00333 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00334 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 1) << ") "
00335 << std::setw(8) << " edge tangent: " << "("
00336 << std::setw(8) << std::right << edgeTangents(pCell, pt, 0) << ", "
00337 << std::setw(8) << std::right << edgeTangents(pCell, pt, 1) << ")\n";
00338 }
00339 std::cout << "\n";
00340 }
00341
00342 std::cout << "\n" \
00343 << "===============================================================================\n"\
00344 << "| EXAMPLE 4: Side normals at Gauss points on a Quadrilateral side workset |\n"\
00345 << "===============================================================================\n";
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355 FieldContainer<double> sideNormals(worksetSize, numCubPoints, pCellDim);
00356
00357
00358 CellTools::getPhysicalSideNormals(sideNormals,
00359 worksetJacobians,
00360 subcellOrd,
00361 quadrilateral_4);
00362
00363
00364 std::cout
00365 << "Side normals computed by CellTools::getPhysicalSideNormals.\n"
00366 << "Local side (edge) ordinal = " << subcellOrd <<"\n";
00367
00368 for(int pCell = 0; pCell < worksetSize; pCell++){
00369
00370 CellTools::printWorksetSubcell(quadNodes, quadrilateral_4, pCell, subcellDim, subcellOrd);
00371
00372 for(int pt = 0; pt < numCubPoints; pt++){
00373 std::cout << "\t 2D Gauss point: ("
00374 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00375 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 1) << ") "
00376 << std::setw(8) << " side normal: " << "("
00377 << std::setw(8) << std::right << sideNormals(pCell, pt, 0) << ", "
00378 << std::setw(8) << std::right << sideNormals(pCell, pt, 1) << ")\n";
00379 }
00380 std::cout << "\n";
00381 }
00382
00383
00384 std::cout << "\n" \
00385 << "===============================================================================\n"\
00386 << "| EXAMPLE 5: Definition of integration points on a Hexahedron edge workset |\n"\
00387 << "===============================================================================\n";
00388
00389
00390
00391
00392
00393
00394
00395 CellTopology hexahedron_8( shards::getCellTopologyData<shards::Hexahedron<8> >() );
00396
00397
00398 worksetSize = 1;
00399 pCellNodeCount = hexahedron_8.getVertexCount();
00400 pCellDim = hexahedron_8.getDimension();
00401
00402 FieldContainer<double> hexNodes(worksetSize, pCellNodeCount, pCellDim);
00403
00404 hexNodes(0, 0, 0) = 0.00; hexNodes(0, 0, 1) = 0.00, hexNodes(0, 0, 2) = 0.00;
00405 hexNodes(0, 1, 0) = 1.00; hexNodes(0, 1, 1) = 0.00, hexNodes(0, 1, 2) = 0.00;
00406 hexNodes(0, 2, 0) = 1.00; hexNodes(0, 2, 1) = 1.00, hexNodes(0, 2, 2) = 0.00;
00407 hexNodes(0, 3, 0) = 0.00; hexNodes(0, 3, 1) = 1.00, hexNodes(0, 3, 2) = 0.00;
00408
00409 hexNodes(0, 4, 0) = 0.00; hexNodes(0, 4, 1) = 0.00, hexNodes(0, 4, 2) = 1.00;
00410 hexNodes(0, 5, 0) = 1.00; hexNodes(0, 5, 1) = 0.00, hexNodes(0, 5, 2) = 1.00;
00411 hexNodes(0, 6, 0) = 1.00; hexNodes(0, 6, 1) = 1.00, hexNodes(0, 6, 2) = 0.75;
00412 hexNodes(0, 7, 0) = 0.00; hexNodes(0, 7, 1) = 1.00, hexNodes(0, 7, 2) = 1.00;
00413
00414
00415
00416
00417 FieldContainer<double> hexNodesAlt(worksetSize, pCellNodeCount, pCellDim);
00418
00419 hexNodesAlt(0, 0, 0) = 0.00; hexNodesAlt(0, 0, 1) = 0.00, hexNodesAlt(0, 0, 2) = 0.00;
00420 hexNodesAlt(0, 1, 0) = 1.00; hexNodesAlt(0, 1, 1) = 0.00, hexNodesAlt(0, 1, 2) = 0.00;
00421 hexNodesAlt(0, 2, 0) = 1.00; hexNodesAlt(0, 2, 1) = 1.00, hexNodesAlt(0, 2, 2) = 0.00;
00422 hexNodesAlt(0, 3, 0) = 0.00; hexNodesAlt(0, 3, 1) = 1.00, hexNodesAlt(0, 3, 2) = 0.00;
00423
00424 hexNodesAlt(0, 4, 0) = 0.00; hexNodesAlt(0, 4, 1) = 0.00, hexNodesAlt(0, 4, 2) = 1.00;
00425 hexNodesAlt(0, 5, 0) = 1.00; hexNodesAlt(0, 5, 1) = 0.00, hexNodesAlt(0, 5, 2) = 0.75;
00426 hexNodesAlt(0, 6, 0) = 1.00; hexNodesAlt(0, 6, 1) = 1.00, hexNodesAlt(0, 6, 2) = 0.50;
00427 hexNodesAlt(0, 7, 0) = 0.00; hexNodesAlt(0, 7, 1) = 1.00, hexNodesAlt(0, 7, 2) = 0.75;
00428
00429
00430 subcellDim = 1;
00431 subcellOrd = 5;
00432
00433
00434
00435
00436
00437
00438 Teuchos::RCP<Cubature<double> > hexEdgeCubature = cubFactory.create(paramEdge, 4);
00439
00440
00441 cubDim = hexEdgeCubature -> getDimension();
00442 numCubPoints = hexEdgeCubature -> getNumPoints();
00443
00444
00445 paramEdgePoints.resize(numCubPoints, cubDim);
00446 paramEdgeWeights.resize(numCubPoints);
00447 hexEdgeCubature -> getCubature(paramEdgePoints, paramEdgeWeights);
00448
00449
00450
00451
00452
00453
00454 refEdgePoints.resize(numCubPoints, pCellDim);
00455
00456
00457 worksetEdgePoints.resize(worksetSize, numCubPoints, pCellDim);
00458
00459
00460 CellTools::mapToReferenceSubcell(refEdgePoints,
00461 paramEdgePoints,
00462 subcellDim,
00463 subcellOrd,
00464 hexahedron_8);
00465
00466
00467 CellTools::mapToPhysicalFrame(worksetEdgePoints,
00468 refEdgePoints,
00469 hexNodes,
00470 hexahedron_8);
00471
00472
00473
00474
00475 for(int pCell = 0; pCell < worksetSize; pCell++){
00476
00477 CellTools::printWorksetSubcell(hexNodes, hexahedron_8, pCell, subcellDim, subcellOrd);
00478
00479 for(int pt = 0; pt < numCubPoints; pt++){
00480 std::cout << "\t 1D Gauss point "
00481 << std::setw(12) << std::right << paramEdgePoints(pt, 0) << std::setw(10) << " --> " << "("
00482 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00483 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 1) << ", "
00484 << std::setw(10) << std::right << worksetEdgePoints(pCell, pt, 2) << ")\n";
00485 }
00486 std::cout << "\n";
00487 }
00488
00489
00490
00491 std::cout << "\n" \
00492 << "===============================================================================\n"\
00493 << "| EXAMPLE 6: Edge tangents at Gauss points on a Hexahedron edge workset |\n"\
00494 << "===============================================================================\n";
00495
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513 worksetJacobians.resize(worksetSize, numCubPoints, pCellDim, pCellDim);
00514
00515
00516 CellTools::setJacobian(worksetJacobians,
00517 refEdgePoints,
00518 hexNodes,
00519 hexahedron_8);
00520
00521
00522
00523
00524
00525 edgeTangents.resize(worksetSize, numCubPoints, pCellDim);
00526
00527
00528 CellTools::getPhysicalEdgeTangents(edgeTangents,
00529 worksetJacobians,
00530 subcellOrd,
00531 hexahedron_8);
00532
00533
00534 std::cout
00535 << "Edge tangents computed by CellTools::getPhysicalEdgeTangents.\n"
00536 << "Local edge ordinal = " << subcellOrd <<"\n";
00537
00538 for(int pCell = 0; pCell < worksetSize; pCell++){
00539
00540 CellTools::printWorksetSubcell(hexNodes, hexahedron_8, pCell, subcellDim, subcellOrd);
00541
00542 for(int pt = 0; pt < numCubPoints; pt++){
00543 std::cout << "\t 3D Gauss point: ("
00544 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 0) << ", "
00545 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 1) << ", "
00546 << std::setw(8) << std::right << worksetEdgePoints(pCell, pt, 2) << ") "
00547 << std::setw(8) << " edge tangent: " << "("
00548 << std::setw(8) << std::right << edgeTangents(pCell, pt, 0) << ", "
00549 << std::setw(8) << std::right << edgeTangents(pCell, pt, 1) << ", "
00550 << std::setw(8) << std::right << edgeTangents(pCell, pt, 2) << ")\n";
00551 }
00552 std::cout << "\n";
00553 }
00554
00555
00556 std::cout << "\n" \
00557 << "===============================================================================\n"\
00558 << "| EXAMPLE 7: Definition of integration points on a Hexahedron face workset |\n"\
00559 << "===============================================================================\n";
00560
00561
00562
00563
00564
00565
00566 CellTopology paramQuadFace(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
00567
00568
00569 FieldContainer<double> paramFaceWeights;
00570 FieldContainer<double> paramFacePoints;
00571
00572
00573 FieldContainer<double> refFacePoints;
00574
00575
00576 FieldContainer<double> worksetFacePoints;
00577
00578
00579
00580
00581
00582
00583
00584
00585 subcellDim = 2;
00586 subcellOrd = 5;
00587
00588
00589
00590
00591
00592
00593 Teuchos::RCP<Cubature<double> > hexFaceCubature = cubFactory.create(paramQuadFace, 3);
00594
00595
00596 cubDim = hexFaceCubature -> getDimension();
00597 numCubPoints = hexFaceCubature -> getNumPoints();
00598
00599
00600 paramFacePoints.resize(numCubPoints, cubDim);
00601 paramFaceWeights.resize(numCubPoints);
00602 hexFaceCubature -> getCubature(paramFacePoints, paramFaceWeights);
00603
00604
00605
00606
00607
00608
00609 refFacePoints.resize(numCubPoints, pCellDim);
00610
00611
00612 worksetFacePoints.resize(worksetSize, numCubPoints, pCellDim);
00613
00614
00615 CellTools::mapToReferenceSubcell(refFacePoints,
00616 paramFacePoints,
00617 subcellDim,
00618 subcellOrd,
00619 hexahedron_8);
00620
00621
00622 CellTools::mapToPhysicalFrame(worksetFacePoints,
00623 refFacePoints,
00624 hexNodesAlt,
00625 hexahedron_8);
00626
00627
00628
00629
00630
00631 for(int pCell = 0; pCell < worksetSize; pCell++){
00632
00633 CellTools::printWorksetSubcell(hexNodesAlt, hexahedron_8, pCell, subcellDim, subcellOrd);
00634
00635 for(int pt = 0; pt < numCubPoints; pt++){
00636 std::cout << "\t 2D Gauss point ("
00637 << std::setw(8) << std::right << paramFacePoints(pt, 0) << ", "
00638 << std::setw(8) << std::right << paramFacePoints(pt, 1) << ") "
00639 << std::setw(8) << " --> " << "("
00640 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 0) << ", "
00641 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 1) << ", "
00642 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 2) << ")\n";
00643 }
00644 std::cout << "\n";
00645 }
00646
00647
00648 std::cout << "\n" \
00649 << "===============================================================================\n"\
00650 << "| EXAMPLE 8: Face normals at Gauss points on a Hexahedron face workset |\n"\
00651 << "===============================================================================\n";
00652
00653
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669 worksetJacobians.resize(worksetSize, numCubPoints, pCellDim, pCellDim);
00670
00671
00672 CellTools::setJacobian(worksetJacobians,
00673 refFacePoints,
00674 hexNodesAlt,
00675 hexahedron_8);
00676
00677
00678
00679
00680
00681
00682 FieldContainer<double> faceNormals(worksetSize, numCubPoints, pCellDim);
00683
00684
00685 CellTools::getPhysicalFaceNormals(faceNormals,
00686 worksetJacobians,
00687 subcellOrd,
00688 hexahedron_8);
00689
00690
00691 std::cout
00692 << "Face normals computed by CellTools::getPhysicalFaceNormals\n"
00693 << "Local face ordinal = " << subcellOrd <<"\n";
00694
00695 for(int pCell = 0; pCell < worksetSize; pCell++){
00696
00697 CellTools::printWorksetSubcell(hexNodesAlt, hexahedron_8, pCell, subcellDim, subcellOrd);
00698
00699 for(int pt = 0; pt < numCubPoints; pt++){
00700 std::cout << "\t 3D Gauss point: ("
00701 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 0) << ", "
00702 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 1) << ", "
00703 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 2) << ") "
00704 << std::setw(8) << " outer normal: " << "("
00705 << std::setw(8) << std::right << faceNormals(pCell, pt, 0) << ", "
00706 << std::setw(8) << std::right << faceNormals(pCell, pt, 1) << ", "
00707 << std::setw(8) << std::right << faceNormals(pCell, pt, 2) << ")\n";
00708 }
00709 std::cout << "\n";
00710 }
00711
00712
00713
00714
00715
00716
00717 FieldContainer<double> uFaceTan(worksetSize, numCubPoints, pCellDim);
00718 FieldContainer<double> vFaceTan(worksetSize, numCubPoints, pCellDim);
00719
00720
00721 CellTools::getPhysicalFaceTangents(uFaceTan,
00722 vFaceTan,
00723 worksetJacobians,
00724 subcellOrd,
00725 hexahedron_8);
00726
00727
00728 RealSpaceTools<double>::vecprod(faceNormals, uFaceTan, vFaceTan);
00729
00730
00731 std::cout << "Face normals computed by CellTools::getPhysicalFaceTangents followed by vecprod\n";
00732 for(int pCell = 0; pCell < worksetSize; pCell++){
00733
00734 CellTools::printWorksetSubcell(hexNodes, hexahedron_8, pCell, subcellDim, subcellOrd);
00735
00736 for(int pt = 0; pt < numCubPoints; pt++){
00737 std::cout << "\t 3D Gauss point: ("
00738 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 0) << ", "
00739 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 1) << ", "
00740 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 2) << ") "
00741 << std::setw(8) << " outer normal: " << "("
00742 << std::setw(8) << std::right << faceNormals(pCell, pt, 0) << ", "
00743 << std::setw(8) << std::right << faceNormals(pCell, pt, 1) << ", "
00744 << std::setw(8) << std::right << faceNormals(pCell, pt, 2) << ")\n";
00745 }
00746 std::cout << "\n";
00747 }
00748
00749
00750 std::cout << "\n" \
00751 << "===============================================================================\n"\
00752 << "| EXAMPLE 8: Side normals at Gauss points on a Hexahedron side workset |\n"\
00753 << "===============================================================================\n";
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764 sideNormals.resize(worksetSize, numCubPoints, pCellDim);
00765
00766
00767 CellTools::getPhysicalSideNormals(sideNormals,
00768 worksetJacobians,
00769 subcellOrd,
00770 hexahedron_8);
00771
00772
00773 std::cout << "Side normals computed by CellTools::getPhysicalSideNormals\n";
00774 for(int pCell = 0; pCell < worksetSize; pCell++){
00775
00776 CellTools::printWorksetSubcell(hexNodes, hexahedron_8, pCell, subcellDim, subcellOrd);
00777
00778 for(int pt = 0; pt < numCubPoints; pt++){
00779 std::cout << "\t 3D Gauss point: ("
00780 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 0) << ", "
00781 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 1) << ", "
00782 << std::setw(8) << std::right << worksetFacePoints(pCell, pt, 2) << ") "
00783 << std::setw(8) << " side normal: " << "("
00784 << std::setw(8) << std::right << sideNormals(pCell, pt, 0) << ", "
00785 << std::setw(8) << std::right << sideNormals(pCell, pt, 1) << ", "
00786 << std::setw(8) << std::right << sideNormals(pCell, pt, 2) << ")\n";
00787 }
00788 std::cout << "\n";
00789 }
00790
00791 return 0;
00792 }