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 
72 int main(int argc, char *argv[]) {
73  Teuchos::oblackholestream oldFormatState;
74  oldFormatState.copyfmt(std::cout);
75 
76 #ifndef KOKKOS_ENABLE_CUDA
77  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
78 
79  // This little trick lets us print to std::cout only if
80  // a (dummy) command-line argument is provided.
81  int iprint = argc - 1;
82  Teuchos::RCP<std::ostream> outStream;
83  Teuchos::oblackholestream bhs; // outputs nothing
84  if (iprint > 0)
85  outStream = Teuchos::rcp(&std::cout, false);
86  else
87  outStream = Teuchos::rcp(&bhs, false);
88 
89  // Save the format state of the original std::cout.
90 
91  *outStream \
92  << "===============================================================================\n" \
93  << "| |\n" \
94  << "| Unit Test (Basis_HGRAD_LINE_C1_FEM) |\n" \
95  << "| |\n" \
96  << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
97  << "| 2) Basis values for VALUE, GRAD, CURL, DIV, and Dk operators |\n" \
98  << "| |\n" \
99  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
100  << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
101  << "| Kara Peterson (dridzal@sandia.gov). |\n" \
102  << "| |\n" \
103  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
104  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
105  << "| |\n" \
106  << "===============================================================================\n"\
107  << "| TEST 1: Basis creation, exception testing |\n"\
108  << "===============================================================================\n";
109 #endif
110 
111  // Define basis and error flag
113  int errorFlag = 0;
114 #ifndef KOKKOS_ENABLE_CUDA
115  // Initialize throw counter for exception testing
116  int nException = 0;
117  int throwCounter = 0;
118 
119  // Define array containing the 2 vertices of the reference Line, its center and another point
120  FieldContainer<double> lineNodes(4, 1);
121  lineNodes(0,0) = -1.0;
122  lineNodes(1,0) = 1.0;
123  lineNodes(2,0) = 0.0;
124  lineNodes(3,0) = 0.5;
125 
126  // Generic array for the output values; needs to be properly resized depending on the operator type
128 
129  try{
130  // resize vals to rank-2 container with dimensions (num. points, num. basis functions)
131  vals.resize(lineBasis.getCardinality(), lineNodes.dimension(0) );
132 
133  // Exceptions 1-5: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
134  // getDofTag() to access invalid array elements thereby causing bounds check exception
135  // exception #1
136  INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(2,0,0), throwCounter, nException );
137  // exception #2
138  INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(1,1,1), throwCounter, nException );
139  // exception #3
140  INTREPID_TEST_COMMAND( lineBasis.getDofOrdinal(0,4,0), throwCounter, nException );
141  // exception #4
142  INTREPID_TEST_COMMAND( lineBasis.getDofTag(5), throwCounter, nException );
143  // exception #5
144  INTREPID_TEST_COMMAND( lineBasis.getDofTag(-1), throwCounter, nException );
145 
146 #ifdef HAVE_INTREPID_DEBUG
147  // Exceptions 6-17 test exception handling with incorrectly dimensioned input/output arrays
148  // exception #6: input points array must be of rank-2
149  FieldContainer<double> badPoints1(4, 5, 3);
150  INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
151 
152  // exception #7 dimension 1 in the input point array must equal space dimension of the cell
153  FieldContainer<double> badPoints2(4, 2);
154  INTREPID_TEST_COMMAND( lineBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
155 
156  // exception #8 output values must be of rank-2 for OPERATOR_VALUE
157  FieldContainer<double> badVals1(4, 3, 1);
158  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals1, lineNodes, OPERATOR_VALUE), throwCounter, nException );
159 
160  // exception #9 output values must be of rank-3 for OPERATOR_GRAD
161  FieldContainer<double> badVals2(4, 3);
162  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_GRAD), throwCounter, nException );
163 
164  // exception #10 output values must be of rank-3 for OPERATOR_DIV
165  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_DIV), throwCounter, nException );
166 
167  // exception #11 output values must be of rank-3 for OPERATOR_CURL
168  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_CURL), throwCounter, nException );
169 
170  // exception #12 output values must be of rank-3 for OPERATOR_D2
171  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals2, lineNodes, OPERATOR_D2), throwCounter, nException );
172 
173  // exception #13 incorrect 1st dimension of output array (must equal number of basis functions)
174  FieldContainer<double> badVals3(lineBasis.getCardinality() + 1, lineNodes.dimension(0));
175  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals3, lineNodes, OPERATOR_VALUE), throwCounter, nException );
176 
177  // exception #14 incorrect 0th dimension of output array (must equal number of points)
178  FieldContainer<double> badVals4(lineBasis.getCardinality(), lineNodes.dimension(0) + 1);
179  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals4, lineNodes, OPERATOR_VALUE), throwCounter, nException );
180 
181  // exception #15: incorrect 2nd dimension of output array (must equal the space dimension)
182  FieldContainer<double> badVals5(lineBasis.getCardinality(), lineNodes.dimension(0), 2);
183  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals5, lineNodes, OPERATOR_GRAD), throwCounter, nException );
184 
185  // exception #16: incorrect 2nd dimension of output array (must equal D2 cardinality in 1D)
186  FieldContainer<double> badVals6(lineBasis.getCardinality(), lineNodes.dimension(0), 40);
187  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals6, lineNodes, OPERATOR_D2), throwCounter, nException );
188 
189  // exception #17: incorrect 2nd dimension of output array (must equal D3 cardinality in 1D)
190  INTREPID_TEST_COMMAND( lineBasis.getValues(badVals6, lineNodes, OPERATOR_D3), throwCounter, nException );
191 #endif
192 
193  }
194  catch (std::logic_error err) {
195  *outStream << "UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
196  *outStream << err.what() << '\n';
197  *outStream << "-------------------------------------------------------------------------------" << "\n\n";
198  errorFlag = -1000;
199  };
200 
201  // Check if number of thrown exceptions matches the one we expect
202  if (throwCounter != nException) {
203  errorFlag++;
204  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
205  }
206 
207  *outStream \
208  << "\n"
209  << "===============================================================================\n"\
210  << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
211  << "===============================================================================\n";
212 
213  try{
214  std::vector<std::vector<int> > allTags = lineBasis.getAllDofTags();
215 
216  // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
217  for (unsigned i = 0; i < allTags.size(); i++) {
218  int bfOrd = lineBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
219 
220  std::vector<int> myTag = lineBasis.getDofTag(bfOrd);
221  if( !( (myTag[0] == allTags[i][0]) &&
222  (myTag[1] == allTags[i][1]) &&
223  (myTag[2] == allTags[i][2]) &&
224  (myTag[3] == allTags[i][3]) ) ) {
225  errorFlag++;
226  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
227  *outStream << " getDofOrdinal( {"
228  << allTags[i][0] << ", "
229  << allTags[i][1] << ", "
230  << allTags[i][2] << ", "
231  << allTags[i][3] << "}) = " << bfOrd <<" but \n";
232  *outStream << " getDofTag(" << bfOrd << ") = { "
233  << myTag[0] << ", "
234  << myTag[1] << ", "
235  << myTag[2] << ", "
236  << myTag[3] << "}\n";
237  }
238  }
239 
240  // Now do the same but loop over basis functions
241  for( int bfOrd = 0; bfOrd < lineBasis.getCardinality(); bfOrd++) {
242  std::vector<int> myTag = lineBasis.getDofTag(bfOrd);
243  int myBfOrd = lineBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
244  if( bfOrd != myBfOrd) {
245  errorFlag++;
246  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
247  *outStream << " getDofTag(" << bfOrd << ") = { "
248  << myTag[0] << ", "
249  << myTag[1] << ", "
250  << myTag[2] << ", "
251  << myTag[3] << "} but getDofOrdinal({"
252  << myTag[0] << ", "
253  << myTag[1] << ", "
254  << myTag[2] << ", "
255  << myTag[3] << "} ) = " << myBfOrd << "\n";
256  }
257  }
258  }
259  catch (std::logic_error err){
260  *outStream << err.what() << "\n\n";
261  errorFlag = -1000;
262  };
263 
264  *outStream \
265  << "\n"
266  << "===============================================================================\n"\
267  << "| TEST 3: correctness of basis function values |\n"\
268  << "===============================================================================\n";
269 
270  outStream -> precision(20);
271 
272  // VALUE: Each row gives the 2 correct basis set values at an evaluation point
273  double basisValues[] = {
274  1.0, 0.0,
275  0.0, 1.0,
276  0.5, 0.5,
277  0.25, 0.75
278  };
279 
280  // GRAD, DIV, CURL and D1: each row gives the 2 correct values of the gradients of the 2 basis functions
281  double basisDerivs[] = {
282  -0.5, 0.5,
283  -0.5, 0.5,
284  -0.5, 0.5,
285  -0.5, 0.5
286  };
287 
288  try{
289 
290  // Dimensions for the output arrays:
291  int numFields = lineBasis.getCardinality();
292  int numPoints = lineNodes.dimension(0);
293  int spaceDim = lineBasis.getBaseCellTopology().getDimension();
294 
295  // Generic array for values, grads, curls, etc. that will be properly sized before each call
297 
298  // Check VALUE of basis functions: resize vals to rank-2 container:
299  vals.resize(numFields, numPoints);
300  lineBasis.getValues(vals, lineNodes, OPERATOR_VALUE);
301  for (int i = 0; i < numFields; i++) {
302  for (int j = 0; j < numPoints; j++) {
303  int l = i + j * numFields;
304  if (std::abs(vals(i,j) - basisValues[l]) > INTREPID_TOL) {
305  errorFlag++;
306  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
307 
308  // Output the multi-index of the value where the error is:
309  *outStream << " At multi-index { ";
310  *outStream << i << " ";*outStream << j << " ";
311  *outStream << "} computed value: " << vals(i,j)
312  << " but reference value: " << basisValues[l] << "\n";
313  }
314  }
315  }
316 
317  // Check derivatives of basis function: resize vals to rank-3 container
318  vals.resize(numFields, numPoints, spaceDim);
319  lineBasis.getValues(vals, lineNodes, OPERATOR_GRAD);
320  for (int i = 0; i < numFields; i++) {
321  for (int j = 0; j < numPoints; j++) {
322  for (int k = 0; k < spaceDim; k++) {
323  int l = k + i * spaceDim + j * spaceDim * numFields;
324  if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
325  errorFlag++;
326  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
327 
328  // Output the multi-index of the value where the error is:
329  *outStream << " At multi-index { ";
330  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
331  *outStream << "} computed grad component: " << vals(i,j,k)
332  << " but reference grad component: " << basisDerivs[l] << "\n";
333  }
334  }
335  }
336  }
337 
338  // Check D1 of basis function (do not resize vals because it has the correct size: D1 = GRAD)
339  lineBasis.getValues(vals, lineNodes, OPERATOR_D1);
340  for (int i = 0; i < numFields; i++) {
341  for (int j = 0; j < numPoints; j++) {
342  for (int k = 0; k < spaceDim; k++) {
343  int l = k + i * spaceDim + j * spaceDim * numFields;
344  if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
345  errorFlag++;
346  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
347 
348  // Output the multi-index of the value where the error is:
349  *outStream << " At multi-index { ";
350  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
351  *outStream << "} computed D1 component: " << vals(i,j,k)
352  << " but reference D1 component: " << basisDerivs[l] << "\n";
353  }
354  }
355  }
356  }
357 
358 
359  // Check CURL of basis function: resize vals just for illustration!
360  vals.resize(numFields, numPoints, spaceDim);
361  lineBasis.getValues(vals, lineNodes, OPERATOR_CURL);
362  for (int i = 0; i < numFields; i++) {
363  for (int j = 0; j < numPoints; j++) {
364  for (int k = 0; k < spaceDim; k++) {
365  int l = k + i * spaceDim + j * spaceDim * numFields;
366  if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
367  errorFlag++;
368  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
369 
370  // Output the multi-index of the value where the error is:
371  *outStream << " At multi-index { ";
372  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
373  *outStream << "} computed curl component: " << vals(i,j,k)
374  << " but reference curl component: " << basisDerivs[l] << "\n";
375  }
376  }
377  }
378  }
379 
380 
381  // Check DIV of basis function (do not resize vals because it has the correct size: DIV = CURL)
382  lineBasis.getValues(vals, lineNodes, OPERATOR_DIV);
383  for (int i = 0; i < numFields; i++) {
384  for (int j = 0; j < numPoints; j++) {
385  for (int k = 0; k < spaceDim; k++) {
386  int l = k + i * spaceDim + j * spaceDim * numFields;
387  if (std::abs(vals(i,j,k) - basisDerivs[l]) > INTREPID_TOL) {
388  errorFlag++;
389  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
390 
391  // Output the multi-index of the value where the error is:
392  *outStream << " At multi-index { ";
393  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
394  *outStream << "} computed D1 component: " << vals(i,j,k)
395  << " but reference D1 component: " << basisDerivs[l] << "\n";
396  }
397  }
398  }
399  }
400 
401 
402  // Check all higher derivatives - must be zero.
403  for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) {
404 
405  // The last dimension is the number of kth derivatives and needs to be resized for every Dk
406  int DkCardin = Intrepid::getDkCardinality(op, spaceDim);
407  vals.resize(numFields, numPoints, DkCardin);
408 
409  lineBasis.getValues(vals, lineNodes, op);
410  for (int i = 0; i < vals.size(); i++) {
411  if (std::abs(vals[i]) > INTREPID_TOL) {
412  errorFlag++;
413  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
414 
415  // Get the multi-index of the value where the error is and the operator order
416  std::vector<int> myIndex;
417  vals.getMultiIndex(myIndex,i);
418  int ord = Intrepid::getOperatorOrder(op);
419  *outStream << " At multi-index { ";
420  for(int j = 0; j < vals.rank(); j++) {
421  *outStream << myIndex[j] << " ";
422  }
423  *outStream << "} computed D"<< ord <<" component: " << vals[i]
424  << " but reference D" << ord << " component: 0 \n";
425  }
426  }
427  }
428  }
429 
430  // Catch unexpected errors
431  catch (std::logic_error err) {
432  *outStream << err.what() << "\n\n";
433  errorFlag = -1000;
434  };
435 
436 #endif
437  if (errorFlag != 0)
438  std::cout << "End Result: TEST FAILED\n";
439  else
440  std::cout << "End Result: TEST PASSED\n";
441 
442  // reset format state of std::cout
443  std::cout.copyfmt(oldFormatState);
444  return errorFlag;
445 }
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
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HGRAD_LINE_C1_FEM class.
int getOperatorOrder(const EOperator operatorType)
Returns order of an operator.
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 default H(grad)-compatible FEM basis of degree 1 on Line cell.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Line 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.
int size() const
Returns size of the FieldContainer defined as the product of its dimensions.
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0.
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.
int rank() const
Return rank of the FieldContainer = number of indices used to tag the multi-indexed value.