Intrepid
http://trilinos.sandia.gov/packages/docs/r10.12/packages/intrepid/example/Drivers/example_15.cpp
Go to the documentation of this file.
00001 // @HEADER
00002 // ************************************************************************
00003 //
00004 //                           Intrepid Package
00005 //                 Copyright (2007) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Pavel Bochev  (pbboche@sandia.gov)
00038 //                    Denis Ridzal  (dridzal@sandia.gov), or
00039 //                    Kara Peterson (kjpeter@sandia.gov)
00040 //
00041 // ************************************************************************
00042 // @HEADER
00043 
00082 // Intrepid includes
00083 #include "Intrepid_FunctionSpaceTools.hpp"
00084 #include "Intrepid_FieldContainer.hpp"
00085 #include "Intrepid_CellTools.hpp"
00086 //#include "Intrepid_ArrayTools.hpp"
00087 #include "Intrepid_HGRAD_HEX_Cn_FEM.hpp"
00088 //#include "Intrepid_RealSpaceTools.hpp"
00089 #include "Intrepid_DefaultCubatureFactory.hpp"
00090 #include "Intrepid_Utils.hpp"
00091 
00092 // Epetra includes
00093 #include "Epetra_Time.h"
00094 #include "Epetra_Map.h"
00095 #include "Epetra_FEVector.h"
00096 #include "Epetra_FECrsMatrix.h"
00097 #include "Epetra_SerialComm.h"
00098 
00099 // Teuchos includes
00100 #include "Teuchos_oblackholestream.hpp"
00101 #include "Teuchos_RCP.hpp"
00102 //#include "Teuchos_BLAS.hpp"
00103 //#include "Teuchos_BLAS_types.hpp"
00104 
00105 // Shards includes
00106 #include "Shards_CellTopology.hpp"
00107 
00108 // EpetraExt includes
00109 #include "EpetraExt_MultiVectorOut.h"
00110 
00111 #include <vector>
00112 #include <map>
00113 
00114 using namespace std;
00115 using namespace Intrepid;
00116 
00117 int main(int argc, char *argv[]) {
00118 
00119   //Check number of arguments
00120   if (argc < 4) {
00121     std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n";
00122     std::cout <<"Usage:\n\n";
00123     std::cout <<"  ./Intrepid_example_Drivers_Example_15.exe deg NX NY NZ verbose\n\n";
00124     std::cout <<" where \n";
00125     std::cout <<"   int deg             - polynomial degree to be used (assumed >= 1) \n";
00126     std::cout <<"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
00127     std::cout <<"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
00128     std::cout <<"   int NZ              - num intervals in y direction (assumed box domain, 0,1) \n";
00129     std::cout <<"   verbose (optional)  - any character, indicates verbose output \n\n";
00130     exit(1);
00131   }
00132   
00133   // This little trick lets us print to std::cout only if
00134   // a (dummy) command-line argument is provided.
00135   int iprint     = argc - 1;
00136   Teuchos::RCP<std::ostream> outStream;
00137   Teuchos::oblackholestream bhs; // outputs nothing
00138   if (iprint > 2)
00139     outStream = Teuchos::rcp(&std::cout, false);
00140   else
00141     outStream = Teuchos::rcp(&bhs, false);
00142   
00143   // Save the format state of the original std::cout.
00144   Teuchos::oblackholestream oldFormatState;
00145   oldFormatState.copyfmt(std::cout);
00146   
00147   *outStream                                                            \
00148     << "===============================================================================\n" \
00149     << "|                                                                             |\n" \
00150     << "|  Example: Build Stiffness Matrix for                                        |\n" \
00151     << "|                   Poisson Equation on Hexahedral Mesh                       |\n" \
00152     << "|                                                                             |\n" \
00153     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00154     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00155     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00156     << "|                                                                             |\n" \
00157     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00158     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00159     << "|                                                                             |\n" \
00160     << "===============================================================================\n";
00161 
00162   
00163   // ************************************ GET INPUTS **************************************
00164   
00165   int deg          = atoi(argv[1]);  // polynomial degree to use
00166   int NX           = atoi(argv[2]);  // num intervals in x direction (assumed box domain, 0,1)
00167   int NY           = atoi(argv[3]);  // num intervals in y direction (assumed box domain, 0,1)
00168   int NZ           = atoi(argv[4]);  // num intervals in y direction (assumed box domain, 0,1)
00169   
00170 
00171   // *********************************** CELL TOPOLOGY **********************************
00172   
00173   // Get cell topology for base hexahedron
00174   typedef shards::CellTopology    CellTopology;
00175   CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
00176   
00177   // Get dimensions 
00178   int numNodesPerElem = hex_8.getNodeCount();
00179   int spaceDim = hex_8.getDimension();
00180   
00181   // *********************************** GENERATE MESH ************************************
00182   
00183   *outStream << "Generating mesh ... \n\n";
00184   
00185   *outStream << "   NX" << "   NY" << "   NZ\n";
00186   *outStream << std::setw(5) << NX <<
00187     std::setw(5) << NY << std::setw(5) << NZ << "\n\n";
00188   
00189   // Print mesh information
00190   int numElems = NX*NY*NZ;
00191   int numNodes = (NX+1)*(NY+1)*(NZ+1);
00192   *outStream << " Number of Elements: " << numElems << " \n";
00193   *outStream << "    Number of Nodes: " << numNodes << " \n\n";
00194   
00195   // Cube
00196   double leftX = 0.0, rightX = 1.0;
00197   double leftY = 0.0, rightY = 1.0;
00198   double leftZ = 0.0, rightZ = 1.0;
00199 
00200   // Mesh spacing
00201   double hx = (rightX-leftX)/((double)NX);
00202   double hy = (rightY-leftY)/((double)NY);
00203   double hz = (rightZ-leftZ)/((double)NZ);
00204 
00205   // Get nodal coordinates
00206   FieldContainer<double> nodeCoord(numNodes, spaceDim);
00207   FieldContainer<int> nodeOnBoundary(numNodes);
00208   int inode = 0;
00209   for (int k=0; k<NZ+1; k++) 
00210     {
00211       for (int j=0; j<NY+1; j++) 
00212         {
00213           for (int i=0; i<NX+1; i++) 
00214             {
00215               nodeCoord(inode,0) = leftX + (double)i*hx;
00216               nodeCoord(inode,1) = leftY + (double)j*hy;
00217               nodeCoord(inode,2) = leftZ + (double)k*hz;
00218               if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
00219                 {
00220                   nodeOnBoundary(inode)=1;
00221                 }
00222               else 
00223                 {
00224                   nodeOnBoundary(inode)=0;
00225                 }
00226               inode++;
00227             }
00228         }
00229     }
00230 #define DUMP_DATA
00231 #ifdef DUMP_DATA
00232   // Print nodal coords
00233   ofstream fcoordout("coords.dat");
00234   for (int i=0; i<numNodes; i++) {
00235     fcoordout << nodeCoord(i,0) <<" ";
00236     fcoordout << nodeCoord(i,1) <<" ";
00237     fcoordout << nodeCoord(i,2) <<"\n";
00238   }
00239   fcoordout.close();
00240 #endif
00241   
00242   
00243   // Element to Node map
00244   // We'll keep it around, but this is only the DOFMap if you are in the lowest order case.
00245   FieldContainer<int> elemToNode(numElems, numNodesPerElem);
00246   int ielem = 0;
00247   for (int k=0; k<NZ; k++) 
00248     {
00249       for (int j=0; j<NY; j++) 
00250         {
00251           for (int i=0; i<NX; i++) 
00252             {
00253               elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
00254               elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
00255               elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
00256               elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
00257               elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
00258               elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
00259               elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
00260               elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
00261               ielem++;
00262             }
00263         }
00264     }
00265 #ifdef DUMP_DATA
00266   // Output connectivity
00267   ofstream fe2nout("elem2node.dat");
00268   for (int k=0;k<NZ;k++)
00269     {
00270       for (int j=0; j<NY; j++) 
00271         {
00272           for (int i=0; i<NX; i++) 
00273             {
00274               int ielem = i + j * NX + k * NY * NY;
00275               for (int m=0; m<numNodesPerElem; m++)
00276                 {
00277                   fe2nout << elemToNode(ielem,m) <<"  ";
00278                 }
00279               fe2nout <<"\n";
00280             }
00281         }
00282     }
00283   fe2nout.close();
00284 #endif
00285   
00286   // ************************************ CUBATURE ************************************** 
00287   *outStream << "Getting cubature ... \n\n";
00288   
00289   // Get numerical integration points and weights
00290   DefaultCubatureFactory<double>  cubFactory;                                   
00291   int cubDegree = 2*deg;
00292   Teuchos::RCP<Cubature<double> > quadCub = cubFactory.create(hex_8, cubDegree); 
00293   
00294   int cubDim       = quadCub->getDimension();
00295   int numCubPoints = quadCub->getNumPoints();
00296   
00297   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00298   FieldContainer<double> cubWeights(numCubPoints);
00299   
00300   quadCub->getCubature(cubPoints, cubWeights);
00301   
00302 
00303   // ************************************** BASIS ***************************************
00304   
00305   *outStream << "Getting basis ... \n\n";
00306   
00307   // Define basis 
00308   Basis_HGRAD_HEX_Cn_FEM<double, FieldContainer<double> > quadHGradBasis(deg,POINTTYPE_SPECTRAL);
00309   int numFieldsG = quadHGradBasis.getCardinality();
00310   FieldContainer<double> quadGVals(numFieldsG, numCubPoints); 
00311   FieldContainer<double> quadGrads(numFieldsG, numCubPoints, spaceDim); 
00312   
00313   // Evaluate basis values and gradients at cubature points
00314   quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
00315   quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
00316 
00317   // create the local-global mapping
00318   FieldContainer<int> ltgMapping(numElems,numFieldsG);
00319   const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
00320   ielem=0;
00321   for (int k=0;k<NZ;k++) 
00322     {
00323       for (int j=0;j<NY;j++) 
00324         {
00325           for (int i=0;i<NX;i++) 
00326             {
00327               const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
00328               // loop over local dof on this cell
00329               int local_dof_cur=0;
00330               for (int kloc=0;kloc<=deg;kloc++) 
00331                 {
00332                   for (int jloc=0;jloc<=deg;jloc++) 
00333                     {
00334                       for (int iloc=0;iloc<=deg;iloc++)
00335                         {
00336                           ltgMapping(ielem,local_dof_cur) = start 
00337                             + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
00338                             + jloc * ( NX * deg + 1 )
00339                             + iloc;
00340                           local_dof_cur++;
00341                         }
00342                     }
00343                 }
00344               ielem++;
00345             }
00346         }
00347     }
00348 #ifdef DUMP_DATA
00349   // Output ltg mapping 
00350   ielem = 0;
00351   ofstream ltgout("ltg.dat");
00352   for (int k=0;k<NZ;k++)  
00353     {
00354       for (int j=0; j<NY; j++) 
00355         {
00356           for (int i=0; i<NX; i++) 
00357             {
00358               int ielem = i + j * NX + k * NX * NY;
00359               for (int m=0; m<numFieldsG; m++)
00360                 {
00361                   ltgout << ltgMapping(ielem,m) <<"  ";
00362                 }
00363               ltgout <<"\n";
00364             }
00365         }
00366     }
00367   ltgout.close();
00368 #endif
00369 
00370   // ********** DECLARE GLOBAL OBJECTS *************
00371   Epetra_SerialComm Comm;
00372   Epetra_Map globalMapG(numDOF, 0, Comm);
00373   Epetra_FEVector u(globalMapG);  u.Random();
00374   Epetra_FEVector Ku(globalMapG);
00375 
00376   // time the instantiation 
00377 //   Epetra_Time instantiateTimer(Comm);
00378 //   Epetra_FECrsMatrix StiffMatrix(Copy,globalMapG,8*numFieldsG); 
00379 //   const double instantiateTime = instantiateTimer.ElapsedTime();
00380 
00381 
00382   // ********** CONSTRUCT AND INSERT LOCAL STIFFNESS MATRICES ***********
00383   *outStream << "Building local stiffness matrices...\n\n";
00384   typedef CellTools<double>  CellTools;
00385   typedef FunctionSpaceTools fst;
00386   int numCells = numElems; 
00387 
00388   // vertices
00389   FieldContainer<double> cellVertices(numCells,numNodesPerElem,spaceDim);
00390 
00391   // jacobian information
00392   FieldContainer<double> cellJacobian(numCells,numCubPoints,spaceDim,spaceDim);
00393   FieldContainer<double> cellJacobInv(numCells,numCubPoints,spaceDim,spaceDim);
00394   FieldContainer<double> cellJacobDet(numCells,numCubPoints);
00395 
00396   // element stiffness matrices and supporting storage space
00397   FieldContainer<double> localStiffMatrices(numCells, numFieldsG, numFieldsG);
00398   FieldContainer<double> transformedBasisGradients(numCells,numFieldsG,numCubPoints,spaceDim);
00399   FieldContainer<double> weightedTransformedBasisGradients(numCells,numFieldsG,numCubPoints,spaceDim);
00400   FieldContainer<double> weightedMeasure(numCells, numCubPoints);
00401 
00402 
00403   // get vertices of cells (for computing Jacobians)
00404   for (int i=0;i<numElems;i++)
00405     {
00406       for (int j=0;j<numNodesPerElem;j++)
00407         {
00408           const int nodeCur = elemToNode(i,j);
00409           for (int k=0;k<spaceDim;k++) 
00410             {
00411               cellVertices(i,j,k) = nodeCoord(nodeCur,k);
00412             }
00413         }
00414     }
00415    
00416   Epetra_Time localConstructTimer( Comm );
00417 
00418   // jacobian evaluation 
00419   CellTools::setJacobian(cellJacobian,cubPoints,cellVertices,hex_8);
00420   CellTools::setJacobianInv(cellJacobInv, cellJacobian );
00421   CellTools::setJacobianDet(cellJacobDet, cellJacobian );
00422 
00423   // transform reference element gradients to each cell
00424   fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
00425       
00426   // compute weighted measure
00427   fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
00428 
00429   // multiply values with weighted measure
00430   fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
00431                                weightedMeasure, transformedBasisGradients);
00432 
00433   // integrate to compute element stiffness matrix
00434   fst::integrate<double>(localStiffMatrices,
00435                          transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
00436 
00437   const double localConstructTime = localConstructTimer.ElapsedTime();
00438 
00439 
00440   Epetra_Time insertionTimer(Comm);
00441 
00442   vector<map<int,double> > mat(numDOF);
00443 
00444 
00445 
00446   // *** Element loop ***
00447   for (int el=0; el<numElems; el++) 
00448     {
00449       for (int i=0;i<numFieldsG;i++) // local rows
00450         {
00451           const int glob_row = ltgMapping(el,i);
00452           map<int,double> & cur_row = mat[glob_row];
00453 
00454           for (int j=0;j<numFieldsG;j++) // local columns
00455             {
00456               const int glob_col = ltgMapping(el,j);
00457               const double cur_val = localStiffMatrices(el,i,j);
00458               map<int,double>::iterator it = cur_row.find( glob_col );
00459               if (it != cur_row.end()) // current column already in row
00460                 {
00461                   it->second += cur_val;
00462                 }
00463               else
00464                 {
00465                   cur_row[glob_col] = cur_val;
00466                 }
00467             }
00468         }
00469     }
00470   //StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
00471   const double insertionTime = insertionTimer.ElapsedTime( );
00472   
00473 //   *outStream << "Time to instantiate global stiffness matrix: " << instantiateTime << "\n";
00474   *outStream << "Time to build local matrices (including Jacobian computation): "<< localConstructTime << "\n";
00475   *outStream << "Time to assemble global matrix from local matrices: " << insertionTime << "\n";
00476   *outStream << "Total construction time: " << localConstructTime + insertionTime << "\n";
00477 
00478 //   Epetra_Time applyTimer(Comm);
00479 //   StiffMatrix.Apply(u,Ku);
00480 //   const double multTime = applyTimer.ElapsedTime();
00481 //   *outStream << "Time to multiply onto a vector: " << multTime << "\n";
00482 
00483   *outStream << "End Result: TEST PASSED\n";
00484   
00485   // reset format state of std::cout
00486   std::cout.copyfmt(oldFormatState);
00487   
00488   return 0;
00489 }
00490