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_PYR_I2_FEM) |\n" \
93 << "| |\n" \
94 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95 << "| 2) Basis values for VALUE, GRAD, CURL, 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
109 Basis_HGRAD_PYR_I2_FEM<double, FieldContainer<double> > pyrBasis;
110 int errorFlag = 0;
111
112 // Initialize throw counter for exception testing
113 int nException = 0;
114 int throwCounter = 0;
115
116 // Define array containing the 13 nodes of the reference PYRAMID and some other points.
117 FieldContainer<double> pyrNodes(18, 3);
118 pyrNodes(0,0) = -1.0; pyrNodes(0,1) = -1.0; pyrNodes(0,2) = 0;
119 pyrNodes(1,0) = 1.0; pyrNodes(1,1) = -1.0; pyrNodes(1,2) = 0;
120 pyrNodes(2,0) = 1.0; pyrNodes(2,1) = 1.0; pyrNodes(2,2) = 0;
121 pyrNodes(3,0) = -1.0; pyrNodes(3,1) = 1.0; pyrNodes(3,2) = 0;
122 pyrNodes(4,0) = 0.0; pyrNodes(4,1) = 0.0; pyrNodes(4,2) = 1.0;
123 pyrNodes(5,0) = 0.0; pyrNodes(5,1) = -1.0; pyrNodes(5,2) = 0.0;
124 pyrNodes(6,0) = 1.0; pyrNodes(6,1) = 0.0; pyrNodes(6,2) = 0.0;
125 pyrNodes(7,0) = 0.0; pyrNodes(7,1) = 1.0; pyrNodes(7,2) = 0.0;
126 pyrNodes(8,0) = -1.0; pyrNodes(8,1) = 0.0; pyrNodes(8,2) = 0.0;
127 pyrNodes(9,0) = -0.5; pyrNodes(9,1) = -0.5; pyrNodes(9,2) = 0.5;
128 pyrNodes(10,0) = 0.5; pyrNodes(10,1) = -0.5; pyrNodes(10,2) = 0.5;
129 pyrNodes(11,0) = 0.5; pyrNodes(11,1) = 0.5; pyrNodes(11,2) = 0.5;
130 pyrNodes(12,0) = -0.5; pyrNodes(12,1) = 0.5; pyrNodes(12,2) = 0.5;
131
132 pyrNodes(13,0) = 0.25; pyrNodes(13,1) = 0.5; pyrNodes(13,2) = 0.2;
133 pyrNodes(14,0) = -0.7 ; pyrNodes(14,1) = 0.3; pyrNodes(14,2) = 0.3;
134 pyrNodes(15,0) = 0.; pyrNodes(15,1) = -0.05; pyrNodes(15,2) = 0.95;
135 pyrNodes(16,0) = -0.15; pyrNodes(16,1) = -0.2; pyrNodes(16,2) = 0.75;
136 pyrNodes(17,0) = -0.4; pyrNodes(17,1) = 0.9; pyrNodes(17,2) = 0.0;
137
138
139 try{
140 // Generic array for the output values; needs to be properly resized depending on the operator type
141 FieldContainer<double> vals;
142
143 // exception #1: CURL cannot be applied to scalar functions
144 // resize vals to rank-3 container with dimensions (num. points, num. basis functions)
145 vals.resize(pyrBasis.getCardinality(), pyrNodes.dimension(0), 3 );
146 INTREPID_TEST_COMMAND( pyrBasis.getValues(vals, pyrNodes, OPERATOR_CURL), throwCounter, nException );
147
148 // exception #2: DIV cannot be applied to scalar functions
149 // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
150 vals.resize(pyrBasis.getCardinality(), pyrNodes.dimension(0) );
151 INTREPID_TEST_COMMAND( pyrBasis.getValues(vals, pyrNodes, OPERATOR_DIV), throwCounter, nException );
152
153 // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
154 // getDofTag() to access invalid array elements thereby causing bounds check exception
155 // exception #3
156 INTREPID_TEST_COMMAND( pyrBasis.getDofOrdinal(3,0,0), throwCounter, nException );
157 // exception #4
158 INTREPID_TEST_COMMAND( pyrBasis.getDofOrdinal(1,1,1), throwCounter, nException );
159 // exception #5
160 INTREPID_TEST_COMMAND( pyrBasis.getDofOrdinal(0,6,0), throwCounter, nException );
161 // exception #6
162 INTREPID_TEST_COMMAND( pyrBasis.getDofTag(14), throwCounter, nException );
163 // exception #7
164 INTREPID_TEST_COMMAND( pyrBasis.getDofTag(-1), throwCounter, nException );
165
166#ifdef HAVE_INTREPID_DEBUG
167 // Exceptions 8-18 test exception handling with incorrectly dimensioned input/output arrays
168 // exception #8: input points array must be of rank-2
169 FieldContainer<double> badPoints1(4, 5, 3);
170 INTREPID_TEST_COMMAND( pyrBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
171
172 // exception #9 dimension 1 in the input point array must equal space dimension of the cell
173 FieldContainer<double> badPoints2(4, 2);
174 INTREPID_TEST_COMMAND( pyrBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
175
176 // exception #10 output values must be of rank-2 for OPERATOR_VALUE
177 FieldContainer<double> badVals1(4, 3, 1);
178 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals1, pyrNodes, OPERATOR_VALUE), throwCounter, nException );
179
180 // exception #11 output values must be of rank-3 for OPERATOR_GRAD
181 FieldContainer<double> badVals2(4, 3);
182 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals2, pyrNodes, OPERATOR_GRAD), throwCounter, nException );
183
184 // exception #12 output values must be of rank-3 for OPERATOR_D1
185 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals2, pyrNodes, OPERATOR_D1), throwCounter, nException );
186
187 // exception #13 output values must be of rank-3 for OPERATOR_D2
188 //INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals2, pyrNodes, OPERATOR_D2), throwCounter, nException );
189
190 // exception #14 incorrect 0th dimension of output array (must equal number of basis functions)
191 FieldContainer<double> badVals3(pyrBasis.getCardinality() + 1, pyrNodes.dimension(0));
192 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals3, pyrNodes, OPERATOR_VALUE), throwCounter, nException );
193
194 // exception #15 incorrect 1st dimension of output array (must equal number of points)
195 FieldContainer<double> badVals4(pyrBasis.getCardinality(), pyrNodes.dimension(0) + 1);
196 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals4, pyrNodes, OPERATOR_VALUE), throwCounter, nException );
197
198 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
199 FieldContainer<double> badVals5(pyrBasis.getCardinality(), pyrNodes.dimension(0), 4);
200 INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals5, pyrNodes, OPERATOR_GRAD), throwCounter, nException );
201
202 // exception #17: incorrect 2nd dimension of output array (must equal D2 cardinality in 3D)
203 //FieldContainer<double> badVals6(pyrBasis.getCardinality(), pyrNodes.dimension(0), 40);
204 //INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals6, pyrNodes, OPERATOR_D2), throwCounter, nException );
205
206 // exception #18: incorrect 2nd dimension of output array (must equal D3 cardinality in 3D)
207 //INTREPID_TEST_COMMAND( pyrBasis.getValues(badVals6, pyrNodes, OPERATOR_D3), throwCounter, nException );
208#endif
209
210 }
211 catch (const std::logic_error & err) {
212 *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
213 *outStream << err.what() << '\n';
214 *outStream << "-------------------------------------------------------------------------------" << "\n\n";
215 errorFlag = -1000;
216 };
217
218 // Check if number of thrown exceptions matches the one we expect - 18
219 if (throwCounter != nException) {
220 errorFlag++;
221 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
222 }
223
224 *outStream \
225 << "\n"
226 << "===============================================================================\n"\
227 << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
228 << "===============================================================================\n";
229
230 try{
231 std::vector<std::vector<int> > allTags = pyrBasis.getAllDofTags();
232
233 // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
234 for (unsigned i = 0; i < allTags.size(); i++) {
235 int bfOrd = pyrBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
236
237 std::vector<int> myTag = pyrBasis.getDofTag(bfOrd);
238 if( !( (myTag[0] == allTags[i][0]) &&
239 (myTag[1] == allTags[i][1]) &&
240 (myTag[2] == allTags[i][2]) &&
241 (myTag[3] == allTags[i][3]) ) ) {
242 errorFlag++;
243 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
244 *outStream << " getDofOrdinal( {"
245 << allTags[i][0] << ", "
246 << allTags[i][1] << ", "
247 << allTags[i][2] << ", "
248 << allTags[i][3] << "}) = " << bfOrd <<" but \n";
249 *outStream << " getDofTag(" << bfOrd << ") = { "
250 << myTag[0] << ", "
251 << myTag[1] << ", "
252 << myTag[2] << ", "
253 << myTag[3] << "}\n";
254 }
255 }
256
257 // Now do the same but loop over basis functions
258 for( int bfOrd = 0; bfOrd < pyrBasis.getCardinality(); bfOrd++) {
259 std::vector<int> myTag = pyrBasis.getDofTag(bfOrd);
260 int myBfOrd = pyrBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
261 if( bfOrd != myBfOrd) {
262 errorFlag++;
263 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
264 *outStream << " getDofTag(" << bfOrd << ") = { "
265 << myTag[0] << ", "
266 << myTag[1] << ", "
267 << myTag[2] << ", "
268 << myTag[3] << "} but getDofOrdinal({"
269 << myTag[0] << ", "
270 << myTag[1] << ", "
271 << myTag[2] << ", "
272 << myTag[3] << "} ) = " << myBfOrd << "\n";
273 }
274 }
275 }
276 catch (const std::logic_error & err){
277 *outStream << err.what() << "\n\n";
278 errorFlag = -1000;
279 };
280
281 *outStream \
282 << "\n"
283 << "===============================================================================\n"\
284 << "| TEST 3: correctness of basis function values |\n"\
285 << "===============================================================================\n";
286
287 outStream -> precision(20);
288
289 // VAL, GRAD, D2 test values are stored in files due to their large size
290 std::string fileName;
291 std::ifstream dataFile;
292
293 // Values are stored in (F,P) format in a data file. Read file and do the test
294 std::vector<double> basisValues; // Flat array for the basis values.
295 fileName = "./testdata/PYR_I2_Vals.dat";
296 dataFile.open(fileName.c_str());
297 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
298 ">>> ERROR (HGRAD_PYR_I2/test01): could not open basis values data file, test aborted.");
299 while (!dataFile.eof() ){
300 double temp;
301 string line; // string for one line of input file
302 std::getline(dataFile, line); // get next line from file
303 stringstream data_line(line); // convert to stringstream
304 while(data_line >> temp){ // extract value from line
305 basisValues.push_back(temp); // push into vector
306 }
307 }
308 // It turns out that just closing and then opening the ifstream variable does not reset it
309 // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
310 // scope the variables.
311 dataFile.close();
312 dataFile.clear();
313
314
315 // GRAD and D1 values are stored in (F,P,D) format in a data file. Read file and do the test
316 std::vector<double> basisGrads; // Flat array for the gradient values.
317 fileName = "./testdata/PYR_I2_GradVals.dat";
318 dataFile.open(fileName.c_str());
319 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
320 ">>> ERROR (HGRAD_PYR_I2/test01): could not open GRAD values data file, test aborted.");
321 while (!dataFile.eof() ){
322 double temp;
323 string line; // string for one line of input file
324 std::getline(dataFile, line); // get next line from file
325 stringstream data_line(line); // convert to stringstream
326 while(data_line >> temp){ // extract value from line
327 basisGrads.push_back(temp); // push into vector
328 }
329 }
330 // It turns out that just closing and then opening the ifstream variable does not reset it
331 // and subsequent open() command fails. One fix is to explicitely clear the ifstream, or
332 // scope the variables.
333 dataFile.close();
334 dataFile.clear();
335
336 //D2: flat array with the values of D2 applied to basis functions. Multi-index is (F,P,D2cardinality)
337 std::vector<double> basisD2;
338
339 fileName = "./testdata/PYR_I2_D2Vals.dat";
340 dataFile.open(fileName.c_str());
341 TEUCHOS_TEST_FOR_EXCEPTION( !dataFile.good(), std::logic_error,
342 ">>> ERROR (HGRAD_PYR_I2/test01): could not open D2 values data file, test aborted.");
343
344 while (!dataFile.eof() ){
345 double temp;
346 string line; // string for one line of input file
347 std::getline(dataFile, line); // get next line from file
348 stringstream data_line(line); // convert to stringstream
349 while(data_line >> temp){ // extract value from line
350 basisD2.push_back(temp); // push into vector
351 }
352 }
353 dataFile.close();
354 dataFile.clear();
355
356 try{
357
358 // Dimensions for the output arrays:
359 int numFields = pyrBasis.getCardinality();
360 int numPoints = pyrNodes.dimension(0);
361 int spaceDim = pyrBasis.getBaseCellTopology().getDimension();
362 int D2Cardin = Intrepid::getDkCardinality(OPERATOR_D2, spaceDim);
363
364 // Generic array for values, grads, curls, etc. that will be properly sized before each call
365 FieldContainer<double> vals;
366
367 // Check VALUE of basis functions: resize vals to rank-2 container:
368 vals.resize(numFields, numPoints);
369 pyrBasis.getValues(vals, pyrNodes, OPERATOR_VALUE);
370 for (int i = 0; i < numFields; i++) {
371 for (int j = 0; j < numPoints; j++) {
372 int l = j + i * numPoints;
373
374 if (std::abs(vals(i,j) - basisValues[l]) > 100.*INTREPID_TOL) {
375 errorFlag++;
376 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
377
378 // Output the multi-index of the value where the error is:
379 *outStream << " At multi-index { ";
380 *outStream << i << " ";*outStream << j << " ";
381 *outStream << "} computed value: " << vals(i,j)
382 << " but reference value: " << basisValues[l] << "\n";
383
384 }
385 }
386 }
387
388
389 // Check GRAD of basis function: resize vals to rank-3 container
390 vals.resize(numFields, numPoints, spaceDim);
391 pyrBasis.getValues(vals, pyrNodes, OPERATOR_GRAD);
392 for (int i = 0; i < numFields; i++) {
393 for (int j = 0; j < numPoints; j++) {
394 for (int k = 0; k < spaceDim; k++) {
395 //int l = k + h * spaceDim + j * spaceDim * numFields;
396 int l = i + j * numFields + k * numFields * numPoints;
397
398 if (std::abs(vals(i,j,k) - basisGrads[l]) > 100.*INTREPID_TOL) {
399 errorFlag++;
400 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
401
402 // Output the multi-index of the value where the error is:
403 *outStream << " At multi-index { ";
404 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
405 *outStream << "} computed grad component: " << vals(i,j,k)
406 << " but reference grad component: " << basisGrads[l] << "\n";
407
408 }
409
410 }
411 }
412 }
413
414 // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
415 pyrBasis.getValues(vals, pyrNodes, OPERATOR_D1);
416 for (int i = 0; i < numFields; i++) {
417 for (int j = 0; j < numPoints; j++) {
418 for (int k = 0; k < spaceDim; k++) {
419 int l = i + j * numFields + k * numFields * numPoints;
420
421 if (std::abs(vals(i,j,k) - basisGrads[l]) > 100.*INTREPID_TOL) {
422 errorFlag++;
423 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
424
425 // Output the multi-index of the value where the error is:
426 *outStream << " At multi-index { ";
427 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
428 *outStream << "} computed D1 component: " << vals(i,j,k)
429 << " but reference D1 component: " << basisGrads[l] << "\n";
430 }
431
432 }
433 }
434 }
435
436
437 // Check D2 of basis function
438 vals.resize(numFields, numPoints, D2Cardin);
439 pyrBasis.getValues(vals, pyrNodes, OPERATOR_D2);
440 //int l=0;
441 for (int i = 0; i < numFields; i++) {
442 for (int j = 0; j < numPoints; j++) {
443 if(j == 4) continue; //derivatives are singular when z = 1.
444 for (int k = 0; k < D2Cardin; k++) {
445 // int l = k + i * D2Cardin + j * D2Cardin * numFields;
446 int l = i + j * numFields + k * numFields * numPoints;
447 //int l = k + j * D2Cardin + i * D2Cardin * numPoints;
448
449
450 if (std::abs(vals(i,j,k) - basisD2[l]) > 100.*INTREPID_TOL) {
451 errorFlag++;
452 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
453
454 // Output the multi-index of the value where the error is:
455 *outStream << " At multi-index { ";
456 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
457 *outStream << "} computed D2 component: " << vals(i,j,k)
458 << " but reference D2 component: " << basisD2[l] << "\n";
459 }
460
461 // l++;
462
463 // *outStream << vals(i,j,k) << "\n";
464
465
466 }
467 }
468 }
469 }
470
471
472 // Catch unexpected errors
473 catch (const std::logic_error & err) {
474 *outStream << err.what() << "\n\n";
475 errorFlag = -1000;
476 };
477
478 if (errorFlag != 0)
479 std::cout << "End Result: TEST FAILED\n";
480 else
481 std::cout << "End Result: TEST PASSED\n";
482
483 // reset format state of std::cout
484 std::cout.copyfmt(oldFormatState);
485
486 return errorFlag;
487}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_PYR_I2_FEM class.
int getDkCardinality(const EOperator operatorType, const int spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.