Intrepid
test_02.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
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*/
66double 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*/
80double computeIntegral(int cubDegree, int polyDegree) {
81
82 DefaultCubatureFactory<double> cubFactory; // create factory
83 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >()); // create cell topology
84 Teuchos::RCP<Cubature<double> > lineCub = cubFactory.create(line, cubDegree); // create default cubature
85 double val = 0.0;
86
87 int cubDim = lineCub->getDimension();
88
89 int numCubPoints = lineCub->getNumPoints();
90
91 FieldContainer<double> point(cubDim);
92 FieldContainer<double> cubPoints(numCubPoints, cubDim);
93 FieldContainer<double> cubWeights(numCubPoints);
94
95 lineCub->getCubature(cubPoints, cubWeights);
96
97 for (int i=0; i<numCubPoints; i++) {
98 for (int j=0; j<cubDim; j++) {
99 point(j) = cubPoints(i,j);
100 }
101 val += computeMonomial(point, polyDegree)*cubWeights(i);
102 }
103
104 return val;
105}
106
107
108int main(int argc, char *argv[]) {
109
110 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
111
112 // This little trick lets us print to std::cout only if
113 // a (dummy) command-line argument is provided.
114 int iprint = argc - 1;
115 Teuchos::RCP<std::ostream> outStream;
116 Teuchos::oblackholestream bhs; // outputs nothing
117 if (iprint > 0)
118 outStream = Teuchos::rcp(&std::cout, false);
119 else
120 outStream = Teuchos::rcp(&bhs, false);
121
122 // Save the format state of the original std::cout.
123 Teuchos::oblackholestream oldFormatState;
124 oldFormatState.copyfmt(std::cout);
125
126 *outStream \
127 << "===============================================================================\n" \
128 << "| |\n" \
129 << "| Unit Test (CubatureDirectLineGauss) |\n" \
130 << "| |\n" \
131 << "| 1) Computing integrals of monomials on reference cells in 1D |\n" \
132 << "| |\n" \
133 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
134 << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
135 << "| |\n" \
136 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
137 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
138 << "| |\n" \
139 << "===============================================================================\n"\
140 << "| TEST 1: integrals of monomials in 1D |\n"\
141 << "===============================================================================\n";
142
143 // internal variables:
144 int errorFlag = 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+01 * INTREPID_TOL;
149 testInt.assign(INTREPID_CUBATURE_LINE_GAUSS_MAX+1, tmparray);
150 analyticInt.assign(INTREPID_CUBATURE_LINE_GAUSS_MAX+1, tmparray);
151
152 // open file with analytic values
153 std::string basedir = "./data";
154 std::stringstream namestream;
155 std::string filename;
156 namestream << basedir << "/EDGE_integrals" << ".dat";
157 namestream >> filename;
158 std::ifstream filecompare(&filename[0]);
159
160 *outStream << "\nIntegrals of monomials on a reference line (edge):\n";
161
162 // compute and compare integrals
163 try {
164 // compute integrals
165 for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_GAUSS_MAX; cubDeg++) {
166 testInt[cubDeg].resize(cubDeg+1);
167 for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
168 testInt[cubDeg][polyDeg] = computeIntegral(cubDeg, polyDeg);
169 }
170 }
171 // get analytic values
172 if (filecompare.is_open()) {
173 getAnalytic(analyticInt, filecompare);
174 // close file
175 filecompare.close();
176 }
177 // perform comparison
178 for (int cubDeg=0; cubDeg <= INTREPID_CUBATURE_LINE_GAUSS_MAX; cubDeg++) {
179 for (int polyDeg=0; polyDeg <= cubDeg; polyDeg++) {
180 double abstol = ( analyticInt[polyDeg][0] == 0.0 ? reltol : std::fabs(reltol*analyticInt[polyDeg][0]) );
181 double absdiff = std::fabs(analyticInt[polyDeg][0] - testInt[cubDeg][polyDeg]);
182 *outStream << "Cubature order " << std::setw(2) << std::left << cubDeg << " integrating "
183 << "x^" << std::setw(2) << std::left << polyDeg << ":" << " "
184 << std::scientific << std::setprecision(16) << testInt[cubDeg][polyDeg] << " " << analyticInt[polyDeg][0] << " "
185 << std::setprecision(4) << absdiff << " " << "<?" << " " << abstol << "\n";
186 if (absdiff > abstol) {
187 errorFlag++;
188 *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
189 }
190 }
191 *outStream << "\n";
192 } // end for cubDeg
193 }
194 catch (const std::logic_error & err) {
195 *outStream << err.what() << "\n";
196 errorFlag = -1;
197 };
198
199
200 if (errorFlag != 0)
201 std::cout << "End Result: TEST FAILED\n";
202 else
203 std::cout << "End Result: TEST PASSED\n";
204
205 // reset format state of std::cout
206 std::cout.copyfmt(oldFormatState);
207
208 return errorFlag;
209}
#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 abstract base class Intrepid::DefaultCubatureFactory.
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...