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_QUAD_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 // get points first
111 const int deg = 1;
112 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
114 FieldContainer<double> openPts(PointTools::getLatticeSize(line,deg+1,1),1);
115 PointTools::getLattice<double,FieldContainer<double> >(closedPts,line,deg);
116 PointTools::getLattice<double,FieldContainer<double> >(openPts,line,deg+1,1);
117
118 Basis_HDIV_QUAD_In_FEM<double, FieldContainer<double> > quadBasis(deg,closedPts,openPts);
119 int errorFlag = 0;
120
121 // Initialize throw counter for exception testing
122 int nException = 0;
123 int throwCounter = 0;
124
125 // Array with a lattice of 9 points on the quad
126 FieldContainer<double> quadNodes(9, 2);
127 quadNodes(0,0) = -1.0; quadNodes(0,1) = -1.0;
128 quadNodes(1,0) = 0.0; quadNodes(1,1) = -1.0;
129 quadNodes(2,0) = 1.0; quadNodes(2,1) = -1.0;
130 quadNodes(3,0) = -1.0; quadNodes(3,1) = 0.0;
131 quadNodes(4,0) = 0.0; quadNodes(4,1) = 0.0;
132 quadNodes(5,0) = 1.0; quadNodes(5,1) = 0.0;
133 quadNodes(6,0) = -1.0; quadNodes(6,1) = 1.0;
134 quadNodes(7,0) = 0.0; quadNodes(7,1) = 1.0;
135 quadNodes(8,0) = 1.0; quadNodes(8,1) = 1.0;
136
137 // Generic array for the output values; needs to be properly resized depending on the operator type
139
140 try{
141
142 int spaceDim = quadBasis.getBaseCellTopology().getDimension();
143
144 // exception #1: GRAD cannot be applied to HDIV functions
145 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
146 vals.resize(quadBasis.getCardinality(), quadNodes.dimension(0), spaceDim );
147 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD), throwCounter, nException );
148
149 // exception #2: CURL cannot be applied to HDIV functions
150 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_CURL), throwCounter, nException );
151
152 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
153 // getDofTag() to access invalid array elements thereby causing bounds check exception
154 // exception #3
155 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException );
156 // exception #4
157 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException );
158 // exception #5
159 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,1), throwCounter, nException );
160 // exception #6
161 INTREPID_TEST_COMMAND( quadBasis.getDofTag(12), throwCounter, nException );
162 // exception #7
163 INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException );
164
165#ifdef HAVE_INTREPID_DEBUG
166 // Exceptions 8- test exception handling with incorrectly dimensioned input/output arrays
167 // exception #8: input points array must be of rank-2
168 FieldContainer<double> badPoints1(4, 5, 3);
169 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
170
171 // exception #9 dimension 1 in the input point array must equal space dimension of the cell
172 FieldContainer<double> badPoints2(4, spaceDim + 1);
173 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
174
175 // exception #10 output values must be of rank-3 for OPERATOR_VALUE
176 FieldContainer<double> badVals1(4, 5);
177 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException );
178
179 // exception #11 output values must be of rank-2 for OPERATOR_DIV
180 FieldContainer<double> badVals2(4, 5, spaceDim);
181 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_DIV), throwCounter, nException );
182
183 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
184 FieldContainer<double> badVals3(quadBasis.getCardinality() + 1, quadNodes.dimension(0), spaceDim);
185 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException );
186
187 // exception #13 incorrect 0th dimension of output array (must equal number of basis functions)
188 FieldContainer<double> badVals4(quadBasis.getCardinality() + 1, quadNodes.dimension(0));
189 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_DIV), throwCounter, nException );
190
191 // exception #14 incorrect 1st dimension of output array (must equal number of points)
192 FieldContainer<double> badVals5(quadBasis.getCardinality(), quadNodes.dimension(0) + 1, spaceDim);
193 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals5, quadNodes, OPERATOR_VALUE), throwCounter, nException );
194
195 // exception #15 incorrect 1st dimension of output array (must equal number of points)
196 FieldContainer<double> badVals6(quadBasis.getCardinality(), quadNodes.dimension(0) + 1);
197 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals6, quadNodes, OPERATOR_DIV), throwCounter, nException );
198
199 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
200 FieldContainer<double> badVals7(quadBasis.getCardinality(), quadNodes.dimension(0), spaceDim + 1);
201 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals7, quadNodes, OPERATOR_VALUE), throwCounter, nException );
202#endif
203
204 }
205 catch (const std::logic_error & err) {
206 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
207 *outStream << err.what() << '\n';
208 *outStream << "-------------------------------------------------------------------------------" << "\n\n";
209 errorFlag = -1000;
210 };
211
212 // Check if number of thrown exceptions matches the one we expect
213 if (throwCounter != nException) {
214 errorFlag++;
215 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
216 }
217
218 *outStream \
219 << "\n"
220 << "===============================================================================\n"\
221 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
222 << "===============================================================================\n";
223
224 try{
225 std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags();
226
227 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
228 for (unsigned i = 0; i < allTags.size(); i++) {
229 int bfOrd = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
230
231 std::vector<int> myTag = quadBasis.getDofTag(bfOrd);
232 if( !( (myTag[0] == allTags[i][0]) &&
233 (myTag[1] == allTags[i][1]) &&
234 (myTag[2] == allTags[i][2]) &&
235 (myTag[3] == allTags[i][3]) ) ) {
236 errorFlag++;
237 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
238 *outStream << " getDofOrdinal( {"
239 << allTags[i][0] << ", "
240 << allTags[i][1] << ", "
241 << allTags[i][2] << ", "
242 << allTags[i][3] << "}) = " << bfOrd <<" but \n";
243 *outStream << " getDofTag(" << bfOrd << ") = { "
244 << myTag[0] << ", "
245 << myTag[1] << ", "
246 << myTag[2] << ", "
247 << myTag[3] << "}\n";
248 }
249 }
250
251 // Now do the same but loop over basis functions
252 for( int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) {
253 std::vector<int> myTag = quadBasis.getDofTag(bfOrd);
254 int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
255 if( bfOrd != myBfOrd) {
256 errorFlag++;
257 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
258 *outStream << " getDofTag(" << bfOrd << ") = { "
259 << myTag[0] << ", "
260 << myTag[1] << ", "
261 << myTag[2] << ", "
262 << myTag[3] << "} but getDofOrdinal({"
263 << myTag[0] << ", "
264 << myTag[1] << ", "
265 << myTag[2] << ", "
266 << myTag[3] << "} ) = " << myBfOrd << "\n";
267 }
268 }
269 }
270 catch (const std::logic_error & err){
271 *outStream << err.what() << "\n\n";
272 errorFlag = -1000;
273 };
274
275 *outStream \
276 << "\n"
277 << "===============================================================================\n"\
278 << "| TEST 3: correctness of basis function values |\n"\
279 << "===============================================================================\n";
280
281 outStream -> precision(20);
282
283 // VALUE:
284 double basisValues[] = {
285 // first bf, first row of points
286 1.0, 0.0, 0.5, 0.0, 0.0, 0.0,
287 // first bf, second row of points
288 1.0, 0.0, 0.5, 0.0, 0.0, 0.0,
289 // first bf, third row of points
290 1.0, 0.0, 0.5, 0.0, 0.0, 0.0,
291 // second bf, first row of points
292 0.0, 0.0, 0.5, 0.0, 1.0, 0.0,
293 // second bf, second row of points
294 0.0, 0.0, 0.5, 0.0, 1.0, 0.0,
295 // second bf, third row of points
296 0.0, 0.0, 0.5, 0.0, 1.0, 0.0,
297 // third bf, first row of points
298 0.0, 1.0, 0.0, 1.0, 0.0, 1.0,
299 // third bf, second row of points
300 0.0, 0.5, 0.0, 0.5, 0.0, 0.5,
301 // third bf, third row of points
302 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
303 // fourth bf, first row of points
304 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
305 // fourth bf, second row of points
306 0.0, 0.5, 0.0, 0.5, 0.0, 0.5,
307 // fourth bf, third row of points
308 0.0, 1.0, 0.0, 1.0, 0.0, 1.0
309 };
310
311 // DIV: each row gives the 6 correct values of the divergence of the 6 basis functions: (P,F) layout
312 double basisDivs[] = {
313 // bf0 has -0.5 divergence
314 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
315 // bf1 has 0.5 divergencece
316 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
317 // bf2 has -0.5 divergence
318 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
319 // bf3 has 0.5 divergence
320 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
321 };
322
323 try{
324
325 // Dimensions for the output arrays:
326 int numPoints = quadNodes.dimension(0);
327 int numFields = quadBasis.getCardinality();
328 int spaceDim = quadBasis.getBaseCellTopology().getDimension();
329
330 // Generic array for values and curls that will be properly sized before each call
332
333 // Check VALUE of basis functions: resize vals to rank-3 container:
334 vals.resize(numFields, numPoints, spaceDim);
335 quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE);
336 for (int i = 0; i < numFields; i++) {
337 for (int j = 0; j < numPoints; j++) {
338 for (int k = 0; k < spaceDim; k++) {
339 // compute offset for (F,P,D) data layout: indices are F->i, P-> j, D->k
340 int l = i * numPoints * spaceDim + j * spaceDim + k;
341 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
342 errorFlag++;
343 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
344
345 // Output the multi-index of the value where the error is:
346 *outStream << " At multi-index { ";
347 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
348 *outStream << "} computed value: " << vals(i,j,k)
349 << " but reference value: " << basisValues[l] << "\n";
350 }
351 }
352 }
353 }
354
355 // Check DIV of basis function: resize vals to rank-2 container
356 vals.resize(numFields, numPoints);
357 quadBasis.getValues(vals, quadNodes, OPERATOR_DIV);
358 for (int i = 0; i < numFields; i++) {
359 for (int j = 0; j < numPoints; j++) {
360 int l = j + i * numPoints;
361 if (std::abs(vals(i,j) - basisDivs[l]) > INTREPID_TOL) {
362 errorFlag++;
363 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
364
365 // Output the multi-index of the value where the error is:
366 *outStream << " At multi-index { ";
367 *outStream << i << " ";*outStream << j << " ";
368 *outStream << "} computed divergence component: " << vals(i,j)
369 << " but reference divergence component: " << basisDivs[l] << "\n";
370 }
371 }
372 }
373
374 }
375
376 // Catch unexpected errors
377 catch (const std::logic_error & err) {
378 *outStream << err.what() << "\n\n";
379 errorFlag = -1000;
380 };
381
382 if (errorFlag != 0)
383 std::cout << "End Result: TEST FAILED\n";
384 else
385 std::cout << "End Result: TEST PASSED\n";
386
387 // reset format state of std::cout
388 std::cout.copyfmt(oldFormatState);
389
390 return errorFlag;
391}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HDIV_QUAD_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 Quadrilateral 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...