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
51#include "Teuchos_oblackholestream.hpp"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_GlobalMPISession.hpp"
54
55using namespace std;
56using namespace Intrepid;
57
58#define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
59{ \
60 ++nException; \
61 try { \
62 S ; \
63 } \
64 catch (const std::logic_error & err) { \
65 ++throwCounter; \
66 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
67 *outStream << err.what() << '\n'; \
68 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
69 }; \
70}
71
72int main(int argc, char *argv[]) {
73
74 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
75
76 // This little trick lets us print to std::cout only if
77 // a (dummy) command-line argument is provided.
78 int iprint = argc - 1;
79 Teuchos::RCP<std::ostream> outStream;
80 Teuchos::oblackholestream bhs; // outputs nothing
81 if (iprint > 0)
82 outStream = Teuchos::rcp(&std::cout, false);
83 else
84 outStream = Teuchos::rcp(&bhs, false);
85
86 // Save the format state of the original std::cout.
87 Teuchos::oblackholestream oldFormatState;
88 oldFormatState.copyfmt(std::cout);
89
90 *outStream \
91 << "===============================================================================\n" \
92 << "| |\n" \
93 << "| Unit Test (Basis_HDIV_HEX_In_FEM) |\n" \
94 << "| |\n" \
95 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
96 << "| 2) Basis values for VALUE and DIV operators |\n" \
97 << "| |\n" \
98 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
99 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
100 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
101 << "| |\n" \
102 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
103 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
104 << "| |\n" \
105 << "===============================================================================\n"\
106 << "| TEST 1: Basis creation, exception testing |\n"\
107 << "===============================================================================\n";
108
109 // Define basis and error flag
110 const int deg = 1;
111 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
113 FieldContainer<double> openPts(PointTools::getLatticeSize(line,deg+1,1),1);
114 PointTools::getLattice<double,FieldContainer<double> >(closedPts,line,deg);
115 PointTools::getLattice<double,FieldContainer<double> >(openPts,line,deg+1,1);
116
117 Basis_HDIV_HEX_In_FEM<double, FieldContainer<double> > hexBasis(deg,closedPts,openPts);
118
119 int errorFlag = 0;
120
121 // Initialize throw counter for exception testing
122 int nException = 0;
123 int throwCounter = 0;
124
125 // compute values at vertices: there are 8 of them
126 FieldContainer<double> hexNodes(8, 3);
127 hexNodes(0,0) = -1.0; hexNodes(0,1) = -1.0; hexNodes(0,2) = -1.0;
128 hexNodes(1,0) = 1.0; hexNodes(1,1) = -1.0; hexNodes(1,2) = -1.0;
129 hexNodes(2,0) = -1.0; hexNodes(2,1) = 1.0; hexNodes(2,2) = -1.0;
130 hexNodes(3,0) = 1.0; hexNodes(3,1) = 1.0; hexNodes(3,2) = -1.0;
131 hexNodes(4,0) = -1.0; hexNodes(4,1) = -1.0; hexNodes(4,2) = 1.0;
132 hexNodes(5,0) = 1.0; hexNodes(5,1) = -1.0; hexNodes(5,2) = 1.0;
133 hexNodes(6,0) = -1.0; hexNodes(6,1) = 1.0; hexNodes(6,2) = 1.0;
134 hexNodes(7,0) = 1.0; hexNodes(7,1) = 1.0; hexNodes(7,2) = 1.0;
135
136
137
138 // Generic array for the output values; needs to be properly resized depending on the operator type
140
141 try{
142 // exception #1: GRAD cannot be applied to HDIV functions
143 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
144 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 3 );
145 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD), throwCounter, nException );
146
147 // exception #2: CURL cannot be applied to HDIV functions
148 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_CURL), throwCounter, nException );
149
150 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
151 // getDofTag() to access invalid array elements thereby causing bounds check exception
152 // exception #3
153 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,0,0), throwCounter, nException );
154 // exception #4
155 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,1,1), throwCounter, nException );
156 // exception #5
157 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
158 // exception #6
159 INTREPID_TEST_COMMAND( hexBasis.getDofTag(12), throwCounter, nException );
160 // exception #7
161 INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException );
162
163#ifdef HAVE_INTREPID_DEBUG
164 // Exceptions 8- test exception handling with incorrectly dimensioned input/output arrays
165 // exception #8: input points array must be of rank-2
166 FieldContainer<double> badPoints1(4, 5, 3);
167 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
168
169 // exception #9 dimension 1 in the input point array must equal space dimension of the cell
170 FieldContainer<double> badPoints2(4, 2);
171 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
172
173 // exception #10 output values must be of rank-3 for OPERATOR_VALUE
174 FieldContainer<double> badVals1(4, 3);
175 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
176
177 // exception #11 output values must be of rank-2 for OPERATOR_DIV
178 FieldContainer<double> badVals2(4, 3, 3);
179 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_DIV), throwCounter, nException );
180
181 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
182 FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0), 3);
183 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException );
184
185 // exception #13 incorrect 0th dimension of output array (must equal number of basis functions)
186 FieldContainer<double> badVals4(hexBasis.getCardinality() + 1, hexNodes.dimension(0));
187 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_DIV), throwCounter, nException );
188
189 // exception #14 incorrect 1st dimension of output array (must equal number of points)
190 FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0) + 1, 3);
191 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_VALUE), throwCounter, nException );
192
193 // exception #15 incorrect 1st dimension of output array (must equal number of points)
194 FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0) + 1);
195 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_DIV), throwCounter, nException );
196
197 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
198 FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 4);
199 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, OPERATOR_VALUE), throwCounter, nException );
200#endif
201
202 }
203 catch (const std::logic_error & err) {
204 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
205 *outStream << err.what() << '\n';
206 *outStream << "-------------------------------------------------------------------------------" << "\n\n";
207 errorFlag = -1000;
208 };
209
210 // Check if number of thrown exceptions matches the one we expect
211 if (throwCounter != nException) {
212 errorFlag++;
213 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
214 }
215
216 *outStream \
217 << "\n"
218 << "===============================================================================\n"\
219 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
220 << "===============================================================================\n";
221
222 try{
223 std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags();
224
225 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
226 for (unsigned i = 0; i < allTags.size(); i++) {
227 int bfOrd = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
228
229 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
230 if( !( (myTag[0] == allTags[i][0]) &&
231 (myTag[1] == allTags[i][1]) &&
232 (myTag[2] == allTags[i][2]) &&
233 (myTag[3] == allTags[i][3]) ) ) {
234 errorFlag++;
235 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
236 *outStream << " getDofOrdinal( {"
237 << allTags[i][0] << ", "
238 << allTags[i][1] << ", "
239 << allTags[i][2] << ", "
240 << allTags[i][3] << "}) = " << bfOrd <<" but \n";
241 *outStream << " getDofTag(" << bfOrd << ") = { "
242 << myTag[0] << ", "
243 << myTag[1] << ", "
244 << myTag[2] << ", "
245 << myTag[3] << "}\n";
246 }
247 }
248
249 // Now do the same but loop over basis functions
250 for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) {
251 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
252 int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
253 if( bfOrd != myBfOrd) {
254 errorFlag++;
255 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
256 *outStream << " getDofTag(" << bfOrd << ") = { "
257 << myTag[0] << ", "
258 << myTag[1] << ", "
259 << myTag[2] << ", "
260 << myTag[3] << "} but getDofOrdinal({"
261 << myTag[0] << ", "
262 << myTag[1] << ", "
263 << myTag[2] << ", "
264 << myTag[3] << "} ) = " << myBfOrd << "\n";
265 }
266 }
267 }
268 catch (const std::logic_error & err){
269 *outStream << err.what() << "\n\n";
270 errorFlag = -1000;
271 };
272
273 *outStream \
274 << "\n"
275 << "===============================================================================\n"\
276 << "| TEST 3: correctness of basis function values |\n"\
277 << "===============================================================================\n";
278
279 outStream -> precision(20);
280
281 // VALUE: Each row pair gives the 6x3 correct basis set values at an evaluation point: (P,F,D) layout
282 double basisValues[] = {
283 // basis function 0 (in from x==-1 plane, y and z are constant functions)
284 1, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0,
285 1, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0,
286 // basis function 1 (out from x==1 plane, y and z are constant functions
287 0, 0, 0, 1, 0, 0, 0, 0, 0, 1, 0, 0,
288 0, 0, 0, 1, 0, 0, 0, 0, 0, 1, 0, 0,
289 // basis function 2 (in from y==-1 plane, x and z are constant functions
290 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0,
291 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0,
292 // basis function 3 (out from y == 1 plane, x and z are constant function
293 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0,
294 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0,
295 // basis function 4 (in from z == -1 plane, x and y are constant function
296 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 1,
297 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
298 // basis function 4 (out from z == 1 plane, x and y are constant function
299 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
300 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 1
301 };
302
303 // DIV: each row gives the 6 correct values of the divergence of the 6 basis functions: (P,F) layout
304 double basisDivs[] = {
305 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
306 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
307 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
308 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
309 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
310 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5
311 };
312
313 try{
314
315 // Dimensions for the output arrays:
316 int numPoints = hexNodes.dimension(0);
317 int numFields = hexBasis.getCardinality();
318 int spaceDim = hexBasis.getBaseCellTopology().getDimension();
319
320 // Generic array for values and curls that will be properly sized before each call
322
323 // Check VALUE of basis functions: resize vals to rank-3 container:
324 vals.resize(numFields, numPoints, spaceDim);
325 hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
326 for (int i = 0; i < numFields; i++) {
327 for (int j = 0; j < numPoints; j++) {
328 for (int k = 0; k < spaceDim; k++) {
329 int l = k + i * numPoints * spaceDim + j * spaceDim;
330 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
331 errorFlag++;
332 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
333
334 // Output the multi-index of the value where the error is:
335 *outStream << " At multi-index { ";
336 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
337 *outStream << "} computed value: " << vals(i,j,k)
338 << " but reference value: " << basisValues[l] << "\n";
339 }
340 }
341 }
342 }
343
344 // Check DIV of basis function: resize vals to rank-2 container
345 vals.resize(numFields, numPoints);
346 hexBasis.getValues(vals, hexNodes, OPERATOR_DIV);
347 for (int i = 0; i < numFields; i++) {
348 for (int j = 0; j < numPoints; j++) {
349 int l = i * numPoints + j;
350 if (std::abs(vals(i,j) - basisDivs[l]) > INTREPID_TOL) {
351 errorFlag++;
352 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
353
354 // Output the multi-index of the value where the error is:
355 *outStream << " At multi-index { ";
356 *outStream << i << " ";*outStream << j << " ";
357 *outStream << "} computed divergence component: " << vals(i,j)
358 << " but reference divergence component: " << basisDivs[l] << "\n";
359 }
360 }
361 }
362
363 }
364
365 // Catch unexpected errors
366 catch (const std::logic_error & err) {
367 *outStream << err.what() << "\n\n";
368 errorFlag = -1000;
369 };
370
371 if (errorFlag != 0)
372 std::cout << "End Result: TEST FAILED\n";
373 else
374 std::cout << "End Result: TEST PASSED\n";
375
376 // reset format state of std::cout
377 std::cout.copyfmt(oldFormatState);
378
379 return errorFlag;
380}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HDIV_HEX_In_FEM class.
Header file for utility class to provide point tools, such as barycentric coordinates,...
Implementation of the default H(div)-compatible FEM basis of degree 1 on Hexahedral cell.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0.
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...