Intrepid
test_09.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, int cubDegree) {
82 
83  CubatureGenSparse<double,3> myCub(cubDegree);
84 
85  int cubDim = myCub.getDimension();
86  int numCubPoints = myCub.getNumPoints();
87  int numPolys = (cubDegree+1)*(cubDegree+2)*(cubDegree+3)/6;
88 
89  FieldContainer<double> point(cubDim);
90  FieldContainer<double> cubPoints(numCubPoints, cubDim);
91  FieldContainer<double> cubWeights(numCubPoints);
92  FieldContainer<double> functValues(numCubPoints, numPolys);
93 
94  myCub.getCubature(cubPoints, cubWeights);
95 
96  int polyCt = 0;
97  for (int xDeg=0; xDeg <= cubDegree; xDeg++) {
98  for (int yDeg=0; yDeg <= cubDegree-xDeg; yDeg++) {
99  for (int zDeg=0; zDeg <= cubDegree-xDeg-yDeg; zDeg++) {
100  for (int i=0; i<numCubPoints; i++) {
101  for (int j=0; j<cubDim; j++) {
102  point(j) = cubPoints(i,j);
103  }
104  functValues(i,polyCt) = computeMonomial(point, xDeg, yDeg, zDeg);
105  }
106  polyCt++;
107  }
108  }
109  }
110 
111  Teuchos::BLAS<int, double> myblas;
112  int inc = 1;
113  double alpha = 1.0;
114  double beta = 0.0;
115  myblas.GEMV(Teuchos::NO_TRANS, numPolys, numCubPoints, alpha, &functValues[0], numPolys,
116  &cubWeights[0], inc, beta, &testIntFixDeg[0], inc);
117 }
118 
119 
120 int main(int argc, char *argv[]) {
121 
122  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
123 
124  // This little trick lets us print to std::cout only if
125  // a (dummy) command-line argument is provided.
126  int iprint = argc - 1;
127  Teuchos::RCP<std::ostream> outStream;
128  Teuchos::oblackholestream bhs; // outputs nothing
129  if (iprint > 0)
130  outStream = Teuchos::rcp(&std::cout, false);
131  else
132  outStream = Teuchos::rcp(&bhs, false);
133 
134  // Save the format state of the original std::cout.
135  Teuchos::oblackholestream oldFormatState;
136  oldFormatState.copyfmt(std::cout);
137 
138  *outStream \
139  << "===============================================================================\n" \
140  << "| |\n" \
141  << "| Unit Test (CubatureGenSparse) |\n" \
142  << "| |\n" \
143  << "| 1) Computing integrals of monomials on reference cells in 3D |\n" \
144  << "| - using Level 2 BLAS - |\n" \
145  << "| |\n" \
146  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
147  << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
148  << "| Matthew Keegan (mskeega@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) using |\n"\
155  << "| Generalized Sparse Grid Construction |\n"\
156  << "===============================================================================\n";
157 
158  // internal variables:
159  int errorFlag = 0;
160  int polyCt = 0;
161  int offset = 0;
162  Teuchos::Array< Teuchos::Array<double> > testInt;
163  Teuchos::Array< Teuchos::Array<double> > analyticInt;
164  Teuchos::Array<double> tmparray(1);
165  double reltol = 1.0e+04 * INTREPID_TOL;
166  int maxDeg = 20; // can be as large as INTREPID_CUBATURE_SPARSE3D_GAUSS_MAX, but runtime is excessive
167  int maxOffset = INTREPID_CUBATURE_LINE_GAUSS_MAX;
168  int numPoly = (maxDeg+1)*(maxDeg+2)*(maxDeg+3)/6;
169  int numAnalytic = (maxOffset+1)*(maxOffset+2)*(maxOffset+3)/6;
170  testInt.assign(numPoly, tmparray);
171  analyticInt.assign(numAnalytic, tmparray);
172 
173  // get names of files with analytic values
174  std::string basedir = "./data";
175  std::stringstream namestream;
176  std::string filename;
177  namestream << basedir << "/HEX_integrals" << ".dat";
178  namestream >> filename;
179 
180  // format of data files with analytic values
181  TypeOfExactData dataFormat = INTREPID_UTILS_FRACTION;
182 
183  // compute and compare integrals
184  try {
185  *outStream << "\nIntegrals of monomials:\n";
186  std::ifstream filecompare(&filename[0]);
187  // compute integrals
188  for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
189  int numMonomials = (cubDeg+1)*(cubDeg+2)*(cubDeg+3)/6;
190  testInt[cubDeg].resize(numMonomials);
191  computeIntegral(testInt[cubDeg], cubDeg);
192  }
193  // get analytic values
194  if (filecompare.is_open()) {
195  getAnalytic(analyticInt, filecompare, dataFormat);
196  // close file
197  filecompare.close();
198  }
199  // perform comparison
200  for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
201  polyCt = 0;
202  offset = 0;
203  int oldErrorFlag = errorFlag;
204  for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
205  for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
206  for (int zDeg=0; zDeg <= cubDeg-xDeg-yDeg; zDeg++) {
207  double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
208  double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
209  if (absdiff > abstol) {
210  *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
211  << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg
212  << " * z^" << std::setw(2) << zDeg << ":" << " "
213  << std::scientific << std::setprecision(16)
214  << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
215  << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
216  errorFlag++;
217  *outStream << std::right << std::setw(118) << "^^^^---FAILURE!\n";
218  }
219  polyCt++;
220  }
221  offset = offset + maxOffset - cubDeg;
222  }
223  offset = offset + (maxOffset - cubDeg)*(maxOffset - cubDeg + 1)/2;
224  }
225  *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg;
226  if (errorFlag == oldErrorFlag)
227  *outStream << " passed.\n";
228  else
229  *outStream << " failed.\n";
230  }
231  *outStream << "\n";
232  }
233  catch (const std::logic_error & err) {
234  *outStream << err.what() << "\n";
235  errorFlag = -1;
236  };
237 
238 
239  if (errorFlag != 0)
240  std::cout << "End Result: TEST FAILED\n";
241  else
242  std::cout << "End Result: TEST PASSED\n";
243 
244  // reset format state of std::cout
245  std::cout.copyfmt(oldFormatState);
246  return errorFlag;
247 }
Header file for the Intrepid::CubatureGenSparse class.
int dimension(const int whichDim) const
Returns the specified dimension.
#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...
Intrepid utilities.