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" 
  101 #include "Teuchos_Time.hpp" 
  102 #include "Teuchos_SerialDenseMatrix.hpp" 
  105 #include "Shards_CellTopology.hpp" 
  108 #include "EpetraExt_MultiVectorOut.h" 
  111 using namespace Intrepid;
 
  114 int main(
int argc, 
char *argv[]) {
 
  118     std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
 
  119     std::cout <<
"Usage:\n\n";
 
  120     std::cout <<
"  ./Intrepid_example_Drivers_Example_14.exe deg NX NY NZ verbose\n\n";
 
  121     std::cout <<
" where \n";
 
  122     std::cout <<
"   int deg             - polynomial degree to be used (assumed >= 1) \n";
 
  123     std::cout <<
"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
 
  124     std::cout <<
"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  125     std::cout <<
"   int NZ              - num intervals in y direction (assumed box domain, 0,1) \n";
 
  126     std::cout <<
"   verbose (optional)  - any character, indicates verbose output \n\n";
 
  132   int iprint     = argc - 1;
 
  133   Teuchos::RCP<std::ostream> outStream;
 
  134   Teuchos::oblackholestream bhs; 
 
  136     outStream = Teuchos::rcp(&std::cout, 
false);
 
  138     outStream = Teuchos::rcp(&bhs, 
false);
 
  141   Teuchos::oblackholestream oldFormatState;
 
  142   oldFormatState.copyfmt(std::cout);
 
  145     << 
"===============================================================================\n" \
 
  147     << 
"|  Example: Apply Stiffness Matrix for                                        |\n" \
 
  148     << 
"|                   Poisson Operator on Hexahedral Mesh with BLAS             |\n" \
 
  150     << 
"|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
 
  151     << 
"|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
 
  152     << 
"|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
 
  154     << 
"|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
 
  155     << 
"|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
 
  157     << 
"===============================================================================\n";
 
  162   int deg          = atoi(argv[1]);  
 
  163   int NX           = atoi(argv[2]);  
 
  164   int NY           = atoi(argv[3]);  
 
  165   int NZ           = atoi(argv[4]);  
 
  171   typedef shards::CellTopology    CellTopology;
 
  172   CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
 
  175   int numNodesPerElem = hex_8.getNodeCount();
 
  176   int spaceDim = hex_8.getDimension();
 
  180   *outStream << 
"Generating mesh ... \n\n";
 
  182   *outStream << 
"   NX" << 
"   NY" << 
"   NZ\n";
 
  183   *outStream << std::setw(5) << NX <<
 
  184     std::setw(5) << NY << std::setw(5) << NZ << 
"\n\n";
 
  187   int numElems = NX*NY*NZ;
 
  188   int numNodes = (NX+1)*(NY+1)*(NZ+1);
 
  189   *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  190   *outStream << 
"    Number of Nodes: " << numNodes << 
" \n\n";
 
  193   double leftX = 0.0, rightX = 1.0;
 
  194   double leftY = 0.0, rightY = 1.0;
 
  195   double leftZ = 0.0, rightZ = 1.0;
 
  198   double hx = (rightX-leftX)/((
double)NX);
 
  199   double hy = (rightY-leftY)/((
double)NY);
 
  200   double hz = (rightZ-leftZ)/((
double)NZ);
 
  206   for (
int k=0; k<NZ+1; k++) 
 
  208       for (
int j=0; j<NY+1; j++) 
 
  210           for (
int i=0; i<NX+1; i++) 
 
  212               nodeCoord(inode,0) = leftX + (double)i*hx;
 
  213               nodeCoord(inode,1) = leftY + (double)j*hy;
 
  214               nodeCoord(inode,2) = leftZ + (double)k*hz;
 
  215               if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
 
  217                   nodeOnBoundary(inode)=1;
 
  221                   nodeOnBoundary(inode)=0;
 
  230   ofstream fcoordout(
"coords.dat");
 
  231   for (
int i=0; i<numNodes; i++) {
 
  232     fcoordout << nodeCoord(i,0) <<
" ";
 
  233     fcoordout << nodeCoord(i,1) <<
" ";
 
  234     fcoordout << nodeCoord(i,2) <<
"\n";
 
  242   *outStream << 
"Getting cubature and basis ... \n\n";
 
  249   const int numCubPoints1D = glcub->getNumPoints();
 
  254   glcub->getCubature(cubPoints1D,cubWeights1D);
 
  258   cub_to_tensor.push_back( glcub );
 
  259   cub_to_tensor.push_back( glcub );
 
  260   cub_to_tensor.push_back( glcub );
 
  265   const int numFields = hexBasis.getCardinality();
 
  268   const int numCubPoints = cubhex.getNumPoints();
 
  272   cubhex.getCubature( cubPoints3D , cubWeights3D );
 
  273   hexBasis.getValues( basisGrads , cubPoints3D , OPERATOR_GRAD );
 
  278   for (
int k=0; k<NZ; k++) 
 
  280       for (
int j=0; j<NY; j++) 
 
  282           for (
int i=0; i<NX; i++) 
 
  284               elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  285               elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  286               elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  287               elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  288               elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
 
  289               elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
 
  290               elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
 
  291               elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
 
  298   ofstream fe2nout(
"elem2node.dat");
 
  299   for (
int k=0;k<NZ;k++)
 
  301       for (
int j=0; j<NY; j++) 
 
  303           for (
int i=0; i<NX; i++) 
 
  305               int ielem = i + j * NX + k * NY * NY;
 
  306               for (
int m=0; m<numNodesPerElem; m++)
 
  308                   fe2nout << elemToNode(ielem,m) <<
"  ";
 
  320   const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
 
  322   for (
int k=0;k<NZ;k++) 
 
  324       for (
int j=0;j<NY;j++) 
 
  326           for (
int i=0;i<NX;i++) 
 
  328               const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
 
  331               for (
int kloc=0;kloc<=deg;kloc++) 
 
  333                   for (
int jloc=0;jloc<=deg;jloc++) 
 
  335                       for (
int iloc=0;iloc<=deg;iloc++)
 
  337                           ltgMapping(ielem,local_dof_cur) = start 
 
  338                             + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
 
  339                             + jloc * ( NX * deg + 1 )
 
  352   ofstream ltgout(
"ltg.dat");
 
  353   for (
int k=0;k<NZ;k++)  
 
  355       for (
int j=0; j<NY; j++) 
 
  357           for (
int i=0; i<NX; i++) 
 
  359               int ielem = i + j * NX + k * NX * NY;
 
  360               for (
int m=0; m<numFields; m++)
 
  362                   ltgout << ltgMapping(ielem,m) <<
"  ";
 
  372   Epetra_SerialComm Comm;
 
  373   Epetra_Map globalMapG(numDOF, 0, Comm);
 
  374   Epetra_FEVector u(globalMapG);  u.Random();
 
  375   Epetra_FEVector Ku(globalMapG);
 
  385   for (
int i=0;i<numElems;i++)
 
  387       for (
int j=0;j<numNodesPerElem;j++)
 
  389           const int nodeCur = elemToNode(i,j);
 
  390           for (
int k=0;k<spaceDim;k++) 
 
  392               cellVertices(i,j,k) = nodeCoord(nodeCur,k);
 
  409   Teuchos::SerialDenseMatrix<int,double> refGradsAsMat( Teuchos::View ,
 
  411                                                      numCubPoints * spaceDim ,
 
  412                                                      numCubPoints * spaceDim ,
 
  414   Teuchos::SerialDenseMatrix<int,double> uScatAsMat( Teuchos::View ,
 
  419   Teuchos::SerialDenseMatrix<int,double> GraduScatAsMat( Teuchos::View ,
 
  421                                                          numCubPoints * spaceDim ,
 
  422                                                          numCubPoints * spaceDim ,
 
  424   Teuchos::SerialDenseMatrix<int,double> KuScatAsMat( Teuchos::View ,
 
  438   double *uVals = u[0];
 
  439   double *KuVals = Ku[0];
 
  441   Teuchos::Time full_timer( 
"Time to apply operator matrix-free:" );
 
  442   Teuchos::Time scatter_timer( 
"Time to scatter dof:" );
 
  443   Teuchos::Time elementwise_timer( 
"Time to do elementwise computation:" ); 
 
  444   Teuchos::Time grad_timer( 
"Time to compute gradients:" );
 
  445   Teuchos::Time pointwise_timer( 
"Time to do pointwise transformations:" );
 
  446   Teuchos::Time moment_timer( 
"Time to compute moments:" );
 
  447   Teuchos::Time gather_timer( 
"Time to gather dof:" );
 
  451   scatter_timer.start();
 
  452   for (
int k=0; k<numElems; k++) 
 
  454       for (
int i=0;i<numFields;i++) 
 
  456           uScattered(k,0,i) = uVals[ltgMapping(k,i)];
 
  459   scatter_timer.stop();
 
  461   elementwise_timer.start();
 
  465   GraduScatAsMat.multiply( Teuchos::NO_TRANS , Teuchos::NO_TRANS ,
 
  466                            1.0 , refGradsAsMat , uScatAsMat , 0.0 );
 
  470   pointwise_timer.start();
 
  472   Intrepid::FunctionSpaceToolsInPlace::HGRADtransformGRAD<double>( gradU , cellJacobian );
 
  473   Intrepid::FunctionSpaceToolsInPlace::HGRADtransformGRADDual<double>( gradU , cellJacobian );
 
  474   Intrepid::FunctionSpaceToolsInPlace::multiplyMeasure<double>( gradU , cellJacobDet );
 
  475   pointwise_timer.stop();
 
  476   moment_timer.start();
 
  478   KuScatAsMat.multiply( Teuchos::TRANS , Teuchos::NO_TRANS , 1.0 , refGradsAsMat , GraduScatAsMat , 0.0 );
 
  481   elementwise_timer.stop();
 
  482   gather_timer.start();
 
  483   for (
int k=0;k<numElems;k++)
 
  485       for (
int i=0;i<numFields;i++)
 
  487           KuVals[ltgMapping(k,i)] += KuScattered(k,0,i);
 
  493   *outStream << full_timer.name() << 
" " << full_timer.totalElapsedTime() << 
" sec\n";
 
  494   *outStream << 
"\t" << scatter_timer.name() << 
" " << scatter_timer.totalElapsedTime() << 
" sec\n";
 
  495   *outStream << 
"\t" << elementwise_timer.name() << 
" " << elementwise_timer.totalElapsedTime() << 
" sec\n";
 
  496   *outStream << 
"\t\t" << grad_timer.name() << 
" " << grad_timer.totalElapsedTime() << 
" sec\n";
 
  497   *outStream << 
"\t\t" << pointwise_timer.name() << 
" " << pointwise_timer.totalElapsedTime() << 
" sec\n";
 
  498   *outStream << 
"\t\t" << moment_timer.name() << 
" " << moment_timer.totalElapsedTime() << 
" sec\n";
 
  499   *outStream << 
"\t" << gather_timer.name() << 
" " << gather_timer.totalElapsedTime() << 
" sec\n";
 
  502   *outStream << 
"End Result: TEST PASSED\n";
 
  505   std::cout.copyfmt(oldFormatState);
 
Header file for utility class to provide multidimensional containers. 
Header file for the Intrepid::CubatureTensor class. 
Header file for the Intrepid::HGRAD_HEX_Cn_FEM class. 
Utilizes cubature (integration) rules contained in the library Polylib (Spencer Sherwin, Aeronautics, Imperial College London) within Intrepid. 
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Header file for the Intrepid::CubaturePolylib class. 
Defines tensor-product cubature (integration) rules in Intrepid.