93 #include "Epetra_Time.h" 
   94 #include "Epetra_Map.h" 
   95 #include "Epetra_FEVector.h" 
   96 #include "Epetra_SerialComm.h" 
   99 #include "Teuchos_oblackholestream.hpp" 
  100 #include "Teuchos_RCP.hpp" 
  105 #include "Shards_CellTopology.hpp" 
  108 #include "EpetraExt_MultiVectorOut.h" 
  111 using namespace Intrepid;
 
  113 int main(
int argc, 
char *argv[]) {
 
  117     std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
 
  118     std::cout <<
"Usage:\n\n";
 
  119     std::cout <<
"  ./Intrepid_example_Drivers_Example_10.exe deg NX NY NZ verbose\n\n";
 
  120     std::cout <<
" where \n";
 
  121     std::cout <<
"   int deg             - polynomial degree to be used (assumed >= 1) \n";
 
  122     std::cout <<
"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
 
  123     std::cout <<
"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  124     std::cout <<
"   int NZ              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  125     std::cout <<
"   verbose (optional)  - any character, indicates verbose output \n\n";
 
  131   int iprint     = argc - 1;
 
  132   Teuchos::RCP<std::ostream> outStream;
 
  133   Teuchos::oblackholestream bhs; 
 
  135     outStream = Teuchos::rcp(&std::cout, 
false);
 
  137     outStream = Teuchos::rcp(&bhs, 
false);
 
  140   Teuchos::oblackholestream oldFormatState;
 
  141   oldFormatState.copyfmt(std::cout);
 
  144     << 
"===============================================================================\n" \
 
  146     << 
"|  Example: Build Stiffness Matrix for                                        |\n" \
 
  147     << 
"|                   Poisson Equation on Hexahedral Mesh                       |\n" \
 
  149     << 
"|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
 
  150     << 
"|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
 
  151     << 
"|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
 
  153     << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  154     << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  156     << 
"===============================================================================\n";
 
  161   int deg          = atoi(argv[1]);  
 
  162   int NX           = atoi(argv[2]);  
 
  163   int NY           = atoi(argv[3]);  
 
  164   int NZ           = atoi(argv[4]);  
 
  170   typedef shards::CellTopology    CellTopology;
 
  171   CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
 
  174   int numNodesPerElem = hex_8.getNodeCount();
 
  175   int spaceDim = hex_8.getDimension();
 
  179   *outStream << 
"Generating mesh ... \n\n";
 
  181   *outStream << 
"   NX" << 
"   NY" << 
"   NZ\n";
 
  182   *outStream << std::setw(5) << NX <<
 
  183     std::setw(5) << NY << std::setw(5) << NZ << 
"\n\n";
 
  186   int numElems = NX*NY*NZ;
 
  187   int numNodes = (NX+1)*(NY+1)*(NZ+1);
 
  188   *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  189   *outStream << 
"    Number of Nodes: " << numNodes << 
" \n\n";
 
  192   double leftX = 0.0, rightX = 1.0;
 
  193   double leftY = 0.0, rightY = 1.0;
 
  194   double leftZ = 0.0, rightZ = 1.0;
 
  197   double hx = (rightX-leftX)/((
double)NX);
 
  198   double hy = (rightY-leftY)/((
double)NY);
 
  199   double hz = (rightZ-leftZ)/((
double)NZ);
 
  205   for (
int k=0; k<NZ+1; k++) 
 
  207       for (
int j=0; j<NY+1; j++) 
 
  209           for (
int i=0; i<NX+1; i++) 
 
  211               nodeCoord(inode,0) = leftX + (double)i*hx;
 
  212               nodeCoord(inode,1) = leftY + (double)j*hy;
 
  213               nodeCoord(inode,2) = leftZ + (double)k*hz;
 
  214               if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
 
  216                   nodeOnBoundary(inode)=1;
 
  220                   nodeOnBoundary(inode)=0;
 
  229   ofstream fcoordout(
"coords.dat");
 
  230   for (
int i=0; i<numNodes; i++) {
 
  231     fcoordout << nodeCoord(i,0) <<
" ";
 
  232     fcoordout << nodeCoord(i,1) <<
" ";
 
  233     fcoordout << nodeCoord(i,2) <<
"\n";
 
  243   for (
int k=0; k<NZ; k++) 
 
  245       for (
int j=0; j<NY; j++) 
 
  247           for (
int i=0; i<NX; i++) 
 
  249               elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  250               elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  251               elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  252               elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  253               elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  254               elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  255               elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  256               elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  263   ofstream fe2nout(
"elem2node.dat");
 
  264   for (
int k=0;k<NZ;k++)
 
  266       for (
int j=0; j<NY; j++) 
 
  268           for (
int i=0; i<NX; i++) 
 
  270               int ielem = i + j * NX + k * NY * NY;
 
  271               for (
int m=0; m<numNodesPerElem; m++)
 
  273                   fe2nout << elemToNode(ielem,m) <<
"  ";
 
  283   *outStream << 
"Getting cubature ... \n\n";
 
  287   int cubDegree = 2*deg;
 
  288   Teuchos::RCP<Cubature<double> > quadCub = cubFactory.
create(hex_8, cubDegree); 
 
  290   int cubDim       = quadCub->getDimension();
 
  291   int numCubPoints = quadCub->getNumPoints();
 
  296   quadCub->getCubature(cubPoints, cubWeights);
 
  300   *outStream << 
"Getting basis ... \n\n";
 
  304   int numFieldsG = quadHGradBasis.getCardinality();
 
  309   quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
 
  310   quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
 
  314   const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
 
  316   for (
int k=0;k<NZ;k++) 
 
  318       for (
int j=0;j<NY;j++) 
 
  320           for (
int i=0;i<NX;i++) 
 
  322               const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
 
  325               for (
int kloc=0;kloc<=deg;kloc++) 
 
  327                   for (
int jloc=0;jloc<=deg;jloc++) 
 
  329                       for (
int iloc=0;iloc<=deg;iloc++)
 
  331                           ltgMapping(ielem,local_dof_cur) = start 
 
  332                             + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
 
  333                             + jloc * ( NX * deg + 1 )
 
  346   ofstream ltgout(
"ltg.dat");
 
  347   for (
int k=0;k<NZ;k++)  
 
  349       for (
int j=0; j<NY; j++) 
 
  351           for (
int i=0; i<NX; i++) 
 
  353               int ielem = i + j * NX + k * NX * NY;
 
  354               for (
int m=0; m<numFieldsG; m++)
 
  356                   ltgout << ltgMapping(ielem,m) <<
"  ";
 
  366   Epetra_SerialComm Comm;
 
  367   Epetra_Map globalMapG(numDOF, 0, Comm);
 
  368   Epetra_FEVector u(globalMapG);  u.Random();
 
  369   Epetra_FEVector Ku(globalMapG);
 
  374   *outStream << 
"Building reference stiffness matrix...\n\n";
 
  390   Epetra_Time localConstructTimer( Comm );
 
  391   refCellNodes(0,0,0) = 0.0;  refCellNodes(0,0,1) = 0.0;  refCellNodes(0,0,2) = 0.0;
 
  392   refCellNodes(0,1,0) = hx;   refCellNodes(0,1,1) = 0.0;  refCellNodes(0,1,2) = 0.0;
 
  393   refCellNodes(0,2,0) = hx;   refCellNodes(0,2,1) = hy;   refCellNodes(0,2,2) = 0.0;
 
  394   refCellNodes(0,3,0) = 0.0;  refCellNodes(0,3,1) = hy;   refCellNodes(0,3,2) = 0.0;
 
  395   refCellNodes(0,4,0) = 0.0;  refCellNodes(0,4,1) = 0.0;  refCellNodes(0,4,2) = hz;
 
  396   refCellNodes(0,5,0) = hx;   refCellNodes(0,5,1) = 0.0;  refCellNodes(0,5,2) = hz;
 
  397   refCellNodes(0,6,0) = hx;   refCellNodes(0,6,1) = hy;   refCellNodes(0,6,2) = hz;
 
  398   refCellNodes(0,7,0) = 0.0;  refCellNodes(0,7,1) = hy;   refCellNodes(0,7,2) = hz;
 
  401   CellTools::setJacobian(cellJacobian,cubPoints,refCellNodes,hex_8);
 
  402   CellTools::setJacobianInv(cellJacobInv, cellJacobian );
 
  403   CellTools::setJacobianDet(cellJacobDet, cellJacobian );
 
  406   fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
 
  409   fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
 
  412   fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
 
  413                                weightedMeasure, transformedBasisGradients);
 
  416   fst::integrate<double>(localStiffMatrix,
 
  417                          transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
 
  419   const double localConstructTime = localConstructTimer.ElapsedTime();
 
  427   Epetra_Time multTimer(Comm);
 
  428   Teuchos::BLAS<int,double> blas;
 
  432   double *uVals = u[0];
 
  433   double *KuVals = Ku[0];
 
  435   Epetra_Time scatterTimer(Comm);
 
  436   std::cout << 
"Scattering\n";
 
  438   for (
int k=0; k<numElems; k++) 
 
  440       for (
int i=0;i<numFieldsG;i++) 
 
  442           uScattered(k,i) = uVals[ltgMapping(k,i)];
 
  445   const double scatterTime = scatterTimer.ElapsedTime();
 
  447   Epetra_Time blasTimer(Comm);
 
  448   blas.GEMM(Teuchos::NO_TRANS , Teuchos::NO_TRANS , 
 
  449             numFieldsG , numElems, numFieldsG  , 
 
  451             &localStiffMatrix(0,0,0) , 
 
  458   const double blasTime = blasTimer.ElapsedTime();
 
  460   Epetra_Time gatherTimer(Comm);
 
  462   for (
int k=0;k<numElems;k++)
 
  464       for (
int i=0;i<numFieldsG;i++)
 
  466           KuVals[ltgMapping(k,i)] += KuScattered(k,i);
 
  470   const double gatherTime = gatherTimer.ElapsedTime();
 
  473   *outStream << 
"Time to build local matrix (including Jacobian computation): "<< localConstructTime << 
"\n";
 
  474   *outStream << 
"Time to scatter " << scatterTime << 
"\n";
 
  475   *outStream << 
"Time for local application " << blasTime << 
"\n";
 
  476   *outStream << 
"Time to gather " << gatherTime << 
"\n";
 
  477   *outStream << 
"Total matrix-free time " << scatterTime + blasTime + gatherTime << 
"\n";
 
  480   *outStream << 
"End Result: TEST PASSED\n";
 
  483   std::cout.copyfmt(oldFormatState);
 
Header file for utility class to provide multidimensional containers. 
Header file for the Intrepid::HGRAD_HEX_Cn_FEM class. 
Header file for the abstract base class Intrepid::DefaultCubatureFactory. 
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron 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.