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