Intrepid
/usr/src/RPM/BUILD/trilinos10-10.6.4/packages/intrepid/example/Drivers/example_05.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 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Pavel Bochev  (pbboche@sandia.gov),
00025 //                    Denis Ridzal  (dridzal@sandia.gov),
00026 //                    Kara Peterson (kjpeter@sandia.gov).
00027 //
00028 // ************************************************************************
00029 // @HEADER
00030 
00070 // Intrepid includes
00071 #include "Intrepid_FunctionSpaceTools.hpp"
00072 #include "Intrepid_FieldContainer.hpp"
00073 #include "Intrepid_CellTools.hpp"
00074 #include "Intrepid_ArrayTools.hpp"
00075 #include "Intrepid_HGRAD_QUAD_Cn_FEM.hpp"
00076 #include "Intrepid_RealSpaceTools.hpp"
00077 #include "Intrepid_DefaultCubatureFactory.hpp"
00078 #include "Intrepid_Utils.hpp"
00079 
00080 // Epetra includes
00081 #include "Epetra_Time.h"
00082 #include "Epetra_Map.h"
00083 #include "Epetra_FECrsMatrix.h"
00084 #include "Epetra_FEVector.h"
00085 #include "Epetra_SerialComm.h"
00086 
00087 // Teuchos includes
00088 #include "Teuchos_oblackholestream.hpp"
00089 #include "Teuchos_RCP.hpp"
00090 #include "Teuchos_BLAS.hpp"
00091 
00092 // Shards includes
00093 #include "Shards_CellTopology.hpp"
00094 
00095 // EpetraExt includes
00096 #include "EpetraExt_RowMatrixOut.h"
00097 #include "EpetraExt_MultiVectorOut.h"
00098 
00099 using namespace std;
00100 using namespace Intrepid;
00101 
00102 // Functions to evaluate exact solution and derivatives
00103 double evalu(double & x, double & y, double & z);
00104 int evalGradu(double & x, double & y, double & z, double & gradu1, double & gradu2, double & gradu3);
00105 double evalDivGradu(double & x, double & y, double & z);
00106 
00107 int main(int argc, char *argv[]) {
00108 
00109   //Check number of arguments
00110    if (argc < 4) {
00111       std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n";
00112       std::cout <<"Usage:\n\n";
00113       std::cout <<"  ./Intrepid_example_Drivers_Example_05.exe deg NX NY verbose\n\n";
00114       std::cout <<" where \n";
00115       std::cout <<"   int deg             - polynomial degree to be used (assumed > 1) \n";
00116       std::cout <<"   int NX              - num intervals in x direction (assumed box domain, 0,1) \n";
00117       std::cout <<"   int NY              - num intervals in y direction (assumed box domain, 0,1) \n";
00118       std::cout <<"   verbose (optional)  - any character, indicates verbose output \n\n";
00119       exit(1);
00120    }
00121   
00122   // This little trick lets us print to std::cout only if
00123   // a (dummy) command-line argument is provided.
00124   int iprint     = argc - 1;
00125   Teuchos::RCP<std::ostream> outStream;
00126   Teuchos::oblackholestream bhs; // outputs nothing
00127   if (iprint > 2)
00128     outStream = Teuchos::rcp(&std::cout, false);
00129   else
00130     outStream = Teuchos::rcp(&bhs, false);
00131   
00132   // Save the format state of the original std::cout.
00133   Teuchos::oblackholestream oldFormatState;
00134   oldFormatState.copyfmt(std::cout);
00135   
00136   *outStream \
00137     << "===============================================================================\n" \
00138     << "|                                                                             |\n" \
00139     << "|  Example: Generate Stiffness Matrix and Right Hand Side Vector for          |\n" \
00140     << "|                   Poisson Equation on Quadrilateral Mesh                    |\n" \
00141     << "|                                                                             |\n" \
00142     << "|  Questions? Contact  Pavel Bochev  (pbboche@sandia.gov),                    |\n" \
00143     << "|                      Denis Ridzal  (dridzal@sandia.gov),                    |\n" \
00144     << "|                      Kara Peterson (kjpeter@sandia.gov).                    |\n" \
00145     << "|                                                                             |\n" \
00146     << "|  Intrepid's website: http://trilinos.sandia.gov/packages/intrepid           |\n" \
00147     << "|  Trilinos website:   http://trilinos.sandia.gov                             |\n" \
00148     << "|                                                                             |\n" \
00149     << "===============================================================================\n";
00150 
00151   
00152   // ************************************ GET INPUTS **************************************
00153   
00154   int deg          = atoi(argv[1]);  // polynomial degree to use
00155   int NX            = atoi(argv[2]);  // num intervals in x direction (assumed box domain, 0,1)
00156   int NY            = atoi(argv[3]);  // num intervals in y direction (assumed box domain, 0,1)
00157   
00158 
00159   // *********************************** CELL TOPOLOGY **********************************
00160   
00161   // Get cell topology for base hexahedron
00162   typedef shards::CellTopology    CellTopology;
00163   CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
00164   
00165   // Get dimensions 
00166   int numNodesPerElem = quad_4.getNodeCount();
00167   int spaceDim = quad_4.getDimension();
00168   
00169   // *********************************** GENERATE MESH ************************************
00170   
00171   *outStream << "Generating mesh ... \n\n";
00172   
00173   *outStream << "   NX" << "   NY\n";
00174   *outStream << std::setw(5) << NX <<
00175     std::setw(5) << NY << "\n\n";
00176   
00177   // Print mesh information
00178   int numElems = NX*NY;
00179   int numNodes = (NX+1)*(NY+1);
00180   *outStream << " Number of Elements: " << numElems << " \n";
00181   *outStream << "    Number of Nodes: " << numNodes << " \n\n";
00182   
00183   // Square
00184   double leftX = 0.0, rightX = 1.0;
00185   double leftY = 0.0, rightY = 1.0;
00186 
00187   // Mesh spacing
00188   double hx = (rightX-leftX)/((double)NX);
00189   double hy = (rightY-leftY)/((double)NY);
00190 
00191   // Get nodal coordinates
00192   FieldContainer<double> nodeCoord(numNodes, spaceDim);
00193   FieldContainer<int> nodeOnBoundary(numNodes);
00194   int inode = 0;
00195   for (int j=0; j<NY+1; j++) {
00196     for (int i=0; i<NX+1; i++) {
00197       nodeCoord(inode,0) = leftX + (double)i*hx;
00198       nodeCoord(inode,1) = leftY + (double)j*hy;
00199       if (j==0 || i==0 || j==NY || i==NX){
00200         nodeOnBoundary(inode)=1;
00201       }
00202       else {
00203         nodeOnBoundary(inode)=0;
00204       }
00205       inode++;
00206     }
00207   }
00208 #define DUMP_DATA
00209 #ifdef DUMP_DATA
00210   // Print nodal coords
00211   ofstream fcoordout("coords.dat");
00212   for (int i=0; i<numNodes; i++) {
00213     fcoordout << nodeCoord(i,0) <<" ";
00214     fcoordout << nodeCoord(i,1) <<"\n";
00215   }
00216   fcoordout.close();
00217 #endif
00218   
00219   
00220   // Element to Node map
00221   // We'll keep it around, but this is only the DOFMap if you are in the lowest order case.
00222   FieldContainer<int> elemToNode(numElems, numNodesPerElem);
00223   int ielem = 0;
00224   for (int j=0; j<NY; j++) {
00225     for (int i=0; i<NX; i++) {
00226       elemToNode(ielem,0) = (NX + 1)*j + i;
00227       elemToNode(ielem,1) = (NX + 1)*j + i + 1;
00228       elemToNode(ielem,2) = (NX + 1)*(j + 1) + i + 1;
00229       elemToNode(ielem,3) = (NX + 1)*(j + 1) + i;
00230       ielem++;
00231     }
00232   }
00233 #ifdef DUMP_DATA
00234   // Output connectivity
00235   ofstream fe2nout("elem2node.dat");
00236   for (int j=0; j<NY; j++) {
00237     for (int i=0; i<NX; i++) {
00238       int ielem = i + j * NX;
00239       for (int m=0; m<numNodesPerElem; m++){
00240         fe2nout << elemToNode(ielem,m) <<"  ";
00241       }
00242       fe2nout <<"\n";
00243     }
00244   }
00245   fe2nout.close();
00246 #endif
00247   
00248 
00249   // ************************************ CUBATURE ************************************** 
00250   *outStream << "Getting cubature ... \n\n";
00251   
00252   // Get numerical integration points and weights
00253   DefaultCubatureFactory<double>  cubFactory;                                   
00254   int cubDegree = 2*deg;
00255   Teuchos::RCP<Cubature<double> > quadCub = cubFactory.create(quad_4, cubDegree); 
00256   
00257   int cubDim       = quadCub->getDimension();
00258   int numCubPoints = quadCub->getNumPoints();
00259   
00260   FieldContainer<double> cubPoints(numCubPoints, cubDim);
00261   FieldContainer<double> cubWeights(numCubPoints);
00262   
00263   quadCub->getCubature(cubPoints, cubWeights);
00264   
00265 
00266   // ************************************** BASIS ***************************************
00267   
00268   *outStream << "Getting basis ... \n\n";
00269   
00270   // Define basis 
00271   Basis_HGRAD_QUAD_Cn_FEM<double, FieldContainer<double> > quadHGradBasis(deg,POINTTYPE_SPECTRAL);
00272   int numFieldsG = quadHGradBasis.getCardinality();
00273   FieldContainer<double> quadGVals(numFieldsG, numCubPoints); 
00274   FieldContainer<double> quadGrads(numFieldsG, numCubPoints, spaceDim); 
00275   
00276   // Evaluate basis values and gradients at cubature points
00277   quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
00278   quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
00279 
00280   // create the local-global mapping for higher order elements
00281   FieldContainer<int> ltgMapping(numElems,numFieldsG);
00282   const int numDOF = (NX*deg+1)*(NY*deg+1);
00283   ielem=0;
00284   for (int j=0;j<NY;j++) {
00285     for (int i=0;i<NX;i++) {
00286       const int start = deg * j * ( NX * deg + 1 ) + i * deg;
00287       // loop over local dof on this cell
00288       int local_dof_cur=0;
00289       for (int vertical=0;vertical<=deg;vertical++) {
00290         for (int horizontal=0;horizontal<=deg;horizontal++) {
00291           ltgMapping(ielem,local_dof_cur) = start + vertical*(NX*deg+1)+horizontal;
00292           local_dof_cur++;
00293         }
00294       }
00295       ielem++;
00296     }
00297   }
00298 #ifdef DUMP_DATA
00299   // Output ltg mapping
00300   ofstream ltgout("ltg.dat");
00301   for (int j=0; j<NY; j++) {
00302     for (int i=0; i<NX; i++) {
00303       int ielem = i + j * NX;
00304       for (int m=0; m<numFieldsG; m++){
00305         ltgout << ltgMapping(ielem,m) <<"  ";
00306       }
00307       ltgout <<"\n";
00308     }
00309   }
00310   ltgout.close();
00311 #endif
00312   
00313   // ******** CREATE A SINGLE STIFFNESS MATRIX, WHICH IS REPLICATED ON ALL ELEMENTS *********
00314   *outStream << "Building stiffness matrix and right hand side ... \n\n";
00315 
00316   // Settings and data structures for mass and stiffness matrices
00317   typedef CellTools<double>  CellTools;
00318   typedef FunctionSpaceTools fst;
00319   int numCells = 1; 
00320 
00321   // Container for nodes
00322   FieldContainer<double> refQuadNodes(numCells, numNodesPerElem, spaceDim);
00323   // Containers for Jacobian
00324   FieldContainer<double> refQuadJacobian(numCells, numCubPoints, spaceDim, spaceDim);
00325   FieldContainer<double> refQuadJacobInv(numCells, numCubPoints, spaceDim, spaceDim);
00326   FieldContainer<double> refQuadJacobDet(numCells, numCubPoints);
00327   // Containers for element HGRAD stiffness matrix
00328   FieldContainer<double> localStiffMatrix(numCells, numFieldsG, numFieldsG);
00329   FieldContainer<double> weightedMeasure(numCells, numCubPoints);
00330   FieldContainer<double> quadGradsTransformed(numCells, numFieldsG, numCubPoints, spaceDim);
00331   FieldContainer<double> quadGradsTransformedWeighted(numCells, numFieldsG, numCubPoints, spaceDim);
00332   // Containers for right hand side vectors
00333   FieldContainer<double> rhsData(numCells, numCubPoints);
00334   FieldContainer<double> localRHS(numCells, numFieldsG);
00335   FieldContainer<double> quadGValsTransformed(numCells, numFieldsG, numCubPoints);
00336   FieldContainer<double> quadGValsTransformedWeighted(numCells, numFieldsG, numCubPoints);
00337   // Container for cubature points in physical space
00338   FieldContainer<double> physCubPoints(numCells, numCubPoints, cubDim);
00339   
00340   // Global arrays in Epetra format 
00341   Epetra_SerialComm Comm;
00342   Epetra_Map globalMapG(numDOF, 0, Comm);
00343   Epetra_Time instantiateTimer(Comm);
00344   Epetra_FECrsMatrix StiffMatrix(Copy, globalMapG, 4*numFieldsG);
00345   const double instantiateTime = instantiateTimer.ElapsedTime();
00346   std::cout << "Time to instantiate sparse matrix " << instantiateTime << "\n";
00347   Epetra_FEVector u(globalMapG);
00348   Epetra_FEVector Ku(globalMapG);
00349 
00350   u.Random();
00351     
00352   // ************************** Compute element HGrad stiffness matrices *******************************  
00353   refQuadNodes(0,0,0) = 0.0;
00354   refQuadNodes(0,0,1) = 0.0;
00355   refQuadNodes(0,1,0) = hx;
00356   refQuadNodes(0,1,1) = 0.0;
00357   refQuadNodes(0,2,0) = hx;
00358   refQuadNodes(0,2,1) = hy;
00359   refQuadNodes(0,3,0) = 0.0;
00360   refQuadNodes(0,3,1) = hy;
00361 
00362   // Compute cell Jacobians, their inverses and their determinants
00363   CellTools::setJacobian(refQuadJacobian, cubPoints, refQuadNodes, quad_4);
00364   CellTools::setJacobianInv(refQuadJacobInv, refQuadJacobian );
00365   CellTools::setJacobianDet(refQuadJacobDet, refQuadJacobian );
00366   
00367   // transform from [-1,1]^2 to [0,hx]x[0,hy]
00368   fst::HGRADtransformGRAD<double>(quadGradsTransformed, refQuadJacobInv, quadGrads);
00369       
00370   // compute weighted measure
00371   fst::computeCellMeasure<double>(weightedMeasure, refQuadJacobDet, cubWeights);
00372 
00373   // multiply values with weighted measure
00374   fst::multiplyMeasure<double>(quadGradsTransformedWeighted,
00375                                weightedMeasure, quadGradsTransformed);
00376 
00377   // integrate to compute element stiffness matrix
00378   fst::integrate<double>(localStiffMatrix,
00379                          quadGradsTransformed, quadGradsTransformedWeighted, COMP_BLAS);
00380 
00381 
00382   Epetra_Time assemblyTimer(Comm);
00383 
00384   // *** Element loop ***
00385    for (int k=0; k<numElems; k++) 
00386      {
00387        // assemble into global matrix
00388        StiffMatrix.InsertGlobalValues(numFieldsG,&ltgMapping(k,0),numFieldsG,&ltgMapping(k,0),&localStiffMatrix(0,0,0));
00389 
00390      }
00391 
00392 
00393   // Assemble global matrices
00394    StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
00395 
00396    double assembleTime = assemblyTimer.ElapsedTime();
00397    std::cout << "Time to insert reference element matrix into global matrix: " << assembleTime << std::endl;
00398    std::cout << "There are " << StiffMatrix.NumGlobalNonzeros() << " nonzeros in the matrix.\n";
00399    std::cout << "There are " << numDOF << " global degrees of freedom.\n";
00400  
00401    Epetra_Time multTimer(Comm);
00402    StiffMatrix.Apply(u,Ku);
00403    double multTime = multTimer.ElapsedTime();
00404    std::cout << "Time to apply: " << multTime << std::endl;
00405 
00406 //    // Adjust stiffness matrix and rhs based on boundary conditions
00407 //    for (int row = 0; row<numNodes; row++){
00408 //        if (nodeOnBoundary(row)) {
00409 //           int rowindex = row;
00410 //           for (int col=0; col<numNodes; col++){
00411 //               double val = 0.0;
00412 //               int colindex = col;
00413 //               StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &colindex, &val);
00414 //           }
00415 //           double val = 1.0;
00416 //           StiffMatrix.ReplaceGlobalValues(1, &rowindex, 1, &rowindex, &val);
00417 //           val = 0.0;
00418 //           rhs.ReplaceGlobalValues(1, &rowindex, &val);
00419 //        }
00420 //     }
00421 
00422 #ifdef DUMP_DATA
00423    // Dump matrices to disk
00424 //    EpetraExt::RowMatrixToMatlabFile("stiff_matrix.dat",StiffMatrix);
00425 //    EpetraExt::MultiVectorToMatrixMarketFile("rhs_vector.dat",rhs,0,0,false);
00426 #endif
00427 
00428    std::cout << "End Result: TEST PASSED\n";   
00429 
00430    // reset format state of std::cout
00431    std::cout.copyfmt(oldFormatState);
00432    
00433    return 0;
00434 }
00435 
00436 
00437 // Calculates value of exact solution u
00438  double evalu(double & x, double & y, double & z)
00439  {
00440  /*
00441    // function1
00442     double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z);
00443  */
00444 
00445    // function2
00446    double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
00447 
00448    return exactu;
00449  }
00450 
00451 // Calculates gradient of exact solution u
00452  int evalGradu(double & x, double & y, double & z, double & gradu1, double & gradu2, double & gradu3)
00453  {
00454  /*
00455    // function 1
00456        gradu1 = M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z);
00457        gradu2 = M_PI*sin(M_PI*x)*cos(M_PI*y)*sin(M_PI*z);
00458        gradu3 = M_PI*sin(M_PI*x)*sin(M_PI*y)*cos(M_PI*z);
00459  */
00460 
00461    // function2
00462        gradu1 = (M_PI*cos(M_PI*x)+sin(M_PI*x))
00463                   *sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
00464        gradu2 = (M_PI*cos(M_PI*y)+sin(M_PI*y))
00465                   *sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z);
00466        gradu3 = (M_PI*cos(M_PI*z)+sin(M_PI*z))
00467                   *sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z);
00468   
00469    return 0;
00470  }
00471 
00472 // Calculates Laplacian of exact solution u
00473  double evalDivGradu(double & x, double & y, double & z)
00474  {
00475  /*
00476    // function 1
00477     double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z);
00478  */
00479 
00480    // function 2
00481    double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
00482                     + 2.0*M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
00483                     + 2.0*M_PI*cos(M_PI*y)*sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z)
00484                     + 2.0*M_PI*cos(M_PI*z)*sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z)
00485                     + 3.0*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
00486    
00487    return divGradu;
00488  }
00489