Intrepid
example_10.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
84 // Intrepid includes
87 #include "Intrepid_CellTools.hpp"
88 //#include "Intrepid_ArrayTools.hpp"
90 //#include "Intrepid_RealSpaceTools.hpp"
92 #include "Intrepid_Utils.hpp"
93 
94 // Epetra includes
95 #include "Epetra_Time.h"
96 #include "Epetra_Map.h"
97 #include "Epetra_FEVector.h"
98 #include "Epetra_FECrsMatrix.h"
99 #include "Epetra_SerialComm.h"
100 
101 // Teuchos includes
102 #include "Teuchos_oblackholestream.hpp"
103 #include "Teuchos_RCP.hpp"
104 //#include "Teuchos_BLAS.hpp"
105 //#include "Teuchos_BLAS_types.hpp"
106 
107 // Shards includes
108 #include "Shards_CellTopology.hpp"
109 
110 // EpetraExt includes
111 #include "EpetraExt_MultiVectorOut.h"
112 
113 using namespace std;
114 using namespace Intrepid;
115 
116 int main(int argc, char *argv[]) {
117 
118  //Check number of arguments
119  if (argc < 4) {
120  std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n";
121  std::cout <<"Usage:\n\n";
122  std::cout <<" ./Intrepid_example_Drivers_Example_10.exe deg NX NY NZ verbose\n\n";
123  std::cout <<" where \n";
124  std::cout <<" int deg - polynomial degree to be used (assumed >= 1) \n";
125  std::cout <<" int NX - num intervals in x direction (assumed box domain, 0,1) \n";
126  std::cout <<" int NY - num intervals in y direction (assumed box domain, 0,1) \n";
127  std::cout <<" int NZ - num intervals in y direction (assumed box domain, 0,1) \n";
128  std::cout <<" verbose (optional) - any character, indicates verbose output \n\n";
129  exit(1);
130  }
131 
132  // This little trick lets us print to std::cout only if
133  // a (dummy) command-line argument is provided.
134  int iprint = argc - 1;
135  Teuchos::RCP<std::ostream> outStream;
136  Teuchos::oblackholestream bhs; // outputs nothing
137  if (iprint > 2)
138  outStream = Teuchos::rcp(&std::cout, false);
139  else
140  outStream = Teuchos::rcp(&bhs, false);
141 
142  // Save the format state of the original std::cout.
143  Teuchos::oblackholestream oldFormatState;
144  oldFormatState.copyfmt(std::cout);
145 
146  *outStream \
147  << "===============================================================================\n" \
148  << "| |\n" \
149  << "| Example: Build Stiffness Matrix for |\n" \
150  << "| Poisson Equation on Hexahedral Mesh |\n" \
151  << "| |\n" \
152  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
153  << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
154  << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
155  << "| |\n" \
156  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
157  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
158  << "| |\n" \
159  << "===============================================================================\n";
160 
161 
162  // ************************************ GET INPUTS **************************************
163 
164  int deg = atoi(argv[1]); // polynomial degree to use
165  int NX = atoi(argv[2]); // num intervals in x direction (assumed box domain, 0,1)
166  int NY = atoi(argv[3]); // num intervals in y direction (assumed box domain, 0,1)
167  int NZ = atoi(argv[4]); // num intervals in y direction (assumed box domain, 0,1)
168 
169 
170  // *********************************** CELL TOPOLOGY **********************************
171 
172  // Get cell topology for base hexahedron
173  typedef shards::CellTopology CellTopology;
174  CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
175 
176  // Get dimensions
177  int numNodesPerElem = hex_8.getNodeCount();
178  int spaceDim = hex_8.getDimension();
179 
180  // *********************************** GENERATE MESH ************************************
181 
182  *outStream << "Generating mesh ... \n\n";
183 
184  *outStream << " NX" << " NY" << " NZ\n";
185  *outStream << std::setw(5) << NX <<
186  std::setw(5) << NY << std::setw(5) << NZ << "\n\n";
187 
188  // Print mesh information
189  int numElems = NX*NY*NZ;
190  int numNodes = (NX+1)*(NY+1)*(NZ+1);
191  *outStream << " Number of Elements: " << numElems << " \n";
192  *outStream << " Number of Nodes: " << numNodes << " \n\n";
193 
194  // Cube
195  double leftX = 0.0, rightX = 1.0;
196  double leftY = 0.0, rightY = 1.0;
197  double leftZ = 0.0, rightZ = 1.0;
198 
199  // Mesh spacing
200  double hx = (rightX-leftX)/((double)NX);
201  double hy = (rightY-leftY)/((double)NY);
202  double hz = (rightZ-leftZ)/((double)NZ);
203 
204  // Get nodal coordinates
205  FieldContainer<double> nodeCoord(numNodes, spaceDim);
206  FieldContainer<int> nodeOnBoundary(numNodes);
207  int inode = 0;
208  for (int k=0; k<NZ+1; k++)
209  {
210  for (int j=0; j<NY+1; j++)
211  {
212  for (int i=0; i<NX+1; i++)
213  {
214  nodeCoord(inode,0) = leftX + (double)i*hx;
215  nodeCoord(inode,1) = leftY + (double)j*hy;
216  nodeCoord(inode,2) = leftZ + (double)k*hz;
217  if (k==0 || k==NZ || j==0 || i==0 || j==NY || i==NX)
218  {
219  nodeOnBoundary(inode)=1;
220  }
221  else
222  {
223  nodeOnBoundary(inode)=0;
224  }
225  inode++;
226  }
227  }
228  }
229 #define DUMP_DATA
230 #ifdef DUMP_DATA
231  // Print nodal coords
232  ofstream fcoordout("coords.dat");
233  for (int i=0; i<numNodes; i++) {
234  fcoordout << nodeCoord(i,0) <<" ";
235  fcoordout << nodeCoord(i,1) <<" ";
236  fcoordout << nodeCoord(i,2) <<"\n";
237  }
238  fcoordout.close();
239 #endif
240 
241 
242  // Element to Node map
243  // We'll keep it around, but this is only the DOFMap if you are in the lowest order case.
244  FieldContainer<int> elemToNode(numElems, numNodesPerElem);
245  int ielem = 0;
246  for (int k=0; k<NZ; k++)
247  {
248  for (int j=0; j<NY; j++)
249  {
250  for (int i=0; i<NX; i++)
251  {
252  elemToNode(ielem,0) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
253  elemToNode(ielem,1) = k * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
254  elemToNode(ielem,2) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
255  elemToNode(ielem,3) = k * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
256  elemToNode(ielem,4) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i;
257  elemToNode(ielem,5) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + j * ( NX + 1 ) + i + 1;
258  elemToNode(ielem,6) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i + 1;
259  elemToNode(ielem,7) = ( k + 1 ) * ( NX + 1 ) * ( NY + 1 ) + ( j + 1 ) * ( NX + 1 ) + i;
260  ielem++;
261  }
262  }
263  }
264 #ifdef DUMP_DATA
265  // Output connectivity
266  ofstream fe2nout("elem2node.dat");
267  for (int k=0;k<NZ;k++)
268  {
269  for (int j=0; j<NY; j++)
270  {
271  for (int i=0; i<NX; i++)
272  {
273  int ielem = i + j * NX + k * NY * NY;
274  for (int m=0; m<numNodesPerElem; m++)
275  {
276  fe2nout << elemToNode(ielem,m) <<" ";
277  }
278  fe2nout <<"\n";
279  }
280  }
281  }
282  fe2nout.close();
283 #endif
284 
285  // ************************************ CUBATURE **************************************
286  *outStream << "Getting cubature ... \n\n";
287 
288  // Get numerical integration points and weights
289  DefaultCubatureFactory<double> cubFactory;
290  int cubDegree = 2*deg;
291  Teuchos::RCP<Cubature<double> > quadCub = cubFactory.create(hex_8, cubDegree);
292 
293  int cubDim = quadCub->getDimension();
294  int numCubPoints = quadCub->getNumPoints();
295 
296  FieldContainer<double> cubPoints(numCubPoints, cubDim);
297  FieldContainer<double> cubWeights(numCubPoints);
298 
299  quadCub->getCubature(cubPoints, cubWeights);
300 
301 
302  // ************************************** BASIS ***************************************
303 
304  *outStream << "Getting basis ... \n\n";
305 
306  // Define basis
307  Basis_HGRAD_HEX_Cn_FEM<double, FieldContainer<double> > quadHGradBasis(deg,POINTTYPE_SPECTRAL);
308  int numFieldsG = quadHGradBasis.getCardinality();
309  FieldContainer<double> quadGVals(numFieldsG, numCubPoints);
310  FieldContainer<double> quadGrads(numFieldsG, numCubPoints, spaceDim);
311 
312  // Evaluate basis values and gradients at cubature points
313  quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
314  quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
315 
316  // create the local-global mapping
317  FieldContainer<int> ltgMapping(numElems,numFieldsG);
318  const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
319  ielem=0;
320  for (int k=0;k<NZ;k++)
321  {
322  for (int j=0;j<NY;j++)
323  {
324  for (int i=0;i<NX;i++)
325  {
326  const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
327  // loop over local dof on this cell
328  int local_dof_cur=0;
329  for (int kloc=0;kloc<=deg;kloc++)
330  {
331  for (int jloc=0;jloc<=deg;jloc++)
332  {
333  for (int iloc=0;iloc<=deg;iloc++)
334  {
335  ltgMapping(ielem,local_dof_cur) = start
336  + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
337  + jloc * ( NX * deg + 1 )
338  + iloc;
339  local_dof_cur++;
340  }
341  }
342  }
343  ielem++;
344  }
345  }
346  }
347 #ifdef DUMP_DATA
348  // Output ltg mapping
349  ielem = 0;
350  ofstream ltgout("ltg.dat");
351  for (int k=0;k<NZ;k++)
352  {
353  for (int j=0; j<NY; j++)
354  {
355  for (int i=0; i<NX; i++)
356  {
357  int ielem = i + j * NX + k * NX * NY;
358  for (int m=0; m<numFieldsG; m++)
359  {
360  ltgout << ltgMapping(ielem,m) <<" ";
361  }
362  ltgout <<"\n";
363  }
364  }
365  }
366  ltgout.close();
367 #endif
368 
369  // ********** DECLARE GLOBAL OBJECTS *************
370  Epetra_SerialComm Comm;
371  Epetra_Map globalMapG(numDOF, 0, Comm);
372  Epetra_FEVector u(globalMapG); u.Random();
373  Epetra_FEVector Ku(globalMapG);
374 
375  // time the instantiation
376  Epetra_Time instantiateTimer(Comm);
377  Epetra_FECrsMatrix StiffMatrix(Copy,globalMapG,8*numFieldsG);
378  const double instantiateTime = instantiateTimer.ElapsedTime();
379 
380 
381  // ********** CONSTRUCT AND INSERT LOCAL STIFFNESS MATRICES ***********
382  *outStream << "Building local stiffness matrices...\n\n";
384  typedef FunctionSpaceTools fst;
385  int numCells = numElems;
386 
387  // vertices
388  FieldContainer<double> cellVertices(numCells,numNodesPerElem,spaceDim);
389 
390  // jacobian information
391  FieldContainer<double> cellJacobian(numCells,numCubPoints,spaceDim,spaceDim);
392  FieldContainer<double> cellJacobInv(numCells,numCubPoints,spaceDim,spaceDim);
393  FieldContainer<double> cellJacobDet(numCells,numCubPoints);
394 
395  // element stiffness matrices and supporting storage space
396  FieldContainer<double> localStiffMatrices(numCells, numFieldsG, numFieldsG);
397  FieldContainer<double> transformedBasisGradients(numCells,numFieldsG,numCubPoints,spaceDim);
398  FieldContainer<double> weightedTransformedBasisGradients(numCells,numFieldsG,numCubPoints,spaceDim);
399  FieldContainer<double> weightedMeasure(numCells, numCubPoints);
400 
401 
402  // get vertices of cells (for computing Jacobians)
403  for (int i=0;i<numElems;i++)
404  {
405  for (int j=0;j<numNodesPerElem;j++)
406  {
407  const int nodeCur = elemToNode(i,j);
408  for (int k=0;k<spaceDim;k++)
409  {
410  cellVertices(i,j,k) = nodeCoord(nodeCur,k);
411  }
412  }
413  }
414 
415  Epetra_Time localConstructTimer( Comm );
416 
417  // jacobian evaluation
418  CellTools::setJacobian(cellJacobian,cubPoints,cellVertices,hex_8);
419  CellTools::setJacobianInv(cellJacobInv, cellJacobian );
420  CellTools::setJacobianDet(cellJacobDet, cellJacobian );
421 
422  // transform reference element gradients to each cell
423  fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
424 
425  // compute weighted measure
426  fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
427 
428  // multiply values with weighted measure
429  fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
430  weightedMeasure, transformedBasisGradients);
431 
432  // integrate to compute element stiffness matrix
433  fst::integrate<double>(localStiffMatrices,
434  transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
435 
436  const double localConstructTime = localConstructTimer.ElapsedTime();
437 
438 
439  Epetra_Time insertionTimer(Comm);
440 
441  // *** Element loop ***
442  for (int k=0; k<numElems; k++)
443  {
444  // assemble into global matrix
445  StiffMatrix.InsertGlobalValues(numFieldsG,&ltgMapping(k,0),numFieldsG,&ltgMapping(k,0),&localStiffMatrices(k,0,0));
446 
447  }
448  StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
449  const double insertionTime = insertionTimer.ElapsedTime( );
450 
451  *outStream << "Time to instantiate global stiffness matrix: " << instantiateTime << "\n";
452  *outStream << "Time to build local matrices (including Jacobian computation): "<< localConstructTime << "\n";
453  *outStream << "Time to assemble global matrix from local matrices: " << insertionTime << "\n";
454  *outStream << "Total construction time: " << instantiateTime + localConstructTime + insertionTime << "\n";
455 
456  Epetra_Time applyTimer(Comm);
457  StiffMatrix.Apply(u,Ku);
458  const double multTime = applyTimer.ElapsedTime();
459  *outStream << "Time to multiply onto a vector: " << multTime << "\n";
460 
461  *outStream << "End Result: TEST PASSED\n";
462 
463  // reset format state of std::cout
464  std::cout.copyfmt(oldFormatState);
465 
466  return 0;
467 }
468 
Header file for the Intrepid::CellTools class.
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_HEX_Cn_FEM class.
Defines expert-level interfaces for the evaluation of functions and operators in physical space (supp...
Intrepid utilities.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Header file for the Intrepid::FunctionSpaceTools 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 > &degree)
Factory method.
A stateless class for operations on cell data. Provides methods for: