Intrepid
example_12.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 
83 // Intrepid includes
86 #include "Intrepid_CellTools.hpp"
87 //#include "Intrepid_ArrayTools.hpp"
89 //#include "Intrepid_RealSpaceTools.hpp"
91 #include "Intrepid_Utils.hpp"
92 
93 // Epetra includes
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"
99 
100 // Teuchos includes
101 #include "Teuchos_oblackholestream.hpp"
102 #include "Teuchos_RCP.hpp"
103 //#include "Teuchos_BLAS.hpp"
104 //#include "Teuchos_BLAS_types.hpp"
105 
106 // Shards includes
107 #include "Shards_CellTopology.hpp"
108 
109 // EpetraExt includes
110 #include "EpetraExt_MultiVectorOut.h"
111 
112 using namespace std;
113 using namespace Intrepid;
114 
115 int main(int argc, char *argv[]) {
116 
117  //Check number of arguments
118  if (argc < 4) {
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";
128  exit(1);
129  }
130 
131  // This little trick lets us print to std::cout only if
132  // a (dummy) command-line argument is provided.
133  int iprint = argc - 1;
134  Teuchos::RCP<std::ostream> outStream;
135  Teuchos::oblackholestream bhs; // outputs nothing
136  if (iprint > 2)
137  outStream = Teuchos::rcp(&std::cout, false);
138  else
139  outStream = Teuchos::rcp(&bhs, false);
140 
141  // Save the format state of the original std::cout.
142  Teuchos::oblackholestream oldFormatState;
143  oldFormatState.copyfmt(std::cout);
144 
145  *outStream \
146  << "===============================================================================\n" \
147  << "| |\n" \
148  << "| Example: Build Stiffness Matrix for |\n" \
149  << "| Poisson Equation on Hexahedral Mesh |\n" \
150  << "| |\n" \
151  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
152  << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
153  << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
154  << "| |\n" \
155  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
156  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
157  << "| |\n" \
158  << "===============================================================================\n";
159 
160 
161  // ************************************ GET INPUTS **************************************
162 
163  int deg = atoi(argv[1]); // polynomial degree to use
164  int NX = atoi(argv[2]); // num intervals in x direction (assumed box domain, 0,1)
165  int NY = atoi(argv[3]); // num intervals in y direction (assumed box domain, 0,1)
166  int NZ = atoi(argv[4]); // num intervals in y direction (assumed box domain, 0,1)
167 
168 
169  // *********************************** CELL TOPOLOGY **********************************
170 
171  // Get cell topology for base hexahedron
172  typedef shards::CellTopology CellTopology;
173  CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() );
174 
175  // Get dimensions
176  int numNodesPerElem = hex_8.getNodeCount();
177  int spaceDim = hex_8.getDimension();
178 
179  // *********************************** GENERATE MESH ************************************
180 
181  *outStream << "Generating mesh ... \n\n";
182 
183  *outStream << " NX" << " NY" << " NZ\n";
184  *outStream << std::setw(5) << NX <<
185  std::setw(5) << NY << std::setw(5) << NZ << "\n\n";
186 
187  // Print mesh information
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";
192 
193  // Cube
194  double leftX = 0.0, rightX = 1.0;
195  double leftY = 0.0, rightY = 1.0;
196  double leftZ = 0.0, rightZ = 1.0;
197 
198  // Mesh spacing
199  double hx = (rightX-leftX)/((double)NX);
200  double hy = (rightY-leftY)/((double)NY);
201  double hz = (rightZ-leftZ)/((double)NZ);
202 
203  // Get nodal coordinates
204  FieldContainer<double> nodeCoord(numNodes, spaceDim);
205  FieldContainer<int> nodeOnBoundary(numNodes);
206  int inode = 0;
207  for (int k=0; k<NZ+1; k++)
208  {
209  for (int j=0; j<NY+1; j++)
210  {
211  for (int i=0; i<NX+1; i++)
212  {
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)
217  {
218  nodeOnBoundary(inode)=1;
219  }
220  else
221  {
222  nodeOnBoundary(inode)=0;
223  }
224  inode++;
225  }
226  }
227  }
228 #define DUMP_DATA
229 #ifdef DUMP_DATA
230  // Print nodal coords
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";
236  }
237  fcoordout.close();
238 #endif
239 
240 
241  // Element to Node map
242  // We'll keep it around, but this is only the DOFMap if you are in the lowest order case.
243  FieldContainer<int> elemToNode(numElems, numNodesPerElem);
244  int ielem = 0;
245  for (int k=0; k<NZ; k++)
246  {
247  for (int j=0; j<NY; j++)
248  {
249  for (int i=0; i<NX; i++)
250  {
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;
259  ielem++;
260  }
261  }
262  }
263 #ifdef DUMP_DATA
264  // Output connectivity
265  ofstream fe2nout("elem2node.dat");
266  for (int k=0;k<NZ;k++)
267  {
268  for (int j=0; j<NY; j++)
269  {
270  for (int i=0; i<NX; i++)
271  {
272  ielem = i + j * NX + k * NY * NY;
273  for (int m=0; m<numNodesPerElem; m++)
274  {
275  fe2nout << elemToNode(ielem,m) <<" ";
276  }
277  fe2nout <<"\n";
278  }
279  }
280  }
281  fe2nout.close();
282 #endif
283 
284  // ************************************ CUBATURE **************************************
285  *outStream << "Getting cubature ... \n\n";
286 
287  // Get numerical integration points and weights
288  DefaultCubatureFactory<double> cubFactory;
289  int cubDegree = 2*deg;
290  Teuchos::RCP<Cubature<double> > quadCub = cubFactory.create(hex_8, cubDegree);
291 
292  int cubDim = quadCub->getDimension();
293  int numCubPoints = quadCub->getNumPoints();
294 
295  FieldContainer<double> cubPoints(numCubPoints, cubDim);
296  FieldContainer<double> cubWeights(numCubPoints);
297 
298  quadCub->getCubature(cubPoints, cubWeights);
299 
300 
301  // ************************************** BASIS ***************************************
302 
303  *outStream << "Getting basis ... \n\n";
304 
305  // Define basis
306  Basis_HGRAD_HEX_Cn_FEM<double, FieldContainer<double> > quadHGradBasis(deg,POINTTYPE_SPECTRAL);
307  int numFieldsG = quadHGradBasis.getCardinality();
308  FieldContainer<double> quadGVals(numFieldsG, numCubPoints);
309  FieldContainer<double> quadGrads(numFieldsG, numCubPoints, spaceDim);
310 
311  // Evaluate basis values and gradients at cubature points
312  quadHGradBasis.getValues(quadGVals, cubPoints, OPERATOR_VALUE);
313  quadHGradBasis.getValues(quadGrads, cubPoints, OPERATOR_GRAD);
314 
315  // create the local-global mapping
316  FieldContainer<int> ltgMapping(numElems,numFieldsG);
317  const int numDOF = (NX*deg+1)*(NY*deg+1)*(NZ*deg+1);
318  ielem=0;
319  for (int k=0;k<NZ;k++)
320  {
321  for (int j=0;j<NY;j++)
322  {
323  for (int i=0;i<NX;i++)
324  {
325  const int start = k * ( NY * deg + 1 ) * ( NX * deg + 1 ) + j * ( NX * deg + 1 ) + i * deg;
326  // loop over local dof on this cell
327  int local_dof_cur=0;
328  for (int kloc=0;kloc<=deg;kloc++)
329  {
330  for (int jloc=0;jloc<=deg;jloc++)
331  {
332  for (int iloc=0;iloc<=deg;iloc++)
333  {
334  ltgMapping(ielem,local_dof_cur) = start
335  + kloc * ( NX * deg + 1 ) * ( NY * deg + 1 )
336  + jloc * ( NX * deg + 1 )
337  + iloc;
338  local_dof_cur++;
339  }
340  }
341  }
342  ielem++;
343  }
344  }
345  }
346 #ifdef DUMP_DATA
347  // Output ltg mapping
348  ielem = 0;
349  ofstream ltgout("ltg.dat");
350  for (int k=0;k<NZ;k++)
351  {
352  for (int j=0; j<NY; j++)
353  {
354  for (int i=0; i<NX; i++)
355  {
356  ielem = i + j * NX + k * NX * NY;
357  for (int m=0; m<numFieldsG; m++)
358  {
359  ltgout << ltgMapping(ielem,m) <<" ";
360  }
361  ltgout <<"\n";
362  }
363  }
364  }
365  ltgout.close();
366 #endif
367 
368  // ********** DECLARE GLOBAL OBJECTS *************
369  Epetra_SerialComm Comm;
370  Epetra_Map globalMapG(numDOF, 0, Comm);
371  Epetra_FEVector u(globalMapG); u.Random();
372  Epetra_FEVector Ku(globalMapG);
373 
374  // Let's preallocate the graph before we instantiate the matrix
375  Epetra_Time graphTimer(Comm);
376  Epetra_CrsGraph grph( Copy , globalMapG , 4 * numFieldsG );
377  for (int k=0;k<numElems;k++)
378  {
379  for (int i=0;i<numFieldsG;i++)
380  {
381  grph.InsertGlobalIndices(ltgMapping(k,i),numFieldsG,&ltgMapping(k,0));
382  }
383  }
384  grph.FillComplete();
385  const double graphTime = graphTimer.ElapsedTime();
386 
387 
388  // time the instantiation
389  Epetra_Time instantiateTimer(Comm);
390  Epetra_FECrsMatrix StiffMatrix(Copy,grph);
391  const double instantiateTime = instantiateTimer.ElapsedTime();
392 
393 
394  // ********** CONSTRUCT AND INSERT LOCAL STIFFNESS MATRICES ***********
395  *outStream << "Building local stiffness matrices...\n\n";
397  typedef FunctionSpaceTools fst;
398 // int numCells = numElems;
399 
400 
401  // jacobian information
402  FieldContainer<double> refCellNodes(1,numNodesPerElem,spaceDim);
403  FieldContainer<double> cellJacobian(1,numCubPoints,spaceDim,spaceDim);
404  FieldContainer<double> cellJacobInv(1,numCubPoints,spaceDim,spaceDim);
405  FieldContainer<double> cellJacobDet(1,numCubPoints);
406 
407  // element stiffness matrices and supporting storage space
408  FieldContainer<double> localStiffMatrix(1, numFieldsG, numFieldsG);
409  FieldContainer<double> transformedBasisGradients(1,numFieldsG,numCubPoints,spaceDim);
410  FieldContainer<double> weightedTransformedBasisGradients(1,numFieldsG,numCubPoints,spaceDim);
411  FieldContainer<double> weightedMeasure(1, numCubPoints);
412 
413 
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;
423 
424 
425  // jacobian evaluation
426  CellTools::setJacobian(cellJacobian,cubPoints,refCellNodes,hex_8);
427  CellTools::setJacobianInv(cellJacobInv, cellJacobian );
428  CellTools::setJacobianDet(cellJacobDet, cellJacobian );
429 
430  // transform reference element gradients to each cell
431  fst::HGRADtransformGRAD<double>(transformedBasisGradients, cellJacobInv, quadGrads);
432 
433  // compute weighted measure
434  fst::computeCellMeasure<double>(weightedMeasure, cellJacobDet, cubWeights);
435 
436  // multiply values with weighted measure
437  fst::multiplyMeasure<double>(weightedTransformedBasisGradients,
438  weightedMeasure, transformedBasisGradients);
439 
440  // integrate to compute element stiffness matrix
441  fst::integrate<double>(localStiffMatrix,
442  transformedBasisGradients, weightedTransformedBasisGradients , COMP_BLAS);
443 
444  const double localConstructTime = localConstructTimer.ElapsedTime();
445 
446 
447  Epetra_Time insertionTimer(Comm);
448 
449  // *** Element loop ***
450  for (int k=0; k<numElems; k++)
451  {
452  // assemble into global matrix
453  StiffMatrix.InsertGlobalValues(numFieldsG,&ltgMapping(k,0),numFieldsG,&ltgMapping(k,0),&localStiffMatrix(0,0,0));
454 
455  }
456  StiffMatrix.GlobalAssemble(); StiffMatrix.FillComplete();
457  const double insertionTime = insertionTimer.ElapsedTime( );
458 
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";
464 
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";
469 
470  *outStream << "End Result: TEST PASSED\n";
471 
472  // reset format state of std::cout
473  std::cout.copyfmt(oldFormatState);
474 
475  return 0;
476 }
477 
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: