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;
115 double evalu(
double & x,
double & y,
double & z);
116 int evalGradu(
double & x,
double & y,
double & z,
double & gradu1,
double & gradu2,
double & gradu3);
117 double evalDivGradu(
double & x,
double & y,
double & z);
119 int main(
int argc,
char *argv[]) {
123 std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
124 std::cout <<
"Usage:\n\n";
125 std::cout <<
" ./Intrepid_example_Drivers_Example_05.exe deg NX NY verbose\n\n";
126 std::cout <<
" where \n";
127 std::cout <<
" int deg - polynomial degree to be used (assumed > 1) \n";
128 std::cout <<
" int NX - num intervals in x direction (assumed box domain, 0,1) \n";
129 std::cout <<
" int NY - num intervals in y direction (assumed box domain, 0,1) \n";
130 std::cout <<
" verbose (optional) - any character, indicates verbose output \n\n";
136 int iprint = argc - 1;
137 Teuchos::RCP<std::ostream> outStream;
138 Teuchos::oblackholestream bhs;
140 outStream = Teuchos::rcp(&std::cout,
false);
142 outStream = Teuchos::rcp(&bhs,
false);
145 Teuchos::oblackholestream oldFormatState;
146 oldFormatState.copyfmt(std::cout);
149 <<
"===============================================================================\n" \
151 <<
"| Example: Generate Stiffness Matrix and Right Hand Side Vector for |\n" \
152 <<
"| Poisson Equation on Quadrilateral Mesh |\n" \
154 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
155 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
156 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
158 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
159 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
161 <<
"===============================================================================\n";
166 int deg = atoi(argv[1]);
167 int NX = atoi(argv[2]);
168 int NY = atoi(argv[3]);
174 typedef shards::CellTopology CellTopology;
175 CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
178 int numNodesPerElem = quad_4.getNodeCount();
179 int spaceDim = quad_4.getDimension();
183 *outStream <<
"Generating mesh ... \n\n";
185 *outStream <<
" NX" <<
" NY\n";
186 *outStream << std::setw(5) << NX <<
187 std::setw(5) << NY <<
"\n\n";
190 int numElems = NX*NY;
191 int numNodes = (NX+1)*(NY+1);
192 *outStream <<
" Number of Elements: " << numElems <<
" \n";
193 *outStream <<
" Number of Nodes: " << numNodes <<
" \n\n";
196 double leftX = 0.0, rightX = 1.0;
197 double leftY = 0.0, rightY = 1.0;
200 double hx = (rightX-leftX)/((
double)NX);
201 double hy = (rightY-leftY)/((
double)NY);
207 for (
int j=0; j<NY+1; j++) {
208 for (
int i=0; i<NX+1; i++) {
209 nodeCoord(inode,0) = leftX + (double)i*hx;
210 nodeCoord(inode,1) = leftY + (double)j*hy;
211 if (j==0 || i==0 || j==NY || i==NX){
212 nodeOnBoundary(inode)=1;
215 nodeOnBoundary(inode)=0;
223 ofstream fcoordout(
"coords.dat");
224 for (
int i=0; i<numNodes; i++) {
225 fcoordout << nodeCoord(i,0) <<
" ";
226 fcoordout << nodeCoord(i,1) <<
"\n";
236 for (
int j=0; j<NY; j++) {
237 for (
int i=0; i<NX; i++) {
238 elemToNode(ielem,0) = (NX + 1)*j + i;
239 elemToNode(ielem,1) = (NX + 1)*j + i + 1;
240 elemToNode(ielem,2) = (NX + 1)*(j + 1) + i + 1;
241 elemToNode(ielem,3) = (NX + 1)*(j + 1) + i;
247 ofstream fe2nout(
"elem2node.dat");
248 for (
int j=0; j<NY; j++) {
249 for (
int i=0; i<NX; i++) {
250 int ielem = i + j * NX;
251 for (
int m=0; m<numNodesPerElem; m++){
252 fe2nout << elemToNode(ielem,m) <<
" ";
262 *outStream <<
"Getting cubature ... \n\n";
266 int cubDegree = 2*deg;
267 Teuchos::RCP<Cubature<double> > quadCub = cubFactory.
create(quad_4, cubDegree);
269 int cubDim = quadCub->getDimension();
270 int numCubPoints = quadCub->getNumPoints();
275 quadCub->getCubature(cubPoints, cubWeights);
280 *outStream <<
"Getting basis ... \n\n";
284 int numFieldsG = quadHGradBasis.getCardinality();
289 quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
290 quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
294 const int numDOF = (NX*deg+1)*(NY*deg+1);
296 for (
int j=0;j<NY;j++) {
297 for (
int i=0;i<NX;i++) {
298 const int start = deg * j * ( NX * deg + 1 ) + i * deg;
301 for (
int vertical=0;vertical<=deg;vertical++) {
302 for (
int horizontal=0;horizontal<=deg;horizontal++) {
303 ltgMapping(ielem,local_dof_cur) = start + vertical*(NX*deg+1)+horizontal;
312 ofstream ltgout(
"ltg.dat");
313 for (
int j=0; j<NY; j++) {
314 for (
int i=0; i<NX; i++) {
315 int ielem = i + j * NX;
316 for (
int m=0; m<numFieldsG; m++){
317 ltgout << ltgMapping(ielem,m) <<
" ";
326 *outStream <<
"Building stiffness matrix and right hand side ... \n\n";
354 Epetra_SerialComm Comm;
355 Epetra_Map globalMapG(numDOF, 0, Comm);
356 Epetra_FEVector u(globalMapG);
357 Epetra_FEVector Ku(globalMapG);
361 refQuadNodes(0,0,0) = 0.0;
362 refQuadNodes(0,0,1) = 0.0;
363 refQuadNodes(0,1,0) = hx;
364 refQuadNodes(0,1,1) = 0.0;
365 refQuadNodes(0,2,0) = hx;
366 refQuadNodes(0,2,1) = hy;
367 refQuadNodes(0,3,0) = 0.0;
368 refQuadNodes(0,3,1) = hy;
371 CellTools::setJacobian(refQuadJacobian, cubPoints, refQuadNodes, quad_4);
372 CellTools::setJacobianInv(refQuadJacobInv, refQuadJacobian );
373 CellTools::setJacobianDet(refQuadJacobDet, refQuadJacobian );
376 fst::HGRADtransformGRAD<double>(quadGradsTransformed, refQuadJacobInv, quadGrads);
379 fst::computeCellMeasure<double>(weightedMeasure, refQuadJacobDet, cubWeights);
382 fst::multiplyMeasure<double>(quadGradsTransformedWeighted,
383 weightedMeasure, quadGradsTransformed);
386 fst::integrate<double>(localStiffMatrix,
387 quadGradsTransformed, quadGradsTransformedWeighted, COMP_BLAS);
389 Epetra_Time graphTimer(Comm);
390 Epetra_CrsGraph grph( Copy , globalMapG , 4 * numFieldsG );
391 for (
int k=0;k<numElems;k++)
393 for (
int i=0;i<numFieldsG;i++)
395 grph.InsertGlobalIndices(ltgMapping(k,i),numFieldsG,<gMapping(k,0));
399 const double graphTime = graphTimer.ElapsedTime();
400 std::cout <<
"Graph computed in " << graphTime <<
"\n";
402 Epetra_Time instantiateTimer( Comm );
403 Epetra_FECrsMatrix StiffMatrix( Copy , grph );
404 const double instantiateTime = instantiateTimer.ElapsedTime( );
405 std::cout <<
"Matrix instantiated in " << instantiateTime <<
"\n";
407 Epetra_Time assemblyTimer(Comm);
410 for (
int k=0; k<numElems; k++)
413 StiffMatrix.InsertGlobalValues(numFieldsG,<gMapping(k,0),numFieldsG,<gMapping(k,0),&localStiffMatrix(0,0,0));
419 StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
421 double assembleTime = assemblyTimer.ElapsedTime();
422 std::cout <<
"Time to insert reference element matrix into global matrix: " << assembleTime << std::endl;
423 std::cout <<
"Total matrix construction time: " << assembleTime + instantiateTime + graphTime <<
"\n";
424 std::cout <<
"There are " << StiffMatrix.NumGlobalNonzeros() <<
" nonzeros in the matrix.\n";
425 std::cout <<
"There are " << numDOF <<
" global degrees of freedom.\n";
427 Epetra_Time multTimer(Comm);
428 StiffMatrix.Apply(u,Ku);
429 double multTime = multTimer.ElapsedTime();
430 std::cout <<
"Time to apply: " << multTime << std::endl;
439 std::cout <<
"End Result: TEST PASSED\n";
442 std::cout.copyfmt(oldFormatState);
449 double evalu(
double & x,
double & y,
double & z)
457 double exactu = sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
463 int evalGradu(
double & x,
double & y,
double & z,
double & gradu1,
double & gradu2,
double & gradu3)
473 gradu1 = (M_PI*cos(M_PI*x)+sin(M_PI*x))
474 *sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
475 gradu2 = (M_PI*cos(M_PI*y)+sin(M_PI*y))
476 *sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z);
477 gradu3 = (M_PI*cos(M_PI*z)+sin(M_PI*z))
478 *sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z);
484 double evalDivGradu(
double & x,
double & y,
double & z)
492 double divGradu = -3.0*M_PI*M_PI*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
493 + 2.0*M_PI*cos(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z)
494 + 2.0*M_PI*cos(M_PI*y)*sin(M_PI*x)*sin(M_PI*z)*exp(x+y+z)
495 + 2.0*M_PI*cos(M_PI*z)*sin(M_PI*x)*sin(M_PI*y)*exp(x+y+z)
496 + 3.0*sin(M_PI*x)*sin(M_PI*y)*sin(M_PI*z)*exp(x+y+z);
Header file for utility class to provide multidimensional containers.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Header file for the Intrepid::HGRAD_QUAD_Cn_FEM class.
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.