93 #include "Epetra_Time.h" 
   94 #include "Epetra_Map.h" 
   95 #include "Epetra_FECrsMatrix.h" 
   96 #include "Epetra_FEVector.h" 
   97 #include "Epetra_SerialComm.h" 
  100 #include "Teuchos_oblackholestream.hpp" 
  101 #include "Teuchos_RCP.hpp" 
  102 #include "Teuchos_BLAS.hpp" 
  105 #include "Shards_CellTopology.hpp" 
  108 #include "EpetraExt_RowMatrixOut.h" 
  109 #include "EpetraExt_MultiVectorOut.h" 
  112 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_09.exe deg NX NY 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 <<
"   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: Apply Stiffness Matrix for                                        |\n" \
 
  147     << 
"|                   Poisson Equation on Quadrilateral 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]);  
 
  169   typedef shards::CellTopology    CellTopology;
 
  170   CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
 
  173   int numNodesPerElem = quad_4.getNodeCount();
 
  174   int spaceDim = quad_4.getDimension();
 
  178   *outStream << 
"Generating mesh ... \n\n";
 
  180   *outStream << 
"   NX" << 
"   NY\n";
 
  181   *outStream << std::setw(5) << NX <<
 
  182     std::setw(5) << NY << 
"\n\n";
 
  185   int numElems = NX*NY;
 
  186   int numNodes = (NX+1)*(NY+1);
 
  187   *outStream << 
" Number of Elements: " << numElems << 
" \n";
 
  188   *outStream << 
