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
50#include "Teuchos_oblackholestream.hpp"
51#include "Teuchos_RCP.hpp"
52#include "Teuchos_GlobalMPISession.hpp"
53
54using namespace std;
55using namespace Intrepid;
56
57#define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
58{ \
59 ++nException; \
60 try { \
61 S ; \
62 } \
63 catch (const std::logic_error & err) { \
64 ++throwCounter; \
65 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
66 *outStream << err.what() << '\n'; \
67 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
68 }; \
69}
70
71int main(int argc, char *argv[]) {
72
73 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
74
75 // This little trick lets us print to std::cout only if
76 // a (dummy) command-line argument is provided.
77 int iprint = argc - 1;
78 Teuchos::RCP<std::ostream> outStream;
79 Teuchos::oblackholestream bhs; // outputs nothing
80 if (iprint > 0)
81 outStream = Teuchos::rcp(&std::cout, false);
82 else
83 outStream = Teuchos::rcp(&bhs, false);
84
85 // Save the format state of the original std::cout.
86 Teuchos::oblackholestream oldFormatState;
87 oldFormatState.copyfmt(std::cout);
88
89 *outStream \
90 << "===============================================================================\n" \
91 << "| |\n" \
92 << "| Unit Test (Basis_HGRAD_HEX_I2_FEM) |\n" \
93 << "| |\n" \
94 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95 << "| 2) Basis values for VALUE, GRAD, and Dk operators |\n" \
96 << "| |\n" \
97 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
98 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
99 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
100 << "| |\n" \
101 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
102 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
103 << "| |\n" \
104 << "===============================================================================\n"\
105 << "| TEST 1: Basis creation, exception testing |\n"\
106 << "===============================================================================\n";
107
108 // Define basis and error flag
110 int errorFlag = 0;
111
112 // Initialize throw counter for exception testing
113 int nException = 0;
114 int throwCounter = 0;
115
116 // Define arrayS containing the 20 nodes of hexahedron<20> topology
117 FieldContainer<double> hexNodes(20, 3);
118 // vertices
119 hexNodes(0, 0) = -1.0; hexNodes(0, 1) = -1.0; hexNodes(0, 2) = -1.0;
120 hexNodes(1, 0) = 1.0; hexNodes(1, 1) = -1.0; hexNodes(1, 2) = -1.0;
121 hexNodes(2, 0) = 1.0; hexNodes(2, 1) = 1.0; hexNodes(2, 2) = -1.0;
122 hexNodes(3, 0) = -1.0; hexNodes(3, 1) = 1.0; hexNodes(3, 2) = -1.0;
123
124 hexNodes(4, 0) = -1.0; hexNodes(4, 1) = -1.0; hexNodes(4, 2) = 1.0;
125 hexNodes(5, 0) = 1.0; hexNodes(5, 1) = -1.0; hexNodes(5, 2) = 1.0;
126 hexNodes(6, 0) = 1.0; hexNodes(6, 1) = 1.0; hexNodes(6, 2) = 1.0;
127 hexNodes(7, 0) = -1.0; hexNodes(7, 1) = 1.0; hexNodes(7, 2) = 1.0;
128
129 // nodes on edges
130 hexNodes(8, 0) = 0.0; hexNodes(8, 1) = -1.0; hexNodes(8, 2) = -1.0;
131 hexNodes(9, 0) = 1.0; hexNodes(9, 1) = 0.0; hexNodes(9, 2) = -1.0;
132 hexNodes(10,0) = 0.0; hexNodes(10,1) = 1.0; hexNodes(10,2) = -1.0;
133 hexNodes(11,0) = -1.0; hexNodes(11,1) = 0.0; hexNodes(11,2) = -1.0;
134 hexNodes(12,0) = -1.0; hexNodes(12,1) = -1.0; hexNodes(12,2) = 0.0;
135 hexNodes(13,0) = 1.0; hexNodes(13,1) = -1.0; hexNodes(13,2) = 0.0;
136 hexNodes(14,0) = 1.0; hexNodes(14,1) = 1.0; hexNodes(14,2) = 0.0;
137 hexNodes(15,0) = -1.0; hexNodes(15,1) = 1.0; hexNodes(15,2) = 0.0;
138 hexNodes(16,0) = 0.0; hexNodes(16,1) = -1.0; hexNodes(16,2) = 1.0;
139 hexNodes(17,0) = 1.0; hexNodes(17,1) = 0.0; hexNodes(17,2) = 1.0;
140 hexNodes(18,0) = 0.0; hexNodes(18,1) = 1.0; hexNodes(18,2) = 1.0;
141 hexNodes(19,0) = -1.0; hexNodes(19,1) = 0.0; hexNodes(19,2) = 1.0;
142
143 // Generic array for the output values; needs to be properly resized depending on the operator type
145
146 try{
147 // exception #1: CURL cannot be applied to scalar functions in 3D
148 // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
149 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0), 4 );
150 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_CURL), throwCounter, nException );
151
152 // exception #2: DIV cannot be applied to scalar functions in 3D
153 // resize vals to rank-2 container with dimensions (num. basis functions, num. points)
154 vals.resize(hexBasis.getCardinality(), hexNodes.dimension(0) );
155 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_DIV), throwCounter, nException );
156
157 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
158 // getDofTag() to access invalid array elements thereby causing bounds check exception
159 // exception #3
160 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(3,10,0), throwCounter, nException );
161 // exception #4
162 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,2,1), throwCounter, nException );
163 // exception #5
164 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
165 // exception #6
166 INTREPID_TEST_COMMAND( hexBasis.getDofTag(21), throwCounter, nException );
167 // exception #7
168 INTREPID_TEST_COMMAND( hexBasis.getDofTag(-1), throwCounter, nException );
169
170#ifdef HAVE_INTREPID_DEBUG
171 // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
172 // exception #8: input points array must be of rank-2
173 FieldContainer<double> badPoints1(4, 5, 3);
174 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
175
176 // exception #9 dimension 1 in the input point array must equal space dimension of the cell
177 FieldContainer<double> badPoints2(4, hexBasis.getBaseCellTopology().getDimension() - 1);
178 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
179
180 // exception #10 output values must be of rank-2 for OPERATOR_VALUE
181 FieldContainer<double> badVals1(4, 3, 1);
182 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
183
184 // exception #11 output values must be of rank-3 for OPERATOR_GRAD
185 FieldContainer<double> badVals2(4, 3);
186 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_GRAD), throwCounter, nException );
187
188 // exception #12 output values must be of rank-3 for OPERATOR_D1
189 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D1), throwCounter, nException );
190
191 // exception #13 output values must be of rank-3 for OPERATOR_D2
192 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_D2), throwCounter, nException );
193
194 // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
195 FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0));
196 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException );
197
198 // exception #15 incorrect 1st dimension of output array (must equal number of points)
199 FieldContainer<double> badVals4(hexBasis.getCardinality(), hexNodes.dimension(0) + 1);
200 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_VALUE), throwCounter, nException );
201
202 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
203 FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0), hexBasis.getBaseCellTopology().getDimension() - 1);
204 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_GRAD), throwCounter, nException );
205
206 // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D)
207 FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0), 40);
208 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_D2), throwCounter, nException );
209
210 // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D)
211 FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 50);
212 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, OPERATOR_D3), throwCounter, nException );
213#endif
214
215 }
216 catch (const std::logic_error & err) {
217 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
218 *outStream << err.what() << '\n';
219 *outStream << "-------------------------------------------------------------------------------" << "\n\n";
220 errorFlag = -1000;
221 };
222
223 // Check if number of thrown exceptions matches the one we expect
224 // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
225 if (throwCounter != nException) {
226 errorFlag++;
227 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
228 }
229
230 *outStream \
231 << "\n"
232 << "===============================================================================\n"\
233 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
234 << "===============================================================================\n";
235
236 try{
237 std::vector<std::vector<int> > allTags = hexBasis.getAllDofTags();
238
239 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
240 for (unsigned i = 0; i < allTags.size(); i++) {
241 int bfOrd = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
242
243 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
244 if( !( (myTag[0] == allTags[i][0]) &&
245 (myTag[1] == allTags[i][1]) &&
246 (myTag[2] == allTags[i][2]) &&
247 (myTag[3] == allTags[i][3]) ) ) {
248 errorFlag++;
249 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
250 *outStream << " getDofOrdinal( {"
251 << allTags[i][0] << ", "
252 << allTags[i][1] << ", "
253 << allTags[i][2] << ", "
254 << allTags[i][3] << "}) = " << bfOrd <<" but \n";
255 *outStream << " getDofTag(" << bfOrd << ") = { "
256 << myTag[0] << ", "
257 << myTag[1] << ", "
258 << myTag[2] << ", "
259 << myTag[3] << "}\n";
260 }
261 }
262
263 // Now do the same but loop over basis functions
264 for( int bfOrd = 0; bfOrd < hexBasis.getCardinality(); bfOrd++) {
265 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
266 int myBfOrd = hexBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
267 if( bfOrd != myBfOrd) {
268 errorFlag++;
269 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
270 *outStream << " getDofTag(" << bfOrd << ") = { "
271 << myTag[0] << ", "
272 << myTag[1] << ", "
273 << myTag[2] << ", "
274 << myTag[3] << "} but getDofOrdinal({"
275 << myTag[0] << ", "
276 << myTag[1] << ", "
277 << myTag[2] << ", "
278 << myTag[3] << "} ) = " << myBfOrd << "\n";
279 }
280 }
281 }
282 catch (const std::logic_error & err){
283 *outStream << err.what() << "\n\n";
284 errorFlag = -1000;
285 };
286
287
288 *outStream \
289 << "\n"
290 << "===============================================================================\n"\
291 << "| TEST 3: correctness of basis function values |\n"\
292 << "===============================================================================\n";
293
294 outStream -> precision(20);
295
296 // VALUE: Each row gives the 8 correct basis set values at an evaluation point
297 double basisValues[] = {
298 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
299 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
300 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
301 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
302 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
303 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
304 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
305 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
306 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
307 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
308 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
309 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, 0, \
310 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, 0, \
311 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, 0, \
312 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, 0, \
313 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, 0, \
314 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, 0, \
315 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, 0, \
316 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0, 0, \
317 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0 };
318
319
320 // GRAD, D1, D2, D3 test values are stored in files due to their large size
321 std::string fileName;
322 std::ifstream dataFile;
323
324
325 // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test
326 std::vector<double> basisGrads; // Flat array for the gradient values.
327
328 fileName = "./testdata/HEX_I2_GradVals.dat";
329 dataFile.open(fileName.c_str());
330 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
331 ">>> ERROR (HGRAD_HEX_I2/test01): could not open GRAD values data file, test aborted.");
332 while (!dataFile.eof() ){
333 double temp;
334 string line; // string for one line of input file
335 std::getline(dataFile, line); // get next line from file
336 stringstream data_line(line); // convert to stringstream
337 while(data_line >> temp){ // extract value from line
338 basisGrads.push_back(temp); // push into vector
339 }
340 }
341 // It turns out that just closing and then opening the ifstream variable does not reset it
342 // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
343 // scope the variables.
344 dataFile.close();
345 dataFile.clear();
346
347
348 //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality)
349 std::vector<double> basisD2;
350 fileName = "./testdata/HEX_I2_D2Vals.dat";
351 dataFile.open(fileName.c_str());
352 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
353 ">>> ERROR (HGRAD_HEX_I2/test01): could not open D2 values data file, test aborted.");
354 while (!dataFile.eof() ){
355 double temp;
356 string line; // string for one line of input file
357 std::getline(dataFile, line); // get next line from file
358 stringstream data_line(line); // convert to stringstream
359 while(data_line >> temp){ // extract value from line
360 basisD2.push_back(temp); // push into vector
361 }
362 }
363 dataFile.close();
364 dataFile.clear();
365
366
367 //D3: flat array with the values of D3 applied to basis functions. Multi-index is (F,P,D3cardinality)
368 std::vector<double> basisD3;
369
370 fileName = "./testdata/HEX_I2_D3Vals.dat";
371 dataFile.open(fileName.c_str());
372 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
373 ">>> ERROR (HGRAD_HEX_I2/test01): could not open D3 values data file, test aborted.");
374
375 while (!dataFile.eof() ){
376 double temp;
377 string line; // string for one line of input file
378 std::getline(dataFile, line); // get next line from file
379 stringstream data_line(line); // convert to stringstream
380 while(data_line >> temp){ // extract value from line
381 basisD3.push_back(temp); // push into vector
382 }
383 }
384 dataFile.close();
385 dataFile.clear();
386
387
388
389 try{
390
391 // Dimensions for the output arrays:
392 int numFields = hexBasis.getCardinality();
393 int numPoints = hexNodes.dimension(0);
394 int spaceDim = hexBasis.getBaseCellTopology().getDimension();
395
396 // Generic array for values, grads, curls, etc. that will be properly sized before each call
398
399 // Check VALUE of basis functions: resize vals to rank-2 container:
400 vals.resize(numFields, numPoints);
401 hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
402 for (int i = 0; i < numFields; i++) {
403 for (int j = 0; j < numPoints; j++) {
404 int l = i + j * numFields;
405 if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
406 errorFlag++;
407 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
408
409 // Output the multi-index of the value where the error is:
410 *outStream << " At multi-index { ";
411 *outStream << i << " ";*outStream << j << " ";
412 *outStream << "} computed value: " << vals(i,j)
413 << " but reference value: " << basisValues[l] << "\n";
414 }
415 }
416 }
417
418
419 // Check GRAD of basis function: resize vals to rank-3 container
420 vals.resize(numFields, numPoints, spaceDim);
421 hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD);
422 for (int i = 0; i < numFields; i++) {
423 for (int j = 0; j < numPoints; j++) {
424 for (int k = 0; k < spaceDim; k++) {
425
426 // basisGrads is (F,P,D), compute offset:
427 int l = k + j * spaceDim + i * spaceDim * numPoints;
428
429 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
430 errorFlag++;
431 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
432
433 // Output the multi-index of the value where the error is:
434 *outStream << " At multi-index { ";
435 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
436 *outStream << "} computed grad component: " << vals(i,j,k)
437 << " but reference grad component: " << basisGrads[l] << "\n";
438 }
439
440 }
441 }
442 }
443
444 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
445 hexBasis.getValues(vals, hexNodes, OPERATOR_D1);
446 for (int i = 0; i < numFields; i++) {
447 for (int j = 0; j < numPoints; j++) {
448 for (int k = 0; k < spaceDim; k++) {
449
450 // basisGrads is (F,P,D), compute offset:
451 int l = k + j * spaceDim + i * spaceDim * numPoints;
452 if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
453 errorFlag++;
454 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
455
456 // Output the multi-index of the value where the error is:
457 *outStream << " At multi-index { ";
458 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
459 *outStream << "} computed D1 component: " << vals(i,j,k)
460 << " but reference D1 component: " << basisGrads[l] << "\n";
461 }
462
463 }
464 }
465 }
466
467
468 // Check D2 of basis function
469 int D2cardinality = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
470 vals.resize(numFields, numPoints, D2cardinality);
471 hexBasis.getValues(vals, hexNodes, OPERATOR_D2);
472 for (int i = 0; i < numFields; i++) {
473 for (int j = 0; j < numPoints; j++) {
474 for (int k = 0; k < D2cardinality; k++) {
475
476 // basisD2 is (F,P,Dk), compute offset:
477 int l = k + j * D2cardinality + i * D2cardinality * numPoints;
478 if (std::abs(vals(i,j,k) - basisD2[l]) > INTREPID_TOL) {
479 errorFlag++;
480 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
481
482 // Output the multi-index of the value where the error is:
483 *outStream << " At multi-index { ";
484 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
485 *outStream << "} computed D2 component: " << vals(i,j,k)
486 << " but reference D2 component: " << basisD2[l] << "\n";
487 }
488 }
489 }
490 }
491
492
493 // Check D3 of basis function
494 int D3cardinality = Intrepid::getDkCardinality(OPERATOR_D3, spaceDim);
495 vals.resize(numFields, numPoints, D3cardinality);
496 hexBasis.getValues(vals, hexNodes, OPERATOR_D3);
497
498 for (int i = 0; i < numFields; i++) {
499 for (int j = 0; j < numPoints; j++) {
500 for (int k = 0; k < D3cardinality; k++) {
501
502 // basisD3 is (F,P,Dk), compute offset:
503 int l = k + j * D3cardinality + i * D3cardinality * numPoints;
504 if (std::abs(vals(i,j,k) - basisD3[l]) > INTREPID_TOL) {
505 errorFlag++;
506 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
507
508 // Output the multi-index of the value where the error is:
509 *outStream << " At multi-index { ";
510 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
511 *outStream << "} computed D3 component: " << vals(i,j,k)
512 << " but reference D3 component: " << basisD3[l] << "\n";
513 }
514 }
515 }
516 }
517
518 }
519
520 // Catch unexpected errors
521 catch (const std::logic_error & err) {
522 *outStream << err.what() << "\n\n";
523 errorFlag = -1000;
524 };
525
526 *outStream \
527 << "\n"
528 << "===============================================================================\n"\
529 << "| TEST 4: correctness of DoF locations |\n"\
530 << "===============================================================================\n";
531
532 try{
533 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis =
534 Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM<double, FieldContainer<double> >);
535 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface =
536 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis);
537
539 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality());
540
541 // Check exceptions.
542#ifdef HAVE_INTREPID_DEBUG
543 cvals.resize(1,2,3);
544 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
545 cvals.resize(3,2);
546 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
547 cvals.resize(20,2);
548 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
549#endif
550 cvals.resize(20,3);
551 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--;
552 // Check if number of thrown exceptions matches the one we expect
553 if (throwCounter != nException) {
554 errorFlag++;
555 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
556 }
557
558 // Check mathematical correctness.
559 basis->getValues(bvals, cvals, OPERATOR_VALUE);
560 char buffer[120];
561 for (int i=0; i<bvals.dimension(0); i++) {
562 for (int j=0; j<bvals.dimension(1); j++) {
563 if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) {
564 errorFlag++;
565 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), bvals(i,j), 0.0);
566 *outStream << buffer;
567 }
568 else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) {
569 errorFlag++;
570 sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), cvals(i,2), bvals(i,j), 1.0);
571 *outStream << buffer;
572 }
573 }
574 }
575
576 }
577 catch (const std::logic_error & err){
578 *outStream << err.what() << "\n\n";
579 errorFlag = -1000;
580 };
581
582 if (errorFlag != 0)
583 std::cout << "End Result: TEST FAILED\n";
584 else
585 std::cout << "End Result: TEST PASSED\n";
586
587 // reset format state of std::cout
588 std::cout.copyfmt(oldFormatState);
589
590 return errorFlag;
591}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_HEX_I2_FEM class.
static const double INTREPID_TOL
General purpose tolerance in, e.g., internal Newton's method to invert ref to phys maps.
int getDkCardinality(const EOperator operatorType, const int spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Implementation of the serendipity-family H(grad)-compatible FEM basis of degree 2 on a Hexahedron cel...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Hexahedron cell.
virtual int getDofOrdinal(const int subcDim, const int subcOrd, const int subcDofOrd)
DoF tag to ordinal lookup.
virtual const std::vector< std::vector< int > > & getAllDofTags()
Retrieves all DoF tags.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual const std::vector< int > & getDofTag(const int dofOrd)
DoF ordinal to DoF tag lookup.
virtual int getCardinality() const
Returns cardinality of the basis.
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.