92 #include "Epetra_Time.h"
93 #include "Epetra_Map.h"
94 #include "Epetra_FECrsMatrix.h"
95 #include "Epetra_FEVector.h"
96 #include "Epetra_SerialComm.h"
99 #include "Teuchos_oblackholestream.hpp"
100 #include "Teuchos_RCP.hpp"
101 #include "Teuchos_BLAS.hpp"
104 #include "Shards_CellTopology.hpp"
107 #include "EpetraExt_RowMatrixOut.h"
108 #include "EpetraExt_MultiVectorOut.h"
111 using namespace Intrepid;
114 double evalu(
double & x,
double & y,
double & z);
115 int evalGradu(
double & x,
double & y,
double & z,
double & gradu1,
double & gradu2,
double & gradu3);
116 double evalDivGradu(
double & x,
double & y,
double & z);
118 int main(
int argc,
char *argv[]) {
122 std::cout <<
"\n>>> ERROR: Invalid number of arguments.\n\n";
123 std::cout <<
"Usage:\n\n";
124 std::cout <<
" ./Intrepid_example_Drivers_Example_05.exe deg NX NY verbose\n\n";
125 std::cout <<
" where \n";
126 std::cout <<
" int deg - polynomial degree to be used (assumed > 1) \n";
127 std::cout <<
" int NX - num intervals in x direction (assumed box domain, 0,1) \n";
128 std::cout <<
" int NY - num intervals in y direction (assumed box domain, 0,1) \n";
129 std::cout <<
" verbose (optional) - any character, indicates verbose output \n\n";
135 int iprint = argc - 1;
136 Teuchos::RCP<std::ostream> outStream;
137 Teuchos::oblackholestream bhs;
139 outStream = Teuchos::rcp(&std::cout,
false);
141 outStream = Teuchos::rcp(&bhs,
false);
144 Teuchos::oblackholestream oldFormatState;
145 oldFormatState.copyfmt(std::cout);
148 <<
"===============================================================================\n" \
150 <<
"| Example: Generate Stiffness Matrix and Right Hand Side Vector for |\n" \
151 <<
"| Poisson Equation on Quadrilateral Mesh |\n" \
153 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
154 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
155 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
157 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
158 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
160 <<
"===============================================================================\n";
165 int deg = atoi(argv[1]);
166 int NX = atoi(argv[2]);
167 int NY = atoi(argv[3]);
173 typedef shards::CellTopology CellTopology;
174 CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
177 int numNodesPerElem = quad_4.getNodeCount();
178 int spaceDim = quad_4.getDimension();
182 *outStream <<
"Generating mesh ... \n\n";
184 *outStream <<
" NX" <<
" NY\n";
185 *outStream << std::setw(5) << NX <<
186 std::setw(5) << NY <<
"\n\n";
189 int numElems = NX*NY;
190 int numNodes = (NX+1)*(NY+1);
191 *outStream <<
" Number of Elements: " << numElems <<
" \n";
192 *outStream <<
" Number of Nodes: " << numNodes <<
" \n\n";
195 double leftX = 0.0, rightX = 1.0;
196 double leftY = 0.0, rightY = 1.0;
199 double hx = (rightX-leftX)/((
double)NX);
200 double hy = (rightY-leftY)/((
double)NY);
206 for (
int j=0; j<NY+1; j++) {
207 for (
int i=0; i<NX+1; i++) {
208 nodeCoord(inode,0) = leftX + (double)i*hx;
209 nodeCoord(inode,1) = leftY + (double)j*hy;
210 if (j==0 || i==0 || j==NY || i==NX){
211 nodeOnBoundary(inode)=1;
214 nodeOnBoundary(inode)=0;
222 ofstream fcoordout(
"coords.dat");
223 for (
int i=0; i<numNodes; i++) {
224 fcoordout << nodeCoord(i,0) <<
" ";
225 fcoordout << nodeCoord(i,1) <<
"\n";
235 for (
int j=0; j<NY; j++) {
236 for (
int i=0; i<NX; i++) {
237 elemToNode(ielem,0) = (NX + 1)*j + i;
238 elemToNode(ielem,1) = (NX + 1)*j + i + 1;
239 elemToNode(ielem,2) = (NX + 1)*(j + 1) + i + 1;
240 elemToNode(ielem,3) = (NX + 1)*(j + 1) + i;
246 ofstream fe2nout(
"elem2node.dat");
247 for (
int j=0; j<NY; j++) {
248 for (
int i=0; i<NX; i++) {
249 int ielem = i + j * NX;
250 for (
int m=0; m<numNodesPerElem; m++){
251 fe2nout << elemToNode(ielem,m) <<
" ";
261 *outStream <<
"Getting cubature ... \n\n";
265 int cubDegree = 2*deg;
266 Teuchos::RCP<Cubature<double> > quadCub = cubFactory.
create(quad_4, cubDegree);
268 int cubDim = quadCub->getDimension();
269 int numCubPoints = quadCub->getNumPoints();
274 quadCub->getCubature(cubPoints, cubWeights);
279 *outStream <<
"Getting basis ... \n\n";
283 int numFieldsG = quadHGradBasis.getCardinality();
288 quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
289 quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
293 const int numDOF = (NX*deg+1)*(NY*deg+1);
295 for (
int j=0;j<NY;j++) {
296 for (
int i=0;i<NX;i++) {
297 const int start = deg * j * ( NX * deg + 1 ) + i * deg;
300 for (
int vertical=0;vertical<=deg;vertical++) {
301 for (
int horizontal=0;horizontal<=deg;horizontal++) {
302 ltgMapping(ielem,local_dof_cur) = start + vertical*(NX*deg+1)+horizontal;
311 ofstream ltgout(
"ltg.dat");
312 for (
int j=0; j<NY; j++) {
313 for (
int i=0; i<NX; i++) {
314 int ielem = i + j * NX;
315 for (
int m=0; m<numFieldsG; m++){
316 ltgout << ltgMapping(ielem,m) <<
" ";
325 *outStream <<
"Building stiffness matrix and right hand side ... \n\n";
330 int numCells = numElems;
348 Epetra_SerialComm Comm;
349 Epetra_Map globalMapG(numDOF, 0, Comm);
351 Epetra_Time graphTimer(Comm);
352 Epetra_CrsGraph grph( Copy , globalMapG , 4 * numFieldsG );
353 for (
int k=0;k<numElems;k++)
355 for (
int i=0;i<numFieldsG;i++)
357 grph.InsertGlobalIndices(ltgMapping(k,i),numFieldsG,<gMapping(k,0));
362 const double graphTime = graphTimer.ElapsedTime();
363 std::cout <<
"Graph computed in " << graphTime <<
"\n";
365 Epetra_Time instantiateTimer( Comm );
366 Epetra_FECrsMatrix StiffMatrix( Copy , grph );
367 const double instantiateTime = instantiateTimer.ElapsedTime( );
368 std::cout <<
"Matrix instantiated in " << instantiateTime <<
"\n";
370 Epetra_FEVector u(globalMapG);
371 Epetra_FEVector Ku(globalMapG);
378 for (
int i=0;i<numElems;i++)
380 for (
int j=0;j<4;j++)
382 const int nodeCur = elemToNode(i,j);
383 for (
int k=0;k<spaceDim;k++)
385 cellVertices(i,j,k) = nodeCoord(nodeCur,k);
390 Epetra_Time localConstructTimer(Comm);
393 CellTools::setJacobian(cellJacobian, cubPoints, cellVertices, quad_4);
394 CellTools::setJacobianInv(cellJacobInv, cellJacobian );
395 CellTools::setJacobianDet(cellJacobDet, cellJacobian );
398 fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
401 fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
404 fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
405 weightedMeasure, transformedBasisGradients);
408 fst::integrate<double>(localStiffMatrices,
409 transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
411 const double localConstructTime = localConstructTimer.ElapsedTime();
412 std::cout <<
"Time to build local matrices (including Jacobian computation): "<< localConstructTime <<
"\n";
414 Epetra_Time insertionTimer(Comm);
417 for (
int k=0; k<numElems; k++)
420 StiffMatrix.InsertGlobalValues(numFieldsG,<gMapping(k,0),numFieldsG,<gMapping(k,0),&localStiffMatrices(k,0,0));
423 StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
424 const double insertionTime = insertionTimer.ElapsedTime( );
426 std::cout <<
"Time to assemble global matrix from local matrices: " << insertionTime <<
"\n";
436 std::cout <<
"End Result: TEST PASSED\n";
439 std::cout.copyfmt(oldFormatState);
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.