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 
55 using namespace std;
56 using namespace Intrepid;
57 
58 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
59 { \
60  ++nException; \
61  try { \
62  S ; \
63  } \
64  catch (std::logic_error err) { \
65  ++throwCounter; \
66  *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
67  *outStream << err.what() << '\n'; \
68  *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
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_HGRAD_TRI_C1_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 (dridzal@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 
109  // Define basis and error flag
111  int errorFlag = 0;
112 
113  // Initialize throw counter for exception testing
114  int nException = 0;
115  int throwCounter = 0;
116 
117  // Define array containing the 3 vertices of the reference Triangle, its center and another point
118  FieldContainer<double> triNodes(5, 2);
119  triNodes(0,0) = 0.0; triNodes(0,1) = 0.0;
120  triNodes(1,0) = 1.0; triNodes(1,1) = 0.0;
121  triNodes(2,0) = 0.0; triNodes(2,1) = 1.0;
122  triNodes(3,0) = 0.5; triNodes(3,1) = 0.5;
123  triNodes(4,0) = 0.0; triNodes(4,1) = 0.75;
124 
125  // Generic array for the output values; needs to be properly resized depending on the operator type
127 
128  try{
129  // exception #1: DIV cannot be applied to scalar functions
130  // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
131  vals.resize(triBasis.getCardinality(), triNodes.dimension(0) );
132  INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_DIV), throwCounter, nException );
133 
134  // Exceptions 2-6: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
135  // getDofTag() to access invalid array elements thereby causing bounds check exception
136  // exception #2
137  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(2,0,0), throwCounter, nException );
138  // exception #3
139  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException );
140  // exception #4
141  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,4,0), throwCounter, nException );
142  // exception #5
143  INTREPID_TEST_COMMAND( triBasis.getDofTag(5), throwCounter, nException );
144  // exception #6
145  INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException );
146 
147 #ifdef HAVE_INTREPID_DEBUG
148  // Exceptions 7-17 test exception handling with incorrectly dimensioned input/output arrays
149  // exception #7: input points array must be of rank-2
150  FieldContainer<double> badPoints1(4, 5, 3);
151  INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
152 
153  // exception #8 dimension 1 in the input point array must equal space dimension of the cell
154  FieldContainer<double> badPoints2(4, 3);
155  INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
156 
157  // exception #9 output values must be of rank-2 for OPERATOR_VALUE
158  FieldContainer<double> badVals1(4, 3, 1);
159  INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException );
160 
161  // exception #10 output values must be of rank-3 for OPERATOR_GRAD
162  FieldContainer<double> badVals2(4, 3);
163  INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_GRAD), throwCounter, nException );
164 
165  // exception #11 output values must be of rank-3 for OPERATOR_CURL
166  INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_CURL), throwCounter, nException );
167 
168  // exception #12 output values must be of rank-3 for OPERATOR_D2
169  INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_D2), throwCounter, nException );
170 
171  // exception #13 incorrect 1st dimension of output array (must equal number of basis functions)
172  FieldContainer<double> badVals3(triBasis.getCardinality() + 1, triNodes.dimension(0));
173  INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException );
174 
175  // exception #14 incorrect 0th dimension of output array (must equal number of points)
176  FieldContainer<double> badVals4(triBasis.getCardinality(), triNodes.dimension(0) + 1);
177  INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_VALUE), throwCounter, nException );
178 
179  // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
180  FieldContainer<double> badVals5(triBasis.getCardinality(), triNodes.dimension(0), 4);
181  INTREPID_TEST_COMMAND( triBasis.getValues(badVals5, triNodes, OPERATOR_GRAD), throwCounter, nException );
182 
183  // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 2D)
184  FieldContainer<double> badVals6(triBasis.getCardinality(), triNodes.dimension(0), 40);
185  INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D2), throwCounter, nException );
186 
187  // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 2D)
188  INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_D3), throwCounter, nException );
189 #endif
190 
191  }
192  catch (std::logic_error err) {
193  *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
194  *outStream << err.what() << '\n';
195  *outStream << "-------------------------------------------------------------------------------" << "\n\n";
196  errorFlag = -1000;
197  };
198 
199  // Check if number of thrown exceptions matches the one we expect
200  if (throwCounter != nException) {
201  errorFlag++;
202  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
203  }
204 
205  *outStream \
206  << "\n"
207  << "===============================================================================\n"\
208  << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
209  << "===============================================================================\n";
210 
211  try{
212  std::vector<std::vector<int> > allTags = triBasis.getAllDofTags();
213 
214  // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
215  for (unsigned i = 0; i < allTags.size(); i++) {
216  int bfOrd = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
217 
218  std::vector<int> myTag = triBasis.getDofTag(bfOrd);
219  if( !( (myTag[0] == allTags[i][0]) &&
220  (myTag[1] == allTags[i][1]) &&
221  (myTag[2] == allTags[i][2]) &&
222  (myTag[3] == allTags[i][3]) ) ) {
223  errorFlag++;
224  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
225  *outStream << " getDofOrdinal( {"
226  << allTags[i][0] << ", "
227  << allTags[i][1] << ", "
228  << allTags[i][2] << ", "
229  << allTags[i][3] << "}) = " << bfOrd <<" but \n";
230  *outStream << " getDofTag(" << bfOrd << ") = { "
231  << myTag[0] << ", "
232  << myTag[1] << ", "
233  << myTag[2] << ", "
234  << myTag[3] << "}\n";
235  }
236  }
237 
238  // Now do the same but loop over basis functions
239  for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) {
240  std::vector<int> myTag = triBasis.getDofTag(bfOrd);
241  int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
242  if( bfOrd != myBfOrd) {
243  errorFlag++;
244  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
245  *outStream << " getDofTag(" << bfOrd << ") = { "
246  << myTag[0] << ", "
247  << myTag[1] << ", "
248  << myTag[2] << ", "
249  << myTag[3] << "} but getDofOrdinal({"
250  << myTag[0] << ", "
251  << myTag[1] << ", "
252  << myTag[2] << ", "
253  << myTag[3] << "} ) = " << myBfOrd << "\n";
254  }
255  }
256  }
257  catch (std::logic_error err){
258  *outStream << err.what() << "\n\n";
259  errorFlag = -1000;
260  };
261 
262  *outStream \
263  << "\n"
264  << "===============================================================================\n"\
265  << "| TEST 3: correctness of basis function values |\n"\
266  << "===============================================================================\n";
267 
268  outStream -> precision(20);
269 
270  // VALUE: Each row gives the 3 correct basis set values at an evaluation point
271  double basisValues[] = {
272  1.0, 0.0, 0.0,
273  0.0, 1.0, 0.0,
274  0.0, 0.0, 1.0,
275  0.0, 0.5, 0.5,
276  0.25,0.0, 0.75
277  };
278 
279  // GRAD and D1: each row gives the 6 correct values of the gradients of the 3 basis functions
280  double basisGrads[] = {
281  -1.0, -1.0, 1.0, 0.0, 0.0, 1.0,
282  -1.0, -1.0, 1.0, 0.0, 0.0, 1.0,
283  -1.0, -1.0, 1.0, 0.0, 0.0, 1.0,
284  -1.0, -1.0, 1.0, 0.0, 0.0, 1.0,
285  -1.0, -1.0, 1.0, 0.0, 0.0, 1.0,
286  };
287 
288  // CURL: each row gives the 6 correct values of the curls of the 3 basis functions
289  double basisCurls[] = {
290  -1.0, 1.0, 0.0, -1.0, 1.0, 0.0,
291  -1.0, 1.0, 0.0, -1.0, 1.0, 0.0,
292  -1.0, 1.0, 0.0, -1.0, 1.0, 0.0,
293  -1.0, 1.0, 0.0, -1.0, 1.0, 0.0,
294  -1.0, 1.0, 0.0, -1.0, 1.0, 0.0
295  };
296 
297  try{
298 
299  // Dimensions for the output arrays:
300  int numFields = triBasis.getCardinality();
301  int numPoints = triNodes.dimension(0);
302  int spaceDim = triBasis.getBaseCellTopology().getDimension();
303 
304  // Generic array for values, grads, curls, etc. that will be properly sized before each call
306 
307  // Check VALUE of basis functions: resize vals to rank-2 container:
308  vals.resize(numFields, numPoints);
309  triBasis.getValues(vals, triNodes, OPERATOR_VALUE);
310  for (int i = 0; i < numFields; i++) {
311  for (int j = 0; j < numPoints; j++) {
312  int l = i + j * numFields;
313  if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
314  errorFlag++;
315  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
316 
317  // Output the multi-index of the value where the error is:
318  *outStream << " At multi-index { ";
319  *outStream << i << " ";*outStream << j << " ";
320  *outStream << "} computed value: " << vals(i,j)
321  << " but reference value: " << basisValues[l] << "\n";
322  }
323  }
324  }
325 
326  // Check GRAD of basis function: resize vals to rank-3 container
327  vals.resize(numFields, numPoints, spaceDim);
328  triBasis.getValues(vals, triNodes, OPERATOR_GRAD);
329  for (int i = 0; i < numFields; i++) {
330  for (int j = 0; j < numPoints; j++) {
331  for (int k = 0; k < spaceDim; k++) {
332  int l = k + i * spaceDim + j * spaceDim * numFields;
333  if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
334  errorFlag++;
335  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
336 
337  // Output the multi-index of the value where the error is:
338  *outStream << " At multi-index { ";
339  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
340  *outStream << "} computed grad component: " << vals(i,j,k)
341  << " but reference grad component: " << basisGrads[l] << "\n";
342  }
343  }
344  }
345  }
346 
347  // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
348  triBasis.getValues(vals, triNodes, OPERATOR_D1);
349  for (int i = 0; i < numFields; i++) {
350  for (int j = 0; j < numPoints; j++) {
351  for (int k = 0; k < spaceDim; k++) {
352  int l = k + i * spaceDim + j * spaceDim * numFields;
353  if (std::abs(vals(i,j,k) - basisGrads[l]) > INTREPID_TOL) {
354  errorFlag++;
355  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
356 
357  // Output the multi-index of the value where the error is:
358  *outStream << " At multi-index { ";
359  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
360  *outStream << "} computed D1 component: " << vals(i,j,k)
361  << " but reference D1 component: " << basisGrads[l] << "\n";
362  }
363  }
364  }
365  }
366 
367 
368  // Check CURL of basis function: resize vals just for illustration!
369  vals.resize(numFields, numPoints, spaceDim);
370  triBasis.getValues(vals, triNodes, OPERATOR_CURL);
371  for (int i = 0; i < numFields; i++) {
372  for (int j = 0; j < numPoints; j++) {
373  for (int k = 0; k < spaceDim; k++) {
374  int l = k + i * spaceDim + j * spaceDim * numFields;
375  if (std::abs(vals(i,j,k) - basisCurls[l]) > INTREPID_TOL) {
376  errorFlag++;
377  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
378 
379  // Output the multi-index of the value where the error is:
380  *outStream << " At multi-index { ";
381  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
382  *outStream << "} computed curl component: " << vals(i,j,k)
383  << " but reference curl component: " << basisCurls[l] << "\n";
384  }
385  }
386  }
387  }
388 
389  // Check all higher derivatives - must be zero.
390  for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) {
391 
392  // The last dimension is the number of kth derivatives and needs to be resized for every Dk
393  int DkCardin = Intrepid::getDkCardinality(op, spaceDim);
394  vals.resize(numFields, numPoints, DkCardin);
395 
396  triBasis.getValues(vals, triNodes, op);
397  for (int i = 0; i < vals.size(); i++) {
398  if (std::abs(vals[i]) > INTREPID_TOL) {
399  errorFlag++;
400  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
401 
402  // Get the multi-index of the value where the error is and the operator order
403  std::vector<int> myIndex;
404  vals.getMultiIndex(myIndex,i);
405  int ord = Intrepid::getOperatorOrder(op);
406  *outStream << " At multi-index { ";
407  for(int j = 0; j < vals.rank(); j++) {
408  *outStream << myIndex[j] << " ";
409  }
410  *outStream << "} computed D"<< ord <<" component: " << vals[i]
411  << " but reference D" << ord << " component: 0 \n";
412  }
413  }
414  }
415  }
416 
417  // Catch unexpected errors
418  catch (std::logic_error err) {
419  *outStream << err.what() << "\n\n";
420  errorFlag = -1000;
421  };
422 
423  *outStream \
424  << "\n"
425  << "===============================================================================\n"\
426  << "| TEST 4: correctness of DoF locations |\n"\
427  << "===============================================================================\n";
428 
429  try{
430  Teuchos::RCP<Basis<double, FieldContainer<double> > > basis =
431  Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM<double, FieldContainer<double> >);
432  Teuchos::RCP<DofCoordsInterface<FieldContainer<double> > > coord_iface =
433  Teuchos::rcp_dynamic_cast<DofCoordsInterface<FieldContainer<double> > >(basis);
434 
436  FieldContainer<double> bvals(basis->getCardinality(), basis->getCardinality());
437 
438  // Check exceptions.
439 #ifdef HAVE_INTREPID_DEBUG
440  cvals.resize(1,2,3);
441  INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
442  cvals.resize(4,2);
443  INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
444  cvals.resize(4,3);
445  INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException );
446 #endif
447  cvals.resize(3,2);
448  INTREPID_TEST_COMMAND( coord_iface->getDofCoords(cvals), throwCounter, nException ); nException--;
449  // Check if number of thrown exceptions matches the one we expect
450  if (throwCounter != nException) {
451  errorFlag++;
452  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
453  }
454 
455  // Check mathematical correctness.
456  basis->getValues(bvals, cvals, OPERATOR_VALUE);
457  char buffer[120];
458  for (int i=0; i<bvals.dimension(0); i++) {
459  for (int j=0; j<bvals.dimension(1); j++) {
460  if ((i != j) && (std::abs(bvals(i,j) - 0.0) > INTREPID_TOL)) {
461  errorFlag++;
462  sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), bvals(i,j), 0.0);
463  *outStream << buffer;
464  }
465  else if ((i == j) && (std::abs(bvals(i,j) - 1.0) > INTREPID_TOL)) {
466  errorFlag++;
467  sprintf(buffer, "\nValue of basis function %d at (%6.4e, %6.4e) is %6.4e but should be %6.4e!\n", i, cvals(i,0), cvals(i,1), bvals(i,j), 1.0);
468  *outStream << buffer;
469  }
470  }
471  }
472 
473  }
474  catch (std::logic_error err){
475  *outStream << err.what() << "\n\n";
476  errorFlag = -1000;
477  };
478 
479  if (errorFlag != 0)
480  std::cout << "End Result: TEST FAILED\n";
481  else
482  std::cout << "End Result: TEST PASSED\n";
483 
484  // reset format state of std::cout
485  std::cout.copyfmt(oldFormatState);
486 
487  return errorFlag;
488 }
Intrepid::FieldContainer::getMultiIndex
void getMultiIndex(int &i0, const int valueEnum) const
Returns the multi-index of a value, based on its enumeration, as a list, for rank-1 containers.
Definition: Intrepid_FieldContainerDef.hpp:863
Intrepid::Basis::getBaseCellTopology
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
Definition: Intrepid_BasisDef.hpp:112
Intrepid::Basis::getCardinality
virtual int getCardinality() const
Returns cardinality of the basis.
Definition: Intrepid_BasisDef.hpp:100
Intrepid::Basis::getAllDofTags
virtual const std::vector< std::vector< int > > & getAllDofTags()
Retrieves all DoF tags.
Definition: Intrepid_BasisDef.hpp:89
Intrepid_FieldContainer.hpp
Header file for utility class to provide multidimensional containers.
Intrepid::FieldContainer::rank
int rank() const
Return rank of the FieldContainer = number of indices used to tag the multi-indexed value.
Definition: Intrepid_FieldContainerDef.hpp:594
Intrepid::FieldContainer::resize
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0.
Definition: Intrepid_FieldContainerDef.hpp:1137
Intrepid::getDkCardinality
int getDkCardinality(const EOperator operatorType, const int spaceDim)
Returns cardinality of Dk, i.e., the number of all derivatives of order k.
Definition: Intrepid_Utils.cpp:400
Intrepid::FieldContainer< double >
Intrepid::getOperatorOrder
int getOperatorOrder(const EOperator operatorType)
Returns order of an operator.
Definition: Intrepid_Utils.cpp:196
Intrepid::Basis_HGRAD_TRI_C1_FEM::getValues
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
Definition: Intrepid_HGRAD_TRI_C1_FEMDef.hpp:95
Intrepid::Basis_HGRAD_TRI_C1_FEM
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell.
Definition: Intrepid_HGRAD_TRI_C1_FEM.hpp:81
main
int main(int argc, char *argv[])
Performs a code-code comparison to FIAT for Nedelec bases on tets (values and curls)
Definition: test_01.cpp:65
Intrepid::EOperator
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
Definition: Intrepid_Types.hpp:206
Intrepid::Basis::getDofOrdinal
virtual int getDofOrdinal(const int subcDim, const int subcOrd, const int subcDofOrd)
DoF tag to ordinal lookup.
Definition: Intrepid_BasisDef.hpp:51
Intrepid::FieldContainer::size
int size() const
Returns size of the FieldContainer defined as the product of its dimensions.
Definition: Intrepid_FieldContainerDef.hpp:601
Intrepid_HGRAD_TRI_C1_FEM.hpp
Header file for the Intrepid::HGRAD_TRI_C1_FEM class.
Intrepid::Basis::getDofTag
virtual const std::vector< int > & getDofTag(const int dofOrd)
DoF ordinal to DoF tag lookup.
Definition: Intrepid_BasisDef.hpp:78
Intrepid::DofCoordsInterface
This is an interface class for bases whose degrees of freedom can be associated with spatial location...
Definition: Intrepid_Basis.hpp:353