"    Number of Nodes: " << numNodes << 
" \n\n";
 
  191   double leftX = 0.0, rightX = 1.0;
 
  192   double leftY = 0.0, rightY = 1.0;
 
  195   double hx = (rightX-leftX)/((
double)NX);
 
  196   double hy = (rightY-leftY)/((
double)NY);
 
  202   for (
int j=0; j<NY+1; j++) {
 
  203     for (
int i=0; i<NX+1; i++) {
 
  204       nodeCoord(inode,0) = leftX + (double)i*hx;
 
  205       nodeCoord(inode,1) = leftY + (double)j*hy;
 
  206       if (j==0 || i==0 || j==NY || i==NX){
 
  207         nodeOnBoundary(inode)=1;
 
  210         nodeOnBoundary(inode)=0;
 
  218   ofstream fcoordout(
"coords.dat");
 
  219   for (
int i=0; i<numNodes; i++) {
 
  220     fcoordout << nodeCoord(i,0) <<
" ";
 
  221     fcoordout << nodeCoord(i,1) <<
"\n";
 
  231   for (
int j=0; j<NY; j++) {
 
  232     for (
int i=0; i<NX; i++) {
 
  233       elemToNode(ielem,0) = (NX + 1)*j + i;
 
  234       elemToNode(ielem,1) = (NX + 1)*j + i + 1;
 
  235       elemToNode(ielem,2) = (NX + 1)*(j + 1) + i + 1;
 
  236       elemToNode(ielem,3) = (NX + 1)*(j + 1) + i;
 
  242   ofstream fe2nout(
"elem2node.dat");
 
  243   for (
int j=0; j<NY; j++) {
 
  244     for (
int i=0; i<NX; i++) {
 
  245       int ielem = i + j * NX;
 
  246       for (
int m=0; m<numNodesPerElem; m++){
 
  247         fe2nout << elemToNode(ielem,m) <<
"  ";
 
  257   *outStream << 
"Getting cubature ... \n\n";
 
  266   const int numCubPoints = glcub->getNumPoints();
 
  271   glcub->getCubature(cubPoints1D,cubWeights1D);
 
  275   *outStream << 
"Getting basis ... \n\n";
 
  279   int numLineFieldsG = lineHGradBasis.getCardinality();
 
  283   lineHGradBasis.getValues(lineGrads, cubPoints1D, OPERATOR_GRAD);
 
  287   const int numDOF = (NX*deg+1)*(NY*deg+1);
 
  289   for (
int j=0;j<NY;j++) {
 
  290     for (
int i=0;i<NX;i++) {
 
  291       const int start = deg * j * ( NX * deg + 1 ) + i * deg;
 
  294       for (
int vertical=0;vertical<=deg;vertical++) {
 
  295         for (
int horizontal=0;horizontal<=deg;horizontal++) {
 
  296           ltgMapping(ielem,local_dof_cur) = start + vertical*(NX*deg+1)+horizontal;
 
  305   ofstream ltgout(
"ltg.dat");
 
  306   for (
int j=0; j<NY; j++) {
 
  307     for (
int i=0; i<NX; i++) {
 
  308       int ielem = i + j * NX;
 
  309       for (
int m=0; m<numLineFieldsG; m++){
 
  310         ltgout << ltgMapping(ielem,m) <<
"  ";
 
  320   Epetra_SerialComm Comm;
 
  321   Epetra_Map globalMapG(numDOF, 0, Comm);
 
  323   Epetra_FEVector u(globalMapG);
 
  324   Epetra_FEVector Ku(globalMapG);
 
  354   double *uVals = u[0];
 
  355   double *KuVals = Ku[0];
 
  356   Epetra_Time scatterTime(Comm);
 
  357   *outStream << 
"Scattering\n";
 
  359   for (
int k=0; k<numElems; k++) 
 
  361       for (
int i=0;i<numLineFieldsG*numLineFieldsG;i++) 
 
  363           uScattered(k,i) = uVals[ltgMapping(k,i)];
 
  366   const double scatTime = scatterTime.ElapsedTime();
 
  367   *outStream << 
"Scattered in time " << scatTime << 
"\n";
 
  369   Epetra_Time applyTimer(Comm);
 
  371   uScattered.resize(numElems,numLineFieldsG,numLineFieldsG);
 
  373   for (
int k=0;k<numElems;k++)
 
  377       for (
int j=0;j<numLineFieldsG;j++)
 
  379           for (
int i=0;i<numLineFieldsG;i++)
 
  382               for (
int q=0;q<numLineFieldsG;q++)
 
  384                   Du(j,i) += uScattered(k,j,i) * lineGrads(i,q,0);
 
  390       for (
int i=0;i<numLineFieldsG*numLineFieldsG;i++)
 
  392           KuScattered(k,i) = 0.0;
 
  397       for (
int j=0;j<numLineFieldsG;j++)
 
  399           for (
int i=0;i<numLineFieldsG;i++)
 
  402               for (
int q1=0;q1<numCubPoints;q1++)
 
  404                   KuScattered(k,cur) += cubWeights1D(j) * cubWeights1D(q1) * Du(j,q1) * lineGrads(i,q1,0);
 
  411       for (
int j=0;j<numLineFieldsG;j++)
 
  413           for (
int i=0;i<numLineFieldsG;i++)
 
  416               for (
int q=0;q<numLineFieldsG;q++)
 
  418                   Du(j,i) += uScattered(k,j,i) * lineGrads(j,q,0);
 
  426       for (
int j=0;j<numLineFieldsG;j++)
 
  428           for (
int i=0;i<numLineFieldsG;i++)
 
  430               for (
int q2=0;q2<numCubPoints;q2++)
 
  432                   KuScattered(k,cur) += cubWeights1D(q2) * Du(q2,i) * lineGrads(j,q2,0) * cubWeights1D(i);
 
  438   uScattered.resize(numElems,numLineFieldsG*numLineFieldsG);
 
  440   const double applyTime = applyTimer.ElapsedTime();
 
  442   *outStream << 
"Local application: " << applyTime << 
"\n";
 
  445   Epetra_Time gatherTimer(Comm);
 
  447   for (
int k=0;k<numElems;k++)
 
  449       for (
int i=0;i<numLineFieldsG*numLineFieldsG;i++)
 
  451           KuVals[ltgMapping(k,i)] += KuScattered(k,i);
 
  455   const double gatherTime = gatherTimer.ElapsedTime();
 
  456   *outStream << 
"Gathered in " << gatherTime << 
"\n";
 
  467   std::cout << 
"End Result: TEST PASSED\n";
 
  470    std::cout.copyfmt(oldFormatState);
 
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
Header file for utility class to provide multidimensional containers. 
Header file for the abstract base class Intrepid::DefaultCubatureFactory. 
Utilizes cubature (integration) rules contained in the library Polylib (Spencer Sherwin, Aeronautics, Imperial College London) within Intrepid. 
Header file for the Intrepid::HGRAD_QUAD_Cn_FEM class. 
A factory class that generates specific instances of cubatures. 
Header file for the Intrepid::CubaturePolylib class.