Intrepid
test_07.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
57using namespace Intrepid;
58
59//Teuchos::RCP<std::ostream> outStream;
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*/
67double 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*/
81double computeIntegral(int cubDegree, int xDeg, int yDeg) {
82
83 CubatureSparse<double,2> myCub(cubDegree);
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
106int 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 (CubatureSparse) |\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), |\n" \
132 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
133 << "| Matthew Keegan (mskeega@sandia.gov) |\n" \
134 << "| |\n" \
135 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
136 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
137 << "| |\n" \
138 << "===============================================================================\n"\
139 << "| TEST 1: integrals of monomials in 2D for Sparse Grid Construction |\n"\
140 << "===============================================================================\n";
141
142 // internal variables:
143 int errorFlag = 0;
144 int polyCt = 0;
145 int offset = 0;
146 Teuchos::Array< Teuchos::Array<double> > testInt;
147 Teuchos::Array< Teuchos::Array<double> > analyticInt;
148 Teuchos::Array<double> tmparray(1);
149 double reltol = 1.0e+03 * INTREPID_TOL;
150 int maxDeg = 30; // can be as large as INTREPID_CUBATURE_SPARSE2D_GAUSS_MAX, but runtime is excessive
151 int maxOffset = INTREPID_CUBATURE_LINE_GAUSS_MAX;
152 int numPoly = (maxDeg+1)*(maxDeg+2)/2;
153 int numAnalytic = (maxOffset+1)*(maxOffset+2)/2;
154 testInt.assign(numPoly, tmparray);
155 analyticInt.assign(numAnalytic, tmparray);
156
157 // get names of files with analytic values
158 std::string basedir = "./data";
159 std::stringstream namestream;
160 std::string filename;
161 namestream << basedir << "/QUAD_integrals" << ".dat";
162 namestream >> filename;
163
164 // compute and compare integrals
165 try {
166 *outStream << "\nIntegrals of monomials:\n";
167
168 std::ifstream filecompare(&filename[0]);
169 // compute integrals
170 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
171 polyCt = 0;
172 testInt[cubDeg].resize((cubDeg+1)*(cubDeg+2)/2);
173 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
174 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
175 testInt[cubDeg][polyCt] = computeIntegral(cubDeg, xDeg, yDeg);
176 polyCt++;
177 }
178 }
179 }
180
181 // get analytic values
182 if (filecompare.is_open()) {
183 getAnalytic(analyticInt, filecompare);
184 // close file
185 filecompare.close();
186 }
187 // perform comparison
188 for (int cubDeg=0; cubDeg <= maxDeg; cubDeg++) {
189 polyCt = 0;
190 offset = 0;
191 for (int xDeg=0; xDeg <= cubDeg; xDeg++) {
192 for (int yDeg=0; yDeg <= cubDeg-xDeg; yDeg++) {
193 double abstol = ( analyticInt[polyCt+offset][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyCt+offset][0]) );
194 double absdiff = std::fabs(analyticInt[polyCt+offset][0] - testInt[cubDeg][polyCt]);
195 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
196 << "x^" << std::setw(2) << std::left << xDeg << " * y^" << std::setw(2) << yDeg << ":" << " "
197 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyCt] << " " << analyticInt[polyCt+offset][0] << " "
198 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
199 if (absdiff > abstol) {
200 errorFlag++;
201 *outStream << std::right << std::setw(111) << "^^^^---FAILURE!\n";
202 }
203 polyCt++;
204 }
205 offset = offset + maxOffset - cubDeg;
206 }
207 *outStream << "\n";
208 }
209 *outStream << "\n";
210 }
211 catch (const std::logic_error & err) {
212 *outStream << err.what() << "\n";
213 errorFlag = -1;
214 };
215
216
217 if (errorFlag != 0)
218 std::cout << "End Result: TEST FAILED\n";
219 else
220 std::cout << "End Result: TEST PASSED\n";
221
222 // reset format state of std::cout
223 std::cout.copyfmt(oldFormatState);
224
225 return errorFlag;
226}
#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...
Header file for the Intrepid::CubatureSparse class.
Intrepid utilities.
void getAnalytic(Teuchos::Array< Teuchos::Array< Scalar > > &testMat, std::ifstream &inputFile, TypeOfExactData analyticDataType=INTREPID_UTILS_FRACTION)
Loads analytic values stored in a file into the matrix testMat, where the output matrix is an array o...