Intrepid
test_15.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 Computes integrals of monomials over a given reference cell.
61*/
62long double evalQuad(int order, int power, EIntrepidBurkardt rule) {
63
64 CubatureLineSorted<long double> lineCub(rule,order,false);
65 CubatureLineSorted<long double> lineCub2(rule,order-1,false);
66 lineCub.update(-1.0,lineCub2,1.0);
67 int size = lineCub.getNumPoints();
68 FieldContainer<long double> cubPoints(size);
69 FieldContainer<long double> cubWeights(size);
70 lineCub.getCubature(cubPoints,cubWeights);
71
72 long double Q = 0.0;
73 for (int i=0; i<size; i++) {
74 Q += cubWeights(i)*powl(cubPoints(i),(long double)power);
75 }
76 return Q;
77
78 /*
79 int mid = size/2;
80 long double Q = 0.0;
81 if (size%2)
82 Q = cubWeights(mid)*powl(cubPoints(mid),power);
83
84 for (int i=0; i<mid; i++) {
85 Q += cubWeights(i)*powl(cubPoints(i),power)+
86 cubWeights(order-i-1)*powl(cubPoints(size-i-1),power);
87 }
88 return Q;
89 */
90}
91
92int main(int argc, char *argv[]) {
93
94 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
95
96 // This little trick lets us print to std::cout only if
97 // a (dummy) command-line argument is provided.
98 int iprint = argc - 1;
99 Teuchos::RCP<std::ostream> outStream;
100 Teuchos::oblackholestream bhs; // outputs nothing
101 if (iprint > 0)
102 outStream = Teuchos::rcp(&std::cout, false);
103 else
104 outStream = Teuchos::rcp(&bhs, false);
105
106 // Save the format state of the original std::cout.
107 Teuchos::oblackholestream oldFormatState;
108 oldFormatState.copyfmt(std::cout);
109
110 *outStream \
111 << "===============================================================================\n" \
112 << "| |\n" \
113 << "| Unit Test (CubatureLineSorted) |\n" \
114 << "| |\n" \
115 << "| 1) Computing differential integrals of monomials in 1D |\n" \
116 << "| |\n" \
117 << "| Questions? Contact Drew Kouri (dpkouri@sandia.gov) or |\n" \
118 << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
119 << "| |\n" \
120 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
121 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
122 << "| |\n" \
123 << "===============================================================================\n"\
124 << "| TEST 15: differential integrals of monomials in 1D |\n"\
125 << "===============================================================================\n";
126
127 // internal variables:
128 int errorFlag = 0;
129 long double abstol = 1.0e+05*INTREPID_TOL;
130 int maxDeg = 0;
131 long double analyticInt = 0;
132 long double testInt = 0;
133 int maxOrder = 60;
134 EIntrepidBurkardt rule = BURK_LEGENDRE;
135
136 *outStream << "\nDifferential integrals of monomials on a reference line:\n";
137 // compute and compare integrals
138 try {
139 *outStream << "Testing " << EIntrepidBurkardtToString(rule) << "\n";
140 // compute integrals
141 for (int i=2; i <= maxOrder; i++) {
142 maxDeg = 2*(i-1)-1;
143 for (int j=0; j <= maxDeg; j++) {
144 testInt = evalQuad(i,j,rule);
145 long double absdiff = std::fabs(testInt);
146 *outStream << "Cubature order " << std::setw(2) << std::left << i
147 << " integrating "
148 << "x^" << std::setw(2) << std::left << j << ":" << " "
149 << std::scientific << std::setprecision(16) << testInt
150 << " " << analyticInt
151 << " " << std::setprecision(4) << absdiff << " "
152 << "<?" << " " << abstol
153 << "\n";
154 if (absdiff > abstol) {
155 errorFlag++;
156 *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
157 }
158 } // end for j
159 *outStream << "\n";
160 } // end for i
161 }
162 catch (const std::logic_error & err) {
163 *outStream << err.what() << "\n";
164 errorFlag = -1;
165 };
166
167
168 if (errorFlag != 0)
169 std::cout << "End Result: TEST FAILED\n";
170 else
171 std::cout << "End Result: TEST PASSED\n";
172
173 // reset format state of std::cout
174 std::cout.copyfmt(oldFormatState);
175
176 return errorFlag;
177}
Header file for the Intrepid::CubatureLineSorted class.
Intrepid utilities.
Utilizes cubature (integration) rules contained in the library sandia_rules (John Burkardt,...
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....