54 #include "Teuchos_oblackholestream.hpp"
55 #include "Teuchos_RCP.hpp"
56 #include "Teuchos_GlobalMPISession.hpp"
58 using namespace Intrepid;
63 long double evalQuad(
int order,
int power, EIntrepidBurkardt rule) {
66 int size = lineCub.getNumPoints();
69 lineCub.getCubature(cubPoints,cubWeights);
74 Q = cubWeights(mid)*powl(cubPoints(mid),power);
76 for (
int i=0; i<mid; i++) {
77 Q += cubWeights(i)*powl(cubPoints(i),power)+
78 cubWeights(size-i-1)*powl(cubPoints(size-i-1),power);
83 long double evalInt(
int power, EIntrepidBurkardt rule) {
86 if (rule==BURK_CLENSHAWCURTIS||rule==BURK_FEJER2||
87 rule==BURK_LEGENDRE||rule==BURK_PATTERSON ||
88 rule==BURK_TRAPEZOIDAL) {
92 I = 2.0/((
long double)power+1.0);
94 else if (rule==BURK_LAGUERRE) {
95 I = tgammal((
long double)(power+1));
97 else if (rule==BURK_CHEBYSHEV1) {
101 for (
int i=2;i<=power;i+=2) {
102 top *= (
long double)(i-1);
103 bot *= (
long double)i;
111 else if (rule==BURK_CHEBYSHEV2) {
112 long double bot, top;
115 for (
int i=2;i<=power;i+=2) {
116 top *= (
long double)(i-1);
117 bot *= (
long double)i;
119 bot *= (
long double)(power+2);
126 else if (rule==BURK_HERMITE||rule==BURK_GENZKEISTER) {
131 long double value = 1.0;
133 int n_copy = power-1;
135 value *= (
long double)n_copy;
139 I = value*sqrt(M_PI)/powl(2.0,(
long double)power/2.0);
145 int main(
int argc,
char *argv[]) {
147 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
151 int iprint = argc - 1;
152 Teuchos::RCP<std::ostream> outStream;
153 Teuchos::oblackholestream bhs;
155 outStream = Teuchos::rcp(&std::cout,
false);
157 outStream = Teuchos::rcp(&bhs,
false);
160 Teuchos::oblackholestream oldFormatState;
161 oldFormatState.copyfmt(std::cout);
164 <<
"===============================================================================\n" \
166 <<
"| Unit Test (CubatureLineSorted) |\n" \
168 <<
"| 1) Computing integrals of monomials in 1D |\n" \
170 <<
"| Questions? Contact Drew Kouri (dpkouri@sandia.gov) or |\n" \
171 <<
"| Denis Ridzal (dridzal@sandia.gov). |\n" \
173 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
174 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
176 <<
"===============================================================================\n"\
177 <<
"| TEST 11: integrals of monomials in 1D |\n"\
178 <<
"===============================================================================\n";
182 long double reltol = 1.0e+05*INTREPID_TOL;
184 long double analyticInt = 0;
185 long double testInt = 0;
188 *outStream <<
"\nIntegrals of monomials on a reference line (edge):\n";
191 for (EIntrepidBurkardt rule=BURK_CHEBYSHEV1; rule <= BURK_LAGUERRE; rule++) {
192 *outStream <<
"Testing " << EIntrepidBurkardtToString(rule) <<
"\n";
194 if (rule==BURK_HERMITE)
196 else if (rule==BURK_TRAPEZOIDAL)
201 if (rule!=BURK_PATTERSON&&rule!=BURK_GENZKEISTER) {
202 for (
int i=1; i <= maxOrder; i++) {
203 if ( rule==BURK_CHEBYSHEV1 ||
204 rule==BURK_CHEBYSHEV2 ||
205 rule==BURK_LEGENDRE ||
206 rule==BURK_LAGUERRE ||
209 else if ( rule==BURK_CLENSHAWCURTIS ||
211 rule==BURK_TRAPEZOIDAL )
214 for (
int j=0; j <= maxDeg; j++) {
215 analyticInt = evalInt(j,rule);
216 testInt = evalQuad(i,j,rule);
218 long double abstol = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
219 long double absdiff = std::fabs(analyticInt - testInt);
220 *outStream <<
"Cubature order " << std::setw(2) << std::left
221 << i <<
" integrating "
222 <<
"x^" << std::setw(2) << std::left << j <<
":"
224 << std::scientific << std::setprecision(16) << testInt
225 <<
" " << analyticInt
226 <<
" " << std::setprecision(4) << absdiff <<
" "
227 <<
"<?" <<
" " << abstol
229 if (absdiff > abstol) {
231 *outStream << std::right << std::setw(104)
232 <<
"^^^^---FAILURE!\n";
238 else if (rule==BURK_PATTERSON) {
239 for (
int i=0; i < 8; i++) {
240 int l = (int)std::pow(2.0,(
double)i+1.0)-1;
244 maxDeg = (int)(1.5*(
double)l+0.5);
245 for (
int j=0; j <= maxDeg; j++) {
246 analyticInt = evalInt(j,rule);
247 testInt = evalQuad(l,j,rule);
249 long double abstol = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
250 long double absdiff = std::fabs(analyticInt - testInt);
251 *outStream <<
"Cubature order " << std::setw(2) << std::left
252 << l <<
" integrating "
253 <<
"x^" << std::setw(2) << std::left << j <<
":"
255 << std::scientific << std::setprecision(16) << testInt
256 <<
" " << analyticInt
257 <<
" " << std::setprecision(4) << absdiff <<
" "
258 <<
"<?" <<
" " << abstol
260 if (absdiff > abstol) {
262 *outStream << std::right << std::setw(104)
263 <<
"^^^^---FAILURE!\n";
269 else if (rule==BURK_GENZKEISTER) {
271 int o_ghk[4] = {1,3, 9,19};
272 int p_ghk[4] = {1,5,15,29};
273 for (
int i=0; i < 4; i++) {
276 for (
int j=0; j <= maxDeg; j++) {
277 analyticInt = evalInt(j,rule);
278 testInt = evalQuad(l,j,rule);
280 long double abstol = (analyticInt == 0.0 ? reltol : std::fabs(reltol*analyticInt) );
281 long double absdiff = std::fabs(analyticInt - testInt);
282 *outStream <<
"Cubature order " << std::setw(2) << std::left
283 << l <<
" integrating "
284 <<
"x^" << std::setw(2) << std::left << j <<
":"
286 << std::scientific << std::setprecision(16) << testInt
287 <<
" " << analyticInt
288 <<
" " << std::setprecision(4) << absdiff <<
" "
289 <<
"<?" <<
" " << abstol
291 if (absdiff > abstol) {
293 *outStream << std::right << std::setw(104) <<
"^^^^---FAILURE!\n";
301 catch (std::logic_error err) {
302 *outStream << err.what() <<
"\n";
308 std::cout <<
"End Result: TEST FAILED\n";
310 std::cout <<
"End Result: TEST PASSED\n";
313 std::cout.copyfmt(oldFormatState);
Utilizes cubature (integration) rules contained in the library sandia_rules (John Burkardt...
Header file for the Intrepid::CubatureLineSorted class.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity...