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_HDIV_HEX_I1_FEM) |\n" \
93 << "| |\n" \
94 << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95 << "| 2) Basis values for VALUE and DIV 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 array containing the 8 vertices of the reference HEX, its center and 6 face centers
117 FieldContainer<double> hexNodes(15, 3);
118 hexNodes(0,0) = -1.0; hexNodes(0,1) = -1.0; hexNodes(0,2) = -1.0;
119 hexNodes(1,0) = 1.0; hexNodes(1,1) = -1.0; hexNodes(1,2) = -1.0;
120 hexNodes(2,0) = 1.0; hexNodes(2,1) = 1.0; hexNodes(2,2) = -1.0;
121 hexNodes(3,0) = -1.0; hexNodes(3,1) = 1.0; hexNodes(3,2) = -1.0;
122
123 hexNodes(4,0) = -1.0; hexNodes(4,1) = -1.0; hexNodes(4,2) = 1.0;
124 hexNodes(5,0) = 1.0; hexNodes(5,1) = -1.0; hexNodes(5,2) = 1.0;
125 hexNodes(6,0) = 1.0; hexNodes(6,1) = 1.0; hexNodes(6,2) = 1.0;
126 hexNodes(7,0) = -1.0; hexNodes(7,1) = 1.0; hexNodes(7,2) = 1.0;
127
128 hexNodes(8,0) = 0.0; hexNodes(8,1) = 0.0; hexNodes(8,2) = 0.0;
129
130 hexNodes(9,0) = 1.0; hexNodes(9,1) = 0.0; hexNodes(9,2) = 0.0;
131 hexNodes(10,0)= -1.0; hexNodes(10,1)= 0.0; hexNodes(10,2)= 0.0;
132
133 hexNodes(11,0)= 0.0; hexNodes(11,1)= 1.0; hexNodes(11,2)= 0.0;
134 hexNodes(12,0)= 0.0; hexNodes(12,1)= -1.0; hexNodes(12,2)= 0.0;
135
136 hexNodes(13,0)= 0.0; hexNodes(13,1)= 0.0; hexNodes(13,2)= 1.0;
137 hexNodes(14,0)= 0.0; hexNodes(14,1)= 0.0; hexNodes(14,2)= -1.0;
138
139
140 // Generic array for the output values; needs to be properly resized depending on the operator type
142
143 try{
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(hexBasis.getCardinality(), hexNodes.dimension(0), 3 );
147 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, OPERATOR_GRAD), throwCounter, nException );
148
149 // exception #2: CURL cannot be applied to HDIV functions
150 INTREPID_TEST_COMMAND( hexBasis.getValues(vals, hexNodes, 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( hexBasis.getDofOrdinal(3,0,0), throwCounter, nException );
156 // exception #4
157 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(1,1,1), throwCounter, nException );
158 // exception #5
159 INTREPID_TEST_COMMAND( hexBasis.getDofOrdinal(0,4,1), throwCounter, nException );
160 // exception #6
161 INTREPID_TEST_COMMAND( hexBasis.getDofTag(12), throwCounter, nException );
162 // exception #7
163 INTREPID_TEST_COMMAND( hexBasis.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( hexBasis.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, 2);
173 INTREPID_TEST_COMMAND( hexBasis.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, 3);
177 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals1, hexNodes, OPERATOR_VALUE), throwCounter, nException );
178
179 // exception #11 output values must be of rank-2 for OPERATOR_DIV
180 FieldContainer<double> badVals2(4, 3, 3);
181 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals2, hexNodes, OPERATOR_DIV), throwCounter, nException );
182
183 // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
184 FieldContainer<double> badVals3(hexBasis.getCardinality() + 1, hexNodes.dimension(0), 3);
185 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals3, hexNodes, OPERATOR_VALUE), throwCounter, nException );
186
187 // exception #13 incorrect 0th dimension of output array (must equal number of basis functions)
188 FieldContainer<double> badVals4(hexBasis.getCardinality() + 1, hexNodes.dimension(0));
189 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals4, hexNodes, OPERATOR_DIV), throwCounter, nException );
190
191 // exception #14 incorrect 1st dimension of output array (must equal number of points)
192 FieldContainer<double> badVals5(hexBasis.getCardinality(), hexNodes.dimension(0) + 1, 3);
193 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals5, hexNodes, OPERATOR_VALUE), throwCounter, nException );
194
195 // exception #15 incorrect 1st dimension of output array (must equal number of points)
196 FieldContainer<double> badVals6(hexBasis.getCardinality(), hexNodes.dimension(0) + 1);
197 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals6, hexNodes, OPERATOR_DIV), throwCounter, nException );
198
199 // exception #16: incorrect 2nd dimension of output array (must equal the space dimension)
200 FieldContainer<double> badVals7(hexBasis.getCardinality(), hexNodes.dimension(0), 4);
201 INTREPID_TEST_COMMAND( hexBasis.getValues(badVals7, hexNodes, 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 = hexBasis.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 = hexBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
230
231 std::vector<int> myTag = hexBasis.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 < hexBasis.getCardinality(); bfOrd++) {
253 std::vector<int> myTag = hexBasis.getDofTag(bfOrd);
254 int myBfOrd = hexBasis.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: Each row pair gives the 6x3 correct basis set values at an evaluation point: (P,F,D) layout
284 double basisValues[] = {
285 // bottom 4 vertices
286 0.,-0.25,0., 0.,0.,0., 0.,0.,0., -0.25,0.,0., 0.,0.,-0.25, 0.,0.,0.,
287 0.,-0.25,0., 0.25,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,-0.25, 0.,0.,0.,
288 0.,0.,0., 0.25,0.,0., 0.,0.25,0., 0.,0.,0., 0.,0.,-0.25, 0.,0.,0.,
289 0.,0.,0., 0.,0.,0., 0.,0.25,0., -0.25,0.,0., 0.,0.,-0.25, 0.,0.,0.,
290 // top 4 vertices
291 0.,-0.25,0., 0.,0.,0., 0.,0.,0., -0.25,0.,0., 0.,0.,0., 0.,0.,0.25,
292 0.,-0.25,0., 0.25,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.25,
293 0.,0.,0., 0.25,0.,0., 0.,0.25,0., 0.,0.,0., 0.,0.,0., 0.,0.,0.25,
294 0.,0.,0., 0.,0.,0., 0.,0.25,0., -0.25,0.,0., 0.,0.,0., 0.,0.,0.25,
295 // center {0, 0, 0}
296 0.,-0.125,0., 0.125,0.,0., 0.,0.125,0., -0.125,0.,0., 0.,0.,-0.125, 0.,0.,0.125,
297 // faces { 1, 0, 0} and {-1, 0, 0}
298 0.,-0.125,0., 0.25,0.,0., 0.,0.125,0., 0.,0.,0., 0.,0.,-0.125, 0.,0.,0.125,
299 0.,-0.125,0., 0.,0.,0., 0.,0.125,0., -0.25,0.,0., 0.,0.,-0.125, 0.,0.,0.125,
300 // faces { 0, 1, 0} and { 0,-1, 0}
301 0.,0.,0., 0.125,0.,0., 0.,0.25,0., -0.125,0.,0., 0.,0.,-0.125, 0.,0.,0.125,
302 0.,-0.25,0., 0.125,0.,0., 0.,0.,0., -0.125,0.,0., 0.,0.,-0.125, 0.,0.,0.125,
303 // faces {0, 0, 1} and {0, 0, -1}
304 0.,-0.125,0., 0.125,0.,0., 0.,0.125,0., -0.125,0.,0., 0.,0.,0., 0.,0.,0.25,
305 0.,-0.125,0., 0.125,0.,0., 0.,0.125,0., -0.125,0.,0., 0.,0.,-0.25, 0.,0.,0.
306 };
307
308 // DIV: each row gives the 6 correct values of the divergence of the 6 basis functions: (P,F) layout
309 double basisDivs[] = {
310 // bottom 4 vertices
311 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
312 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
313 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
314 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
315 // top 4 vertices
316 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
317 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
318 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
319 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
320 // center {0, 0, 0}
321 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
322 // faces { 1, 0, 0} and {-1, 0, 0}
323 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
324 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
325 // faces { 0, 1, 0} and { 0,-1, 0}
326 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
327 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
328 // faces {0, 0, 1} and {0, 0, -1}
329 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
330 0.125, 0.125, 0.125, 0.125, 0.125, 0.125,
331 };
332
333 try{
334
335 // Dimensions for the output arrays:
336 int numPoints = hexNodes.dimension(0);
337 int numFields = hexBasis.getCardinality();
338 int spaceDim = hexBasis.getBaseCellTopology().getDimension();
339
340 // Generic array for values and curls that will be properly sized before each call
342
343 // Check VALUE of basis functions: resize vals to rank-3 container:
344 vals.resize(numFields, numPoints, spaceDim);
345 hexBasis.getValues(vals, hexNodes, OPERATOR_VALUE);
346 for (int i = 0; i < numFields; i++) {
347 for (int j = 0; j < numPoints; j++) {
348 for (int k = 0; k < spaceDim; k++) {
349
350 // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
351 int l = k + i * spaceDim + j * spaceDim * numFields;
352 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
353 errorFlag++;
354 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
355
356 // Output the multi-index of the value where the error is:
357 *outStream << " At multi-index { ";
358 *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
359 *outStream << "} computed value: " << vals(i,j,k)
360 << " but reference value: " << basisValues[l] << "\n";
361 }
362 }
363 }
364 }
365
366 // Check DIV of basis function: resize vals to rank-2 container
367 vals.resize(numFields, numPoints);
368 hexBasis.getValues(vals, hexNodes, OPERATOR_DIV);
369 for (int i = 0; i < numFields; i++) {
370 for (int j = 0; j < numPoints; j++) {
371 int l = i + j * numFields;
372 if (std::abs(vals(i,j) - basisDivs[l]) > INTREPID_TOL) {
373 errorFlag++;
374 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
375
376 // Output the multi-index of the value where the error is:
377 *outStream << " At multi-index { ";
378 *outStream << i << " ";*outStream << j << " ";
379 *outStream << "} computed divergence component: " << vals(i,j)
380 << " but reference divergence component: " << basisDivs[l] << "\n";
381 }
382 }
383 }
384
385 }
386
387 // Catch unexpected errors
388 catch (const std::logic_error & err) {
389 *outStream << err.what() << "\n\n";
390 errorFlag = -1000;
391 };
392
393 *outStream \
394 << "\n"
395 << "===============================================================================\n"\
396 << "| TEST 4: correctness of DoF locations |\n"\
397 << "===============================================================================\n";
398
399 try{
400 Teuchos::RCP<Basis<double, FieldContainer<double> > > basis =
401 Teuchos::rcp(new Basis_HDIV_HEX_I1_FEM<double, FieldContainer<double> >);
402 Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface =
403 Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis);
404
405 int spaceDim = 3;
407 FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality(),spaceDim); // last dimension is spatial dim
408
409 // Check exceptions.
410#ifdef HAVE_INTREPID_DEBUG
411 cvals.resize(1,2,3);
412 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
413 cvals.resize(3,2);
414 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
415 cvals.resize(4,2);
416 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
417#endif
418 cvals.resize(6,spaceDim);
419 INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--;
420 // Check if number of thrown exceptions matches the one we expect
421 if (throwCounter != nException) {
422 errorFlag++;
423 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
424 }
425
426 // Check mathematical correctness
427 FieldContainer<double> normals(basis->getCardinality(),spaceDim); // normals at each point basis point
428 normals(0,0) = 0.0; normals(0,1) = -4.0; normals(0,2) = 0.0;
429 normals(1,0) = 4.0; normals(1,1) = 0.0; normals(1,2) = 0.0;
430 normals(2,0) = 0.0; normals(2,1) = 4.0; normals(2,2) = 0.0;
431 normals(3,0) = -4.0; normals(3,1) = 0.0; normals(3,2) = 0.0;
432 normals(4,0) = 0.0; normals(4,1) = 0.0; normals(4,2) = -4.0;
433 normals(5,0) = 0.0; normals(5,1) = 0.0; normals(5,2) = 4.0;
434
435 basis->getValues(bvals, cvals, OPERATOR_VALUE);
436 char buffer[120];
437 for (int i=0; i<bvals.dimension(0); i++) {
438 for (int j=0; j<bvals.dimension(1); j++) {
439
440 double normal = 0.0;
441 for(int d=0;d<spaceDim;d++)
442 normal += bvals(i,j,d)*normals(j,d);
443
444 if ((i != j) && (std::abs(normal - 0.0) > INTREPID_TOL)) {
445 errorFlag++;
446 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), normal, 0.0);
447 *outStream << buffer;
448 }
449 else if ((i == j) && (std::abs(normal - 1.0) > INTREPID_TOL)) {
450 errorFlag++;
451 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), normal, 1.0);
452 *outStream << buffer;
453 }
454 }
455 }
456
457 }
458 catch (const std::logic_error & err){
459 *outStream << err.what() << "\n\n";
460 errorFlag = -1000;
461 };
462
463 if (errorFlag != 0)
464 std::cout << "End Result: TEST FAILED\n";
465 else
466 std::cout << "End Result: TEST PASSED\n";
467
468 // reset format state of std::cout
469 std::cout.copyfmt(oldFormatState);
470
471 return errorFlag;
472}
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HDIV_HEX_I1_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on Hexahedron cell.
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.