95 #include "Epetra_Time.h" 
   96 #include "Epetra_Map.h" 
   97 #include "Epetra_FEVector.h" 
   98 #include "Epetra_FECrsMatrix.h" 
   99 #include "Epetra_SerialComm.h" 
  102 #include "Teuchos_oblackholestream.hpp" 
  103 #include "Teuchos_RCP.hpp" 
  108 #include "Shards_CellTopology.hpp" 
  111 #include "EpetraExt_MultiVectorOut.h" 
  114 using namespace Intrepid;
 
  116 int main(
int argc, 
char *argv[]) {
 
  120     std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
 
  121     std::cout <<
"Usage:\n\n";
 
  122     std::cout <<
"  ./Intrepid_example_Drivers_Example_10.exe deg NX NY NZ verbose\n\n";
 
  123     std::cout <<
" where \n";
 
  124     std::cout <<
"   int deg             - polynomial degree to be used (assumed >= 1) \n";
 
  125     std::cout <<
"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
 
  126     std::cout <<
"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  127     std::cout <<
"   int NZ              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  128     std::cout <<
"   verbose (optional)  - any character, indicates verbose output \n\n";
 
  134   int iprint     = argc - 1;
 
  135   Teuchos::RCP<std::ostream> outStream;
 
  136   Teuchos::oblackholestream bhs; 
 
  138     outStream = Teuchos::rcp(&std::cout, 
false);
 
  140     outStream = Teuchos::rcp(&bhs, 
false);
 
  143   Teuchos::oblackholestream oldFormatState;
 
  144   oldFormatState.copyfmt(std::cout);
 
  147     << 
"===============================================================================\n" \
 
  149     << 
"|  Example: Build Stiffness Matrix for                                        |\n" \
 
  150     << 
"|                   Poisson Equation on Hexahedral Mesh                       |\n" \
 
  152     << 
"|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
 
  153     << 
"|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
 
  154     << 
"|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
 
  156     << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  157     << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  159     << 
"===============================================================================\n";
 
  164   int deg          = atoi(argv[1]);  
 
  165   int NX           = atoi(argv[2]);  
 
  166   int NY           = atoi(argv[3]);  
 
  167   int NZ           = atoi(argv[4]);  
 
  173   typedef shards::CellTopology    CellTopology;
 
  174   CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
 
  177   int numNodesPerElem = hex_8.getNodeCount();
 
  178   int spaceDim = hex_8.getDimension();
 
  182   *outStream << 
"Generating mesh ... \n\n";
 
  184   *outStream << 
"   NX" << 
"   NY" << 
"   NZ\n";
 
  185   *outStream << std::setw(5) << NX <<
 
  186     std::setw(5) << NY << std::setw(5) << NZ << 
"\n\n";
 
  189   int numElems = NX*NY*NZ;
 
  190   int numNodes = (NX+1)*(NY+1)*(NZ+1);
 
  191   *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  192   *outStream << 
"    Number of Nodes: " << numNodes << 
" \n\n";
 
  195   double leftX = 0.0, rightX = 1.0;
 
  196   double leftY = 0.0, rightY = 1.0;
 
  197   double leftZ = 0.0, rightZ = 1.0;
 
  200   double hx = (rightX-leftX)/((
double)NX);
 
  201   double hy = (rightY-leftY)/((
double)NY);
 
  202   double hz = (rightZ-leftZ)/((
double)NZ);
 
  208   for (
int k=0; k<NZ+1; k++) 
 
  210       for (
int j=0; j<NY+1; j++) 
 
  212           for (
int i=0; i<NX+1; i++) 
 
  214               nodeCoord(inode,0) = leftX + (double)i*hx;
 
  215               nodeCoord(inode,1) = leftY + (double)j*hy;
 
  216               nodeCoord(inode,2) = leftZ + (double)k*hz;
 
  217               if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
 
  219                   nodeOnBoundary(inode)=1;
 
  223                   nodeOnBoundary(inode)=0;
 
  232   ofstream fcoordout(
"coords.dat");
 
  233   for (
int i=0; i<numNodes; i++) {
 
  234     fcoordout << nodeCoord(i,0) <<
" ";
 
  235     fcoordout << nodeCoord(i,1) <<
" ";
 
  236     fcoordout << nodeCoord(i,2) <<
"\n";
 
  246   for (
int k=0; k<NZ; k++) 
 
  248       for (
int j=0; j<NY; j++) 
 
  250           for (
int i=0; i<NX; i++) 
 
  252               elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  253               elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  254               elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  255               elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  256               elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  257               elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  258               elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  259               elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  266   ofstream fe2nout(
"elem2node.dat");
 
  267   for (
int k=0;k<NZ;k++)
 
  269       for (
int j=0; j<NY; j++) 
 
  271           for (
int i=0; i<NX; i++) 
 
  273               int ielem = i + j * NX + k * NY * NY;
 
  274               for (
int m=0; m<numNodesPerElem; m++)
 
  276                   fe2nout << elemToNode(ielem,m) <<
"  ";
 
  286   *outStream << 
"Getting cubature ... \n\n";
 
  290   int cubDegree = 2*deg;
 
  291   Teuchos::RCP<Cubature<double> > quadCub = cubFactory.
create(hex_8, cubDegree); 
 
  293   int cubDim       = quadCub->getDimension();
 
  294   int numCubPoints = quadCub->getNumPoints();
 
  299   quadCub->getCubature(cubPoints, cubWeights);
 
  304   *outStream << 
"Getting basis ... \n\n";
 
  308   int numFieldsG = quadHGradBasis.getCardinality();
 
  313   quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
 
  314   quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
 
  318   const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
 
  320   for (
int k=0;k<NZ;k++) 
 
  322       for (
int j=0;j<NY;j++) 
 
  324           for (
int i=0;i<NX;i++) 
 
  326               const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
 
  329               for (
int kloc=0;kloc<=deg;kloc++) 
 
  331                   for (
int jloc=0;jloc<=deg;jloc++) 
 
  333                       for (
int iloc=0;iloc<=deg;iloc++)
 
  335                           ltgMapping(ielem,local_dof_cur) = start 
 
  336                             + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
 
  337                             + jloc * ( NX * deg + 1 )
 
  350   ofstream ltgout(
"ltg.dat");
 
  351   for (
int k=0;k<NZ;k++)  
 
  353       for (
int j=0; j<NY; j++) 
 
  355           for (
int i=0; i<NX; i++) 
 
  357               int ielem = i + j * NX + k * NX * NY;
 
  358               for (
int m=0; m<numFieldsG; m++)
 
  360                   ltgout << ltgMapping(ielem,m) <<
"  ";
 
  370   Epetra_SerialComm Comm;
 
  371   Epetra_Map globalMapG(numDOF, 0, Comm);
 
  372   Epetra_FEVector u(globalMapG);  u.Random();
 
  373   Epetra_FEVector Ku(globalMapG);
 
  376   Epetra_Time instantiateTimer(Comm);
 
  377   Epetra_FECrsMatrix StiffMatrix(Copy,globalMapG,8*numFieldsG); 
 
  378   const double instantiateTime = instantiateTimer.ElapsedTime();
 
  382   *outStream << 
"Building local stiffness matrices...\n\n";
 
  385   int numCells = numElems; 
 
  403   for (
int i=0;i<numElems;i++)
 
  405       for (
int j=0;j<numNodesPerElem;j++)
 
  407           const int nodeCur = elemToNode(i,j);
 
  408           for (
int k=0;k<spaceDim;k++) 
 
  410               cellVertices(i,j,k) = nodeCoord(nodeCur,k);
 
  415   Epetra_Time localConstructTimer( Comm );
 
  418   CellTools::setJacobian(cellJacobian,cubPoints,cellVertices,hex_8);
 
  419   CellTools::setJacobianInv(cellJacobInv, cellJacobian );
 
  420   CellTools::setJacobianDet(cellJacobDet, cellJacobian );
 
  423   fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
 
  426   fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
 
  429   fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
 
  430                                weightedMeasure, transformedBasisGradients);
 
  433   fst::integrate<double>(localStiffMatrices,
 
  434                          transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
 
  436   const double localConstructTime = localConstructTimer.ElapsedTime();
 
  439   Epetra_Time insertionTimer(Comm);
 
  442   for (
int k=0; k<numElems; k++) 
 
  445       StiffMatrix.InsertGlobalValues(numFieldsG,<gMapping(k,0),numFieldsG,<gMapping(k,0),&localStiffMatrices(k,0,0));
 
  448   StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
 
  449   const double insertionTime = insertionTimer.ElapsedTime( );
 
  451   *outStream << 
"Time to instantiate global stiffness matrix: " << instantiateTime << 
"\n";
 
  452   *outStream << 
"Time to build local matrices (including Jacobian computation): "<< localConstructTime << 
"\n";
 
  453   *outStream << 
"Time to assemble global matrix from local matrices: " << insertionTime << 
"\n";
 
  454   *outStream << 
"Total construction time: " << instantiateTime + localConstructTime + insertionTime << 
"\n";
 
  456   Epetra_Time applyTimer(Comm);
 
  457   StiffMatrix.Apply(u,Ku);
 
  458   const double multTime = applyTimer.ElapsedTime();
 
  459   *outStream << 
"Time to multiply onto a vector: " << multTime << 
"\n";
 
  461   *outStream << 
"End Result: TEST PASSED\n";
 
  464   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.