94 #include "Epetra_Time.h"
95 #include "Epetra_Map.h"
96 #include "Epetra_FEVector.h"
97 #include "Epetra_FECrsMatrix.h"
98 #include "Epetra_SerialComm.h"
101 #include "Teuchos_oblackholestream.hpp"
102 #include "Teuchos_RCP.hpp"
107 #include "Shards_CellTopology.hpp"
110 #include "EpetraExt_MultiVectorOut.h"
113 using namespace Intrepid;
115 int main(
int argc,
char *argv[]) {
119 std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
120 std::cout <<
"Usage:\n\n";
121 std::cout <<
" ./Intrepid_example_Drivers_Example_10.exe deg NX NY NZ verbose\n\n";
122 std::cout <<
" where \n";
123 std::cout <<
" int deg - polynomial degree to be used (assumed >= 1) \n";
124 std::cout <<
" int NX - num intervals in x direction (assumed box domain, 0,1) \n";
125 std::cout <<
" int NY - num intervals in y direction (assumed box domain, 0,1) \n";
126 std::cout <<
" int NZ - num intervals in y direction (assumed box domain, 0,1) \n";
127 std::cout <<
" verbose (optional) - any character, indicates verbose output \n\n";
133 int iprint = argc - 1;
134 Teuchos::RCP<std::ostream> outStream;
135 Teuchos::oblackholestream bhs;
137 outStream = Teuchos::rcp(&std::cout,
false);
139 outStream = Teuchos::rcp(&bhs,
false);
142 Teuchos::oblackholestream oldFormatState;
143 oldFormatState.copyfmt(std::cout);
146 <<
"===============================================================================\n" \
148 <<
"| Example: Build Stiffness Matrix for |\n" \
149 <<
"| Poisson Equation on Hexahedral Mesh |\n" \
151 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
152 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
153 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
155 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
156 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
158 <<
"===============================================================================\n";
163 int deg = atoi(argv[1]);
164 int NX = atoi(argv[2]);
165 int NY = atoi(argv[3]);
166 int NZ = atoi(argv[4]);
172 typedef shards::CellTopology CellTopology;
173 CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
176 int numNodesPerElem = hex_8.getNodeCount();
177 int spaceDim = hex_8.getDimension();
181 *outStream <<
"Generating mesh ... \n\n";
183 *outStream <<
" NX" <<
" NY" <<
" NZ\n";
184 *outStream << std::setw(5) << NX <<
185 std::setw(5) << NY << std::setw(5) << NZ <<
"\n\n";
188 int numElems = NX*NY*NZ;
189 int numNodes = (NX+1)*(NY+1)*(NZ+1);
190 *outStream <<
" Number of Elements: " << numElems <<
" \n";
191 *outStream <<
" Number of Nodes: " << numNodes <<
" \n\n";
194 double leftX = 0.0, rightX = 1.0;
195 double leftY = 0.0, rightY = 1.0;
196 double leftZ = 0.0, rightZ = 1.0;
199 double hx = (rightX-leftX)/((
double)NX);
200 double hy = (rightY-leftY)/((
double)NY);
201 double hz = (rightZ-leftZ)/((
double)NZ);
207 for (
int k=0; k<NZ+1; k++)
209 for (
int j=0; j<NY+1; j++)
211 for (
int i=0; i<NX+1; i++)
213 nodeCoord(inode,0) = leftX + (double)i*hx;
214 nodeCoord(inode,1) = leftY + (double)j*hy;
215 nodeCoord(inode,2) = leftZ + (double)k*hz;
216 if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
218 nodeOnBoundary(inode)=1;
222 nodeOnBoundary(inode)=0;
231 ofstream fcoordout(
"coords.dat");
232 for (
int i=0; i<numNodes; i++) {
233 fcoordout << nodeCoord(i,0) <<
" ";
234 fcoordout << nodeCoord(i,1) <<
" ";
235 fcoordout << nodeCoord(i,2) <<
"\n";
245 for (
int k=0; k<NZ; k++)
247 for (
int j=0; j<NY; j++)
249 for (
int i=0; i<NX; i++)
251 elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
252 elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
253 elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
254 elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
255 elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
256 elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
257 elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
258 elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
265 ofstream fe2nout(
"elem2node.dat");
266 for (
int k=0;k<NZ;k++)
268 for (
int j=0; j<NY; j++)
270 for (
int i=0; i<NX; i++)
272 int ielem = i + j * NX + k * NY * NY;
273 for (
int m=0; m<numNodesPerElem; m++)
275 fe2nout << elemToNode(ielem,m) <<
" ";
285 *outStream <<
"Getting cubature ... \n\n";
289 int cubDegree = 2*deg;
290 Teuchos::RCP<Cubature<double> > quadCub = cubFactory.
create(hex_8, cubDegree);
292 int cubDim = quadCub->getDimension();
293 int numCubPoints = quadCub->getNumPoints();
298 quadCub->getCubature(cubPoints, cubWeights);
303 *outStream <<
"Getting basis ... \n\n";
307 int numFieldsG = quadHGradBasis.getCardinality();
312 quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
313 quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
317 const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
319 for (
int k=0;k<NZ;k++)
321 for (
int j=0;j<NY;j++)
323 for (
int i=0;i<NX;i++)
325 const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
328 for (
int kloc=0;kloc<=deg;kloc++)
330 for (
int jloc=0;jloc<=deg;jloc++)
332 for (
int iloc=0;iloc<=deg;iloc++)
334 ltgMapping(ielem,local_dof_cur) = start
335 + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
336 + jloc * ( NX * deg + 1 )
349 ofstream ltgout(
"ltg.dat");
350 for (
int k=0;k<NZ;k++)
352 for (
int j=0; j<NY; j++)
354 for (
int i=0; i<NX; i++)
356 int ielem = i + j * NX + k * NX * NY;
357 for (
int m=0; m<numFieldsG; m++)
359 ltgout << ltgMapping(ielem,m) <<
" ";
369 Epetra_SerialComm Comm;
370 Epetra_Map globalMapG(numDOF, 0, Comm);
371 Epetra_FEVector u(globalMapG); u.Random();
372 Epetra_FEVector Ku(globalMapG);
375 Epetra_Time graphTimer(Comm);
376 Epetra_CrsGraph grph( Copy , globalMapG , 4 * numFieldsG );
377 for (
int k=0;k<numElems;k++)
379 for (
int i=0;i<numFieldsG;i++)
381 grph.InsertGlobalIndices(ltgMapping(k,i),numFieldsG,<gMapping(k,0));
385 const double graphTime = graphTimer.ElapsedTime();
389 Epetra_Time instantiateTimer(Comm);
390 Epetra_FECrsMatrix StiffMatrix(Copy,grph);
391 const double instantiateTime = instantiateTimer.ElapsedTime();
395 *outStream <<
"Building local stiffness matrices...\n\n";
414 Epetra_Time localConstructTimer( Comm );
415 refCellNodes(0,0,0) = 0.0; refCellNodes(0,0,1) = 0.0; refCellNodes(0,0,2) = 0.0;
416 refCellNodes(0,1,0) = hx; refCellNodes(0,1,1) = 0.0; refCellNodes(0,1,2) = 0.0;
417 refCellNodes(0,2,0) = hx; refCellNodes(0,2,1) = hy; refCellNodes(0,2,2) = 0.0;
418 refCellNodes(0,3,0) = 0.0; refCellNodes(0,3,1) = hy; refCellNodes(0,3,2) = 0.0;
419 refCellNodes(0,4,0) = 0.0; refCellNodes(0,4,1) = 0.0; refCellNodes(0,4,2) = hz;
420 refCellNodes(0,5,0) = hx; refCellNodes(0,5,1) = 0.0; refCellNodes(0,5,2) = hz;
421 refCellNodes(0,6,0) = hx; refCellNodes(0,6,1) = hy; refCellNodes(0,6,2) = hz;
422 refCellNodes(0,7,0) = 0.0; refCellNodes(0,7,1) = hy; refCellNodes(0,7,2) = hz;
426 CellTools::setJacobian(cellJacobian,cubPoints,refCellNodes,hex_8);
427 CellTools::setJacobianInv(cellJacobInv, cellJacobian );
428 CellTools::setJacobianDet(cellJacobDet, cellJacobian );
431 fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
434 fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
437 fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
438 weightedMeasure, transformedBasisGradients);
441 fst::integrate<double>(localStiffMatrix,
442 transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
444 const double localConstructTime = localConstructTimer.ElapsedTime();
447 Epetra_Time insertionTimer(Comm);
450 for (
int k=0; k<numElems; k++)
453 StiffMatrix.InsertGlobalValues(numFieldsG,<gMapping(k,0),numFieldsG,<gMapping(k,0),&localStiffMatrix(0,0,0));
456 StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
457 const double insertionTime = insertionTimer.ElapsedTime( );
459 *outStream <<
"Time to construct matrix graph: " << graphTime <<
"\n";
460 *outStream <<
"Time to instantiate global stiffness matrix: " << instantiateTime <<
"\n";
461 *outStream <<
"Time to build local matrices (including Jacobian computation): "<< localConstructTime <<
"\n";
462 *outStream <<
"Time to assemble global matrix from local matrices: " << insertionTime <<
"\n";
463 *outStream <<
"Total construction time: " << graphTime + instantiateTime + localConstructTime + insertionTime <<
"\n";
465 Epetra_Time applyTimer(Comm);
466 StiffMatrix.Apply(u,Ku);
467 const double multTime = applyTimer.ElapsedTime();
468 *outStream <<
"Time to multiply onto a vector: " << multTime <<
"\n";
470 *outStream <<
"End Result: TEST PASSED\n";
473 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.