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