Intrepid
test_03.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 /*
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) {
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 
93  myCub->getCubature(cubPoints, cubWeights);
94 
95  for (int i=0; i<numCubPoints; i++) {
96  for (int j=0; j<cubDim; j++) {
97  point(j) = cubPoints(i,j);
98  }
99  val += computeMonomial(point, xDeg, yDeg)*cubWeights(i);
100  }
101 
102  return val;
103 }
104 
105 
106 int main(int argc, char *argv[]) {
107 
108  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
109 
110  // This little trick lets us print to std::cout only if
111  // a (dummy) command-line argument is provided.
112  int iprint = argc - 1;
113  Teuchos::RCP<std::ostream> outStream;
114  Teuchos::oblackholestream bhs; // outputs nothing
115  if (iprint > 0)
116  outStream = Teuchos::rcp(&std::cout, false);
117  else
118  outStream = Teuchos::rcp(&bhs, false);
119 
120  // Save the format state of the original std::cout.
121  Teuchos::oblackholestream oldFormatState;
122  oldFormatState.copyfmt(std::cout);
123 
124  *outStream \
125  << "===============================================================================\n" \
126  << "| |\n" \
127  << "| Unit Test (CubatureDirect,CubatureTensor,DefaultCubatureFactory) |\n" \
128  << "| |\n" \
129  << "| 1) Computing integrals of monomials on reference cells in 2D |\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 2D |\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+03 * INTREPID_TOL;
149  int maxDeg[2];
150  int numPoly[2];
155 
156  // get names of files with analytic values
157  std::string basedir = "./data";
158  std::stringstream namestream[2];
159  std::string filename[2];
160  namestream[0] << basedir << "/TRI_integrals" << ".dat";
161  namestream[0] >> filename[0];
162  namestream[1] << basedir << "/QUAD_integrals" << ".dat";
163  namestream[1] >> filename[1];
164 
165  shards::CellTopology cellType[] = {shards::getCellTopologyData< shards::Triangle<> >(),
166  shards::getCellTopologyData< shards::Quadrilateral<> >()};
167 
168  // compute and compare integrals
169  try {
170  for (int cellCt=0; cellCt < 2; cellCt++) {
171  testInt.assign(numPoly[cellCt], tmparray);
172  analyticInt.assign(numPoly[cellCt], tmparray);
173  *outStream << "\nIntegrals of monomials on a reference " << cellType[cellCt].getBaseCellTopologyData()->name << ":\n";
174  std::ifstream filecompare(&filename[cellCt][0]);
175  // compute integrals
176  for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
177  polyCt = 0;
178  testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
179  for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
180  for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
181  testInt[cubDeg][polyCt] = computeIntegral(cellType[cellCt], cubDeg, xDeg, yDeg);
182  polyCt++;
183  }
184  }
185  }
186  // get analytic values
187  if (filecompare.is_open()) {
188  getAnalytic(analyticInt, filecompare);
189  // close file
190  filecompare.close();
191  }
192  // perform comparison
193  for (int cubDeg=0; cubDeg <= maxDeg[cellCt]; cubDeg++) {
194  polyCt = 0;
195  offset = 0;
196  for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
197  for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
198  double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
199  double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
200  *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
201  << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << " "
202  << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
203  << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
204  if (absdiff > abstol) {
205  errorFlag++;
206  *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n";
207  }
208  polyCt++;
209  }
210  offset = offset + maxDeg[cellCt] - cubDeg;
211  }
212  *outStream << "\n";
213  }
214  *outStream << "\n";
215  } // end for cellCt
216  }
217  catch (const std::logic_error & err) {
218  *outStream << err.what() << "\n";
219  errorFlag = -1;
220  };
221 
222 
223  if (errorFlag != 0)
224  std::cout << "End Result: TEST FAILED\n";
225  else
226  std::cout << "End Result: TEST PASSED\n";
227 
228  // reset format state of std::cout
229  std::cout.copyfmt(oldFormatState);
230 
231  return errorFlag;
232 }
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...
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.