Intrepid
test_06.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 
44 
52 #include "Intrepid_Utils.hpp"
53 #include "Teuchos_oblackholestream.hpp"
54 #include "Teuchos_RCP.hpp"
55 #include "Teuchos_BLAS.hpp"
56 #include "Teuchos_GlobalMPISession.hpp"
57 
58 using namespace Intrepid;
59 
60 
61 /*
62  Monomial evaluation.
63  in 1D, for point p(x) : x^xDeg
64  in 2D, for point p(x,y) : x^xDeg * y^yDeg
65  in 3D, for point p(x,y,z): x^xDeg * y^yDeg * z^zDeg
66 */
67 double computeMonomial(FieldContainer<double> & p, int xDeg, int yDeg=0, int zDeg=0) {
68  double val = 1.0;
69  int polydeg[3];
70  polydeg[0] = xDeg; polydeg[1] = yDeg; polydeg[2] = zDeg;
71  for (int i=0; i<p.dimension(0); i++) {
72  val *= std::pow(p(i),polydeg[i]);
73  }
74  return val;
75 }
76 
77 
78 /*
79  Computes integrals of monomials over a given reference cell.
80 */
81 void computeIntegral(Teuchos::Array<double>& testIntFixDeg, shards::CellTopology & cellTopology, int cubDegree) {
82 
83  DefaultCubatureFactory<double> cubFactory; // create factory
84  Teuchos::RCP<Cubature<double> > myCub = cubFactory.create(cellTopology, cubDegree); // create default cubature
85 
86  int cubDim = myCub->getDimension();
87  int numCubPoints = myCub->getNumPoints();
88  int numPolys = (cubDegree+1)*(cubDegree+2)*(cubDegree+3)/6;
89 
90  FieldContainer<double> point(cubDim);
91  FieldContainer<double> cubPoints(numCubPoints, cubDim);
92  FieldContainer<double> cubWeights(numCubPoints);
93  FieldContainer<double> functValues(numCubPoints, numPolys);
94 
95  myCub->getCubature(cubPoints, cubWeights);
96 
97  int polyCt = 0;
98  for (int xDeg=0; xDeg <= cubDegree; xDeg++) {
99  for (int yDeg=0; yDeg <= cubDegree-xDeg; yDeg++) {
100  for (int zDeg=0; zDeg <= cubDegree-xDeg-yDeg; zDeg++) {
101  for (int i=0; i<numCubPoints; i++) {
102  for (int j=0; j<cubDim; j++) {
103  point(j) = cubPoints(i,j);
104  }
105  functValues(i,polyCt) = computeMonomial(point, xDeg, yDeg, zDeg);
106  }
107  polyCt++;
108  }
109  }
110  }
111 
112  Teuchos::BLAS<int, double> myblas;
113  int inc = 1;
114  double alpha = 1.0;
115  double beta = 0.0;
116  myblas.GEMV(Teuchos::NO_TRANS, numPolys, numCubPoints, alpha, &functValues[0], numPolys,
117  &cubWeights[0], inc, beta, &testIntFixDeg[0], inc);
118 }
119 
120 
121 int main(int argc, char *argv[]) {
122 
123  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
124 
125  // This little trick lets us print to std::cout only if
126  // a (dummy) command-line argument is provided.
127  int iprint = argc - 1;
128  Teuchos::RCP<std::ostream> outStream;
129  Teuchos::oblackholestream bhs; // outputs nothing
130  if (iprint > 0)
131  outStream = Teuchos::rcp(&std::cout, false);
132  else
133  outStream = Teuchos::rcp(&bhs, false);
134 
135  // Save the format state of the original std::cout.
136  Teuchos::oblackholestream oldFormatState;
137  oldFormatState.copyfmt(std::cout);
138 
139  *outStream \
140  << "===============================================================================\n" \
141  << "| |\n" \
142  << "| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
143  << "| |\n" \
144  << "| 1) Computing integrals of monomials on reference cells in 3D |\n" \
145  << "| - using Level 2 BLAS - |\n" \
146  << "| |\n" \
147  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
148  << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
149  << "| |\n" \
150  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
151  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
152  << "| |\n" \
153  << "===============================================================================\n"\
154  << "| TEST 1: integrals of monomials in 3D (Level 2 BLAS version) |\n"\
155  << "===============================================================================\n";
156 
157  // internal variables:
158  int errorFlag = 0;
159  int polyCt = 0;
160  int offset = 0;
161  Teuchos::Array< Teuchos::Array<double> > testInt;
162  Teuchos::Array< Teuchos::Array<double> > analyticInt;
163  Teuchos::Array<double> tmparray(1);
164  double reltol = 1.0e+04 * INTREPID_TOL;
165  int maxDeg[4];
166  int maxOffset[4];
167  int numPoly[4];
168  int numAnalytic[4];
169  // max polynomial degree tested, per cell type:
171  maxDeg[1] = 20; // can be as large as INTREPID_CUBATURE_LINE_GAUSS_MAX, but runtime is excessive
174  // max polynomial degree recorded in analytic comparison files, per cell type:
175  maxOffset[0] = INTREPID_CUBATURE_TET_DEFAULT_MAX;
176  maxOffset[1] = INTREPID_CUBATURE_LINE_GAUSS_MAX;
179  for (int i=0; i<4; i++) {
180  numPoly[i] = (maxDeg[i]+1)*(maxDeg[i]+2)*(maxDeg[i]+3)/6;
181  }
182  for (int i=0; i<4; i++) {
183  numAnalytic[i] = (maxOffset[i]+1)*(maxOffset[i]+2)*(maxOffset[i]+3)/6;
184  }
185 
186 
187  // get names of files with analytic values
188  std::string basedir = "./data";
189  std::stringstream namestream[4];
190  std::string filename[4];
191  namestream[0] << basedir << "/TET_integrals" << ".dat";
192  namestream[0] >> filename[0];
193  namestream[1] << basedir << "/HEX_integrals" << ".dat";
194  namestream[1] >> filename[1];
195  namestream[2] << basedir << "/TRIPRISM_integrals" << ".dat";
196  namestream[2] >> filename[2];
197  namestream[3] << basedir << "/PYR_integrals" << ".dat";
198  namestream[3] >> filename[3];
199 
200  // reference cells tested
201  shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Tetrahedron<> >(),
202  shards::getCellTopologyData< shards::Hexahedron<> >(),
203  shards::getCellTopologyData< shards::Wedge<> >(),
204  shards::getCellTopologyData< shards::Pyramid<> >() };
205  // format of data files with analytic values
206  TypeOfExactData dataFormat[] = {INTREPID_UTILS_SCALAR, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION, INTREPID_UTILS_FRACTION};
207 
208  // compute and compare integrals
209  try {
210  for (int cellCt=0; cellCt < 4; cellCt++) {
211  testInt.assign(numPoly[cellCt], tmparray);
212  analyticInt.assign(numAnalytic[cellCt], tmparray);
213 
214  *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name << ":\n";
215  std::ifstream filecompare(&filename[cellCt][0]);
216  // compute integrals
217  for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
218  int numMonomials = (cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6;
219  testInt[cubDeg].resize(numMonomials);
220  computeIntegral(testInt[cubDeg], cellType[cellCt], cubDeg);
221  }
222  // get analytic values
223  if (filecompare.is_open()) {
224  getAnalytic(analyticInt, filecompare, dataFormat[cellCt]);
225  // close file
226  filecompare.close();
227  }
228  // perform comparison
229  for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
230  polyCt = 0;
231  offset = 0;
232  int oldErrorFlag = errorFlag;
233  for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
234  for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
235  for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
236  double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
237  double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
238  if (absdiff > abstol) {
239  *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
240  << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg
241  << " * z^" << std::setw(2) << zDeg << ":" << " "
242  << std::scientific << std::setprecision(16)
243  << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
244  << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
245  errorFlag++;
246  *outStream << std::right << std::setw(118) << "^^^^---FAILURE!\n";
247  }
248  polyCt++;
249  }
250  offset = offset + maxOffset[cellCt] - cubDeg;
251  }
252  offset = offset + (maxOffset[cellCt] - cubDeg)*(maxOffset[cellCt] - cubDeg + 1)/2;
253  }
254  *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg;
255  if (errorFlag == oldErrorFlag)
256  *outStream << " passed.\n";
257  else
258  *outStream << " failed.\n";
259  }
260  *outStream << "\n";
261  } // end for cellCt
262  }
263  catch (std::logic_error err) {
264  *outStream << err.what() << "\n";
265  errorFlag = -1;
266  };
267 
268 
269  if (errorFlag != 0)
270  std::cout << "End Result: TEST FAILED\n";
271  else
272  std::cout << "End Result: TEST PASSED\n";
273 
274  // reset format state of std::cout
275  std::cout.copyfmt(oldFormatState);
276 
277  return errorFlag;
278 }
int dimension(const int whichDim) const
Returns the specified dimension.
#define INTREPID_CUBATURE_TRI_DEFAULT_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct triangle rule of the ...
#define INTREPID_CUBATURE_LINE_GAUSS_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct line rule of the Gaus...
#define INTREPID_CUBATURE_LINE_GAUSSJACOBI20_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct line rule of the Gaus...
#define INTREPID_CUBATURE_TET_DEFAULT_MAX
The maximum degree of the polynomial that can be integrated exactly by a direct tetrahedron rule of t...
Intrepid utilities.
Header file for the abstract base class Intrepid::DefaultCubatureFactory.
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.