Intrepid
test_01.cpp
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 
54 using namespace std;
55 using 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 
71 int 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_HCURL_WEDGE_I1_FEM) |\n" \
93  << "| |\n" \
94  << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95  << "| 2) Basis values for VALUE and CURL 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 6 vertices of the reference WEDGE and 6 other points.
117  FieldContainer<double> wedgeNodes(12, 3);
118  wedgeNodes(0,0) = 0.0; wedgeNodes(0,1) = 0.0; wedgeNodes(0,2) = -1.0;
119  wedgeNodes(1,0) = 1.0; wedgeNodes(1,1) = 0.0; wedgeNodes(1,2) = -1.0;
120  wedgeNodes(2,0) = 0.0; wedgeNodes(2,1) = 1.0; wedgeNodes(2,2) = -1.0;
121  wedgeNodes(3,0) = 0.0; wedgeNodes(3,1) = 0.0; wedgeNodes(3,2) = 1.0;
122  wedgeNodes(4,0) = 1.0; wedgeNodes(4,1) = 0.0; wedgeNodes(4,2) = 1.0;
123  wedgeNodes(5,0) = 0.0; wedgeNodes(5,1) = 1.0; wedgeNodes(5,2) = 1.0;
124 
125  wedgeNodes(6,0) = 0.25; wedgeNodes(6,1) = 0.5; wedgeNodes(6,2) = -1.0;
126  wedgeNodes(7,0) = 0.5; wedgeNodes(7,1) = 0.25; wedgeNodes(7,2) = 0.0;
127  wedgeNodes(8,0) = 0.25; wedgeNodes(8,1) = 0.25; wedgeNodes(8,2) = 1.0;
128  wedgeNodes(9,0) = 0.25; wedgeNodes(9,1) = 0.0; wedgeNodes(9,2) = 0.75;
129  wedgeNodes(10,0)= 0.0; wedgeNodes(10,1)= 0.5; wedgeNodes(10,2)= -0.25;
130  wedgeNodes(11,0)= 0.5; wedgeNodes(11,1)= 0.5; wedgeNodes(11,2)= 0.0;
131 
132 
133 
134  // Generic array for the output values; needs to be properly resized depending on the operator type
136 
137  try{
138  // exception #1: GRAD cannot be applied to HCURL functions
139  // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
140  vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 3 );
141  INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_GRAD), throwCounter, nException );
142 
143  // exception #2: DIV cannot be applied to HCURL functions
144  // resize vals to rank-2 container with dimensions (num. basis functions, num. points)
145  vals.resize(wedgeBasis.getCardinality(), wedgeNodes.dimension(0));
146  INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_DIV), throwCounter, nException );
147 
148  // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
149  // getDofTag() to access invalid array elements thereby causing bounds check exception
150  // exception #3
151  INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(3,0,0), throwCounter, nException );
152  // exception #4
153  INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(1,1,1), throwCounter, nException );
154  // exception #5
155  INTREPID_TEST_COMMAND( wedgeBasis.getDofOrdinal(0,4,1), throwCounter, nException );
156  // exception #6
157  INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(10), throwCounter, nException );
158  // exception #7
159  INTREPID_TEST_COMMAND( wedgeBasis.getDofTag(-1), throwCounter, nException );
160 
161 #ifdef HAVE_INTREPID_DEBUG
162  // Exceptions 8- test exception handling with incorrectly dimensioned input/output arrays
163  // exception #8: input points array must be of rank-2
164  FieldContainer<double> badPoints1(4, 5, 3);
165  INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
166 
167  // exception #9 dimension 1 in the input point array must equal space dimension of the cell
168  FieldContainer<double> badPoints2(4, 2);
169  INTREPID_TEST_COMMAND( wedgeBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
170 
171  // exception #10 output values must be of rank-3 for OPERATOR_VALUE
172  FieldContainer<double> badVals1(4, 3);
173  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals1, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
174 
175  // exception #11 output values must be of rank-3 for OPERATOR_CURL
176  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals1, wedgeNodes, OPERATOR_CURL), throwCounter, nException );
177 
178  // exception #12 incorrect 0th dimension of output array (must equal number of basis functions)
179  FieldContainer<double> badVals2(wedgeBasis.getCardinality() + 1, wedgeNodes.dimension(0), 3);
180  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals2, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
181 
182  // exception #13 incorrect 1st dimension of output array (must equal number of points)
183  FieldContainer<double> badVals3(wedgeBasis.getCardinality(), wedgeNodes.dimension(0) + 1, 3);
184  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals3, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
185 
186  // exception #14: incorrect 2nd dimension of output array (must equal the space dimension)
187  FieldContainer<double> badVals4(wedgeBasis.getCardinality(), wedgeNodes.dimension(0), 4);
188  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals4, wedgeNodes, OPERATOR_VALUE), throwCounter, nException );
189 
190  // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
191  INTREPID_TEST_COMMAND( wedgeBasis.getValues(badVals4, wedgeNodes, OPERATOR_CURL), throwCounter, nException );
192 #endif
193 
194  }
195  catch (const std::logic_error & err) {
196  *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
197  *outStream << err.what() << '\n';
198  *outStream << "-------------------------------------------------------------------------------" << "\n\n";
199  errorFlag = -1000;
200  };
201 
202  // Check if number of thrown exceptions matches the one we expect
203  // Note Teuchos throw number will not pick up exceptions 3-7 and therefore will not match.
204  if (throwCounter != nException) {
205  errorFlag++;
206  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
207  }
208 
209  *outStream \
210  << "\n"
211  << "===============================================================================\n"\
212  << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
213  << "===============================================================================\n";
214 
215  try{
216  std::vector<std::vector<int> > allTags = wedgeBasis.getAllDofTags();
217 
218  // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
219  for (unsigned i = 0; i < allTags.size(); i++) {
220  int bfOrd = wedgeBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
221 
222  std::vector<int> myTag = wedgeBasis.getDofTag(bfOrd);
223  if( !( (myTag[0] == allTags[i][0]) &&
224  (myTag[1] == allTags[i][1]) &&
225  (myTag[2] == allTags[i][2]) &&
226  (myTag[3] == allTags[i][3]) ) ) {
227  errorFlag++;
228  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
229  *outStream << " getDofOrdinal( {"
230  << allTags[i][0] << ", "
231  << allTags[i][1] << ", "
232  << allTags[i][2] << ", "
233  << allTags[i][3] << "}) = " << bfOrd <<" but \n";
234  *outStream << " getDofTag(" << bfOrd << ") = { "
235  << myTag[0] << ", "
236  << myTag[1] << ", "
237  << myTag[2] << ", "
238  << myTag[3] << "}\n";
239  }
240  }
241 
242  // Now do the same but loop over basis functions
243  for( int bfOrd = 0; bfOrd < wedgeBasis.getCardinality(); bfOrd++) {
244  std::vector<int> myTag = wedgeBasis.getDofTag(bfOrd);
245  int myBfOrd = wedgeBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
246  if( bfOrd != myBfOrd) {
247  errorFlag++;
248  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
249  *outStream << " getDofTag(" << bfOrd << ") = { "
250  << myTag[0] << ", "
251  << myTag[1] << ", "
252  << myTag[2] << ", "
253  << myTag[3] << "} but getDofOrdinal({"
254  << myTag[0] << ", "
255  << myTag[1] << ", "
256  << myTag[2] << ", "
257  << myTag[3] << "} ) = " << myBfOrd << "\n";
258  }
259  }
260  }
261  catch (const std::logic_error & err){
262  *outStream << err.what() << "\n\n";
263  errorFlag = -1000;
264  };
265 
266  *outStream \
267  << "\n"
268  << "===============================================================================\n"\
269  << "| TEST 3: correctness of basis function values |\n"\
270  << "===============================================================================\n";
271 
272  outStream -> precision(20);
273 
274  // VALUE: Each row pair gives the 9x3 correct basis set values at an evaluation point: (P,F,D) layout
275  double basisValues[] = {
276  1.00000, 0, 0, 0, 0, 0, 0, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
277  0, 0.500000, 0, 0, 0, 0, 0, 0, 1.00000, 1.00000, 0, 0, 1.00000, 0, 0, \
278  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, \
279  0, 0, -1.00000, 0, 0, -1.00000, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
280  0, 0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
281  1.00000, 0, 0, 0, 0, 0, 0, -1.00000, 0, 0, 0, 0.500000, 0, 0, 0, 0, \
282  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.00000, 1.00000, 0, 0, 1.00000, 0, \
283  0, 0, 0, 0, 0, 0, 0, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, \
284  0, 0, 0, -1.00000, 0, 0, -1.00000, -1.00000, 0, 0, 0, 0, 0, 0, 0, 0, \
285  0, 0.500000, 0.500000, 0.250000, 0, -0.500000, 0.250000, 0, \
286  -0.500000, -0.750000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.125000, \
287  0, 0, 0.125000, 0, 0, 0.250000, 0.375000, 0.250000, 0, -0.125000, \
288  0.250000, 0, -0.125000, -0.250000, 0, 0.375000, 0.250000, 0, \
289  -0.125000, 0.250000, 0, -0.125000, -0.250000, 0, 0, 0, 0.125000, 0, \
290  0, 0.250000, 0, 0, 0.125000, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.750000, \
291  0.250000, 0, -0.250000, 0.250000, 0, -0.250000, -0.750000, 0, 0, 0, \
292  0.250000, 0, 0, 0.125000, 0, 0, 0.125000, 0.125000, 0.0312500, 0, 0, \
293  0.0312500, 0, 0, -0.0937500, 0, 0.875000, 0.218750, 0, 0, 0.218750, \
294  0, 0, -0.656250, 0, 0, 0, 0.375000, 0, 0, 0.125000, 0, 0, 0, \
295  0.312500, 0, 0, -0.312500, 0, 0, -0.312500, -0.625000, 0, 0.187500, \
296  0, 0, -0.187500, 0, 0, -0.187500, -0.375000, 0, 0, 0, 0.250000, 0, 0, \
297  0, 0, 0, 0.250000, 0.250000, 0.250000, 0, -0.250000, 0.250000, 0, \
298  -0.250000, -0.250000, 0, 0.250000, 0.250000, 0, -0.250000, 0.250000, \
299  0, -0.250000, -0.250000, 0, 0, 0, 0, 0, 0, 0.250000, 0, 0, 0.250000
300  };
301 
302  // CURL: each row pair gives the 9x3 correct values of the curls of the 9 basis functions: (P,F,D) layout
303  double basisCurls[] = {
304  0, -0.500000, 2.00000, 0, 0, 2.00000, -0.500000, 0, 2.00000, 0, \
305  0.500000, 0, 0, 0, 0, 0.500000, 0, 0, -0.500000, 0.500000, 0, 0, \
306  -0.500000, 0, 0.500000, 0, 0, 0.500000, -0.500000, 2.00000, 0.500000, \
307  0, 2.00000, 0, 0, 2.00000, -0.500000, 0.500000, 0, -0.500000, 0, 0, \
308  0, 0, 0, -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0, \
309  0, 2.00000, 0, 0.500000, 2.00000, -0.500000, 0.500000, 2.00000, 0, 0, \
310  0, 0, -0.500000, 0, 0.500000, -0.500000, 0, -0.500000, 0.500000, 0, \
311  0, -0.500000, 0, 0.500000, 0, 0, 0, -0.500000, 0, 0, 0, 0, -0.500000, \
312  0, 0, 0, 0.500000, 2.00000, 0, 0, 2.00000, 0.500000, 0, 2.00000, \
313  -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0.500000, \
314  -0.500000, 0, 0.500000, 0, 0, 0, 0, 0, -0.500000, 0.500000, 2.00000, \
315  -0.500000, 0, 2.00000, 0, 0, 2.00000, -0.500000, 0.500000, 0, 0, \
316  -0.500000, 0, 0.500000, 0, 0, 0, 0, 0, 0, 0.500000, 0, -0.500000, \
317  0.500000, 0, 0, 0, 2.00000, 0, -0.500000, 2.00000, 0.500000, \
318  -0.500000, 2.00000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \
319  0.500000, 0, 0, 0.125000, -0.250000, 2.00000, 0.125000, 0.250000, \
320  2.00000, -0.375000, 0.250000, 2.00000, -0.125000, 0.250000, 0, \
321  -0.125000, -0.250000, 0, 0.375000, -0.250000, 0, -0.500000, 0.500000, \
322  0, 0, -0.500000, 0, 0.500000, 0, 0, 0.250000, -0.375000, 1.00000, \
323  0.250000, 0.125000, 1.00000, -0.250000, 0.125000, 1.00000, -0.250000, \
324  0.375000, 1.00000, -0.250000, -0.125000, 1.00000, 0.250000, \
325  -0.125000, 1.00000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \
326  0.500000, 0, 0, 0.125000, -0.375000, 0, 0.125000, 0.125000, 0, \
327  -0.375000, 0.125000, 0, -0.125000, 0.375000, 2.00000, -0.125000, \
328  -0.125000, 2.00000, 0.375000, -0.125000, 2.00000, -0.500000, \
329  0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0.125000, -0.500000, \
330  0.250000, 0.125000, 0, 0.250000, -0.375000, 0, 0.250000, -0.125000, \
331  0.500000, 1.75000, -0.125000, 0, 1.75000, 0.375000, 0, 1.75000, \
332  -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0, 0, \
333  -0.250000, 1.25000, 0, 0.250000, 1.25000, -0.500000, 0.250000, \
334  1.25000, 0, 0.250000, 0.750000, 0, -0.250000, 0.750000, 0.500000, \
335  -0.250000, 0.750000, -0.500000, 0.500000, 0, 0, -0.500000, 0, \
336  0.500000, 0, 0, 0.250000, -0.250000, 1.00000, 0.250000, 0.250000, \
337  1.00000, -0.250000, 0.250000, 1.00000, -0.250000, 0.250000, 1.00000, \
338  -0.250000, -0.250000, 1.00000, 0.250000, -0.250000, 1.00000, \
339  -0.500000, 0.500000, 0, 0, -0.500000, 0, 0.500000, 0, 0
340  };
341 
342  try{
343 
344  // Dimensions for the output arrays:
345  int numFields = wedgeBasis.getCardinality();
346  int numPoints = wedgeNodes.dimension(0);
347  int spaceDim = wedgeBasis.getBaseCellTopology().getDimension();
348 
349  // Generic array for values and curls that will be properly sized before each call
351 
352  // Check VALUE of basis functions: resize vals to rank-3 container:
353  vals.resize(numFields, numPoints, spaceDim);
354  wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_VALUE);
355  for (int i = 0; i < numFields; i++) {
356  for (int j = 0; j < numPoints; j++) {
357  for (int k = 0; k < spaceDim; k++) {
358 
359  // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
360  int l = k + i * spaceDim + j * spaceDim * numFields;
361  if (std::abs(vals(i,j,k) - basisValues[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 << " ";*outStream << k << " ";
368  *outStream << "} computed value: " << vals(i,j,k)
369  << " but reference value: " << basisValues[l] << "\n";
370  }
371  }
372  }
373  }
374 
375  // Check CURL of basis function: resize vals to rank-3 container
376  vals.resize(numFields, numPoints, spaceDim);
377  wedgeBasis.getValues(vals, wedgeNodes, OPERATOR_CURL);
378  for (int i = 0; i < numFields; i++) {
379  for (int j = 0; j < numPoints; j++) {
380  for (int k = 0; k < spaceDim; k++) {
381 
382  // compute offset for (P,F,D) data layout: indices are P->j, F->i, D->k
383  int l = k + i * spaceDim + j * spaceDim * numFields;
384  if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) {
385  errorFlag++;
386  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
387 
388  // Output the multi-index of the value where the error is:
389  *outStream << " At multi-index { ";
390  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
391  *outStream << "} computed curl component: " << vals(i,j,k)
392  << " but reference curl component: " << basisCurls[l] << "\n";
393  }
394  }
395  }
396  }
397 
398  }
399 
400  // Catch unexpected errors
401  catch (const std::logic_error & err) {
402  *outStream << err.what() << "\n\n";
403  errorFlag = -1000;
404  };
405 
406  if (errorFlag != 0)
407  std::cout << "End Result: TEST FAILED\n";
408  else
409  std::cout << "End Result: TEST PASSED\n";
410 
411  // reset format state of std::cout
412  std::cout.copyfmt(oldFormatState);
413 
414  return errorFlag;
415 }
virtual const std::vector< int > & getDofTag(const int dofOrd)
DoF ordinal to DoF tag lookup.
virtual int getCardinality() const
Returns cardinality of the basis.
virtual const std::vector< std::vector< int > > & getAllDofTags()
Retrieves all DoF tags.
Header file for utility class to provide multidimensional containers.
Implementation of the default H(curl)-compatible FEM basis of degree 1 on Wedge cell.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Wedge cell.
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
virtual int getDofOrdinal(const int subcDim, const int subcOrd, const int subcDofOrd)
DoF tag to ordinal lookup.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
Header file for the Intrepid::HCURL_WEDGE_I1_FEM class.