Intrepid
test_01.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
51#include "Teuchos_oblackholestream.hpp"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_GlobalMPISession.hpp"
56#include "Shards_CellTopology.hpp"
57
58#include <iostream>
59using namespace Intrepid;
60
65int main(int argc, char *argv[]) {
66
67 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
68
69 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
70 int iprint = argc - 1;
71
72 Teuchos::RCP<std::ostream> outStream;
73 Teuchos::oblackholestream bhs; // outputs nothing
74
75 if (iprint > 0)
76 outStream = Teuchos::rcp(&std::cout, false);
77 else
78 outStream = Teuchos::rcp(&bhs, false);
79
80 // Save the format state of the original std::cout.
81 Teuchos::oblackholestream oldFormatState;
82 oldFormatState.copyfmt(std::cout);
83
84 *outStream \
85 << "===============================================================================\n" \
86 << "| |\n" \
87 << "| Unit Test HDIV_TRI_In_FEM |\n" \
88 << "| |\n" \
89 << "| 1) Tests triangular Raviart-Thomas basis |\n" \
90 << "| |\n" \
91 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
92 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
93 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
94 << "| |\n" \
95 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
96 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
97 << "| |\n" \
98 << "===============================================================================\n";
99
100 int errorFlag = 0;
101
102 // test for basis values, compare against fiat
103 try {
104 const int deg = 2;
105 Basis_HDIV_TRI_In_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_EQUISPACED );
106
107 // Get a lattice
108 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 );
109 FieldContainer<double> lattice( np_lattice , 2 );
110 FieldContainer<double> myBasisValues( myBasis.getCardinality() , np_lattice , 2 );
111 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
112 myBasis.getBaseCellTopology() ,
113 deg ,
114 0 ,
115 POINTTYPE_EQUISPACED );
116
117 myBasis.getValues( myBasisValues , lattice , OPERATOR_VALUE );
118
119 const double fiat_vals[] = {
120 0.000000000000000e+00, -2.000000000000000e+00,
121 2.500000000000000e-01, -5.000000000000000e-01,
122 -1.000000000000000e+00, 1.000000000000000e+00,
123 0.000000000000000e+00, -2.500000000000000e-01,
124 -5.000000000000000e-01, 5.000000000000000e-01,
125 0.000000000000000e+00, 0.000000000000000e+00,
126 0.000000000000000e+00, 1.000000000000000e+00,
127 2.500000000000000e-01, -5.000000000000000e-01,
128 2.000000000000000e+00, -2.000000000000000e+00,
129 0.000000000000000e+00, 5.000000000000000e-01,
130 2.500000000000000e-01, -2.500000000000000e-01,
131 0.000000000000000e+00, 0.000000000000000e+00,
132 0.000000000000000e+00, 0.000000000000000e+00,
133 2.500000000000000e-01, 0.000000000000000e+00,
134 2.000000000000000e+00, 0.000000000000000e+00,
135 0.000000000000000e+00, -5.000000000000000e-01,
136 2.500000000000000e-01, 2.500000000000000e-01,
137 0.000000000000000e+00, -1.000000000000000e+00,
138 0.000000000000000e+00, 0.000000000000000e+00,
139 -5.000000000000000e-01, 0.000000000000000e+00,
140 -1.000000000000000e+00, 0.000000000000000e+00,
141 0.000000000000000e+00, 2.500000000000000e-01,
142 2.500000000000000e-01, 2.500000000000000e-01,
143 0.000000000000000e+00, 2.000000000000000e+00,
144 1.000000000000000e+00, 0.000000000000000e+00,
145 5.000000000000000e-01, 0.000000000000000e+00,
146 0.000000000000000e+00, 0.000000000000000e+00,
147 -5.000000000000000e-01, 2.500000000000000e-01,
148 -2.500000000000000e-01, 2.500000000000000e-01,
149 -2.000000000000000e+00, 2.000000000000000e+00,
150 -2.000000000000000e+00, 0.000000000000000e+00,
151 -2.500000000000000e-01, 0.000000000000000e+00,
152 0.000000000000000e+00, 0.000000000000000e+00,
153 -5.000000000000000e-01, 2.500000000000000e-01,
154 5.000000000000000e-01, -5.000000000000000e-01,
155 1.000000000000000e+00, -1.000000000000000e+00,
156 0.000000000000000e+00, 0.000000000000000e+00,
157 1.500000000000000e+00, 0.000000000000000e+00,
158 0.000000000000000e+00, 0.000000000000000e+00,
159 0.000000000000000e+00, 7.500000000000000e-01,
160 7.500000000000000e-01, -7.500000000000000e-01,
161 0.000000000000000e+00, 0.000000000000000e+00,
162 0.000000000000000e+00, 0.000000000000000e+00,
163 7.500000000000000e-01, 0.000000000000000e+00,
164 0.000000000000000e+00, 0.000000000000000e+00,
165 0.000000000000000e+00, 1.500000000000000e+00,
166 -7.500000000000000e-01, 7.500000000000000e-01,
167 0.000000000000000e+00, 0.000000000000000e+00
168 };
169
170 int cur=0;
171 for (int i=0;i<myBasisValues.dimension(0);i++) {
172 for (int j=0;j<myBasisValues.dimension(1);j++) {
173 for (int k=0;k<myBasisValues.dimension(2);k++) {
174 if (std::abs( myBasisValues(i,j,k) - fiat_vals[cur] ) > INTREPID_TOL ) {
175 errorFlag++;
176 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
177
178 // Output the multi-index of the value where the error is:
179 *outStream << " At multi-index { ";
180 *outStream << i << " " << j << " " << k;
181 *outStream << "} computed value: " << myBasisValues(i,j,k)
182 << " but correct value: " << fiat_vals[cur] << "\n";
183 *outStream << "Difference: " << std::abs( myBasisValues(i,j,k) - fiat_vals[cur] ) << "\n";
184 }
185 cur++;
186 }
187 }
188 }
189 }
190 catch (const std::exception & err) {
191 *outStream << err.what() << "\n\n";
192 errorFlag = -1000;
193 }
194 try {
195 const int deg = 2;
196 Basis_HDIV_TRI_In_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_EQUISPACED );
197
198 // Get a lattice
199 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 );
200 FieldContainer<double> lattice( np_lattice , 2 );
201 FieldContainer<double> myBasisDivs( myBasis.getCardinality() , np_lattice );
202 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
203 myBasis.getBaseCellTopology() ,
204 deg ,
205 0 ,
206 POINTTYPE_EQUISPACED );
207
208 myBasis.getValues( myBasisDivs , lattice , OPERATOR_DIV );
209
210
211 const double fiat_divs[] = {
212 7.000000000000000e+00,
213 2.500000000000000e+00,
214 -2.000000000000000e+00,
215 2.500000000000000e+00,
216 -2.000000000000000e+00,
217 -2.000000000000000e+00,
218 -2.000000000000000e+00,
219 2.500000000000000e+00,
220 7.000000000000000e+00,
221 -2.000000000000000e+00,
222 2.500000000000000e+00,
223 -2.000000000000000e+00,
224 -2.000000000000000e+00,
225 2.500000000000000e+00,
226 7.000000000000000e+00,
227 -2.000000000000000e+00,
228 2.500000000000000e+00,
229 -2.000000000000000e+00,
230 -2.000000000000000e+00,
231 -2.000000000000000e+00,
232 -2.000000000000000e+00,
233 2.500000000000000e+00,
234 2.500000000000000e+00,
235 7.000000000000000e+00,
236 -2.000000000000000e+00,
237 -2.000000000000000e+00,
238 -2.000000000000000e+00,
239 2.500000000000000e+00,
240 2.500000000000000e+00,
241 7.000000000000000e+00,
242 7.000000000000000e+00,
243 2.500000000000000e+00,
244 -2.000000000000000e+00,
245 2.500000000000000e+00,
246 -2.000000000000000e+00,
247 -2.000000000000000e+00,
248 9.000000000000000e+00,
249 0.000000000000000e+00,
250 -9.000000000000000e+00,
251 4.500000000000000e+00,
252 -4.500000000000000e+00,
253 0.000000000000000e+00,
254 9.000000000000000e+00,
255 4.500000000000000e+00,
256 0.000000000000000e+00,
257 0.000000000000000e+00,
258 -4.500000000000000e+00,
259 -9.000000000000000e+00
260 };
261
262 int cur=0;
263 for (int i=0;i<myBasisDivs.dimension(0);i++) {
264 for (int j=0;j<myBasisDivs.dimension(1);j++) {
265 if (std::abs( myBasisDivs(i,j) - fiat_divs[cur] ) > INTREPID_TOL ) {
266 errorFlag++;
267 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
268
269 // Output the multi-index of the value where the error is:
270 *outStream << " At multi-index { ";
271 *outStream << i << " " << j;
272 *outStream << "} computed value: " << myBasisDivs(i,j)
273 << " but correct value: " << fiat_divs[cur] << "\n";
274 *outStream << "Difference: " << std::abs( myBasisDivs(i,j) - fiat_divs[cur] ) << "\n";
275 }
276 cur++;
277 }
278 }
279 }
280 catch (const std::exception & err) {
281 *outStream << err.what() << "\n\n";
282 errorFlag = -1000;
283 }
284
285
286 if (errorFlag != 0)
287 std::cout << "End Result: TEST FAILED\n";
288 else
289 std::cout << "End Result: TEST PASSED\n";
290
291 // reset format state of std::cout
292 std::cout.copyfmt(oldFormatState);
293
294 return errorFlag;
295}
int main(int argc, char *argv[])
Performs a code-code comparison to FIAT for Raviart-Thomas bases on triangles (values and divs)
Definition: test_01.cpp:65
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HDIV_TRI_In_FEM class.
Header file for utility class to provide point tools, such as barycentric coordinates,...
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual int getCardinality() const
Returns cardinality of the basis.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
int dimension(const int whichDim) const
Returns the specified dimension.
static int getLatticeSize(const shards::CellTopology &cellType, const int order, const int offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...