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_HDIV_TRI_I1_FEM) |\n" \
93  << "| |\n" \
94  << "| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
95  << "| 3) Exception tests on input arguments and invalid operators |\n" \
96  << "| 2) Basis values for VALUE and DIV operators |\n" \
97  << "| |\n" \
98  << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
99  << "| Denis Ridzal (dridzal@sandia.gov), |\n" \
100  << "| Kara Peterson (kjpeter@sandia.gov). |\n" \
101  << "| |\n" \
102  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
103  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
104  << "| |\n" \
105  << "===============================================================================\n"\
106  << "| TEST 1: Basis creation, exception testing |\n"\
107  << "===============================================================================\n";
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 TRI and its 3 edge midpoints.
118  FieldContainer<double> triNodes(6, 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  // edge midpoints
123  triNodes(3,0) = 0.5; triNodes(3,1) = 0.0;
124  triNodes(4,0) = 0.5; triNodes(4,1) = 0.5;
125  triNodes(5,0) = 0.0; triNodes(5,1) = 0.5;
126 
127 
128  // Generic array for the output values; needs to be properly resized depending on the operator type
130 
131  try{
132  // exception #1: GRAD cannot be applied to HDIV functions
133  // resize vals to rank-3 container with dimensions (num. basis functions, num. points, arbitrary)
134  vals.resize(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() );
135  INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_GRAD), throwCounter, nException );
136 
137  // exception #2: CURL cannot be applied to HDIV functions
138  INTREPID_TEST_COMMAND( triBasis.getValues(vals, triNodes, OPERATOR_CURL), throwCounter, nException );
139 
140  // Exceptions 3-7: all bf tags/bf Ids below are wrong and should cause getDofOrdinal() and
141  // getDofTag() to access invalid array elements thereby causing bounds check exception
142  // exception #3
143  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(3,0,0), throwCounter, nException );
144  // exception #4
145  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(1,1,1), throwCounter, nException );
146  // exception #5
147  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(0,2,1), throwCounter, nException );
148  // exception #6
149  INTREPID_TEST_COMMAND( triBasis.getDofTag(3), throwCounter, nException );
150  // exception #7
151  INTREPID_TEST_COMMAND( triBasis.getDofTag(-1), throwCounter, nException );
152  // exception #8
153  INTREPID_TEST_COMMAND( triBasis.getDofOrdinal(2,0,0), throwCounter, nException );
154 
155 #ifdef HAVE_INTREPID_DEBUG
156  // Exceptions 9-16 test exception handling with incorrectly dimensioned input/output arrays
157  // exception #9: input points array must be of rank-2
158  FieldContainer<double> badPoints1(4, 5, 3);
159  INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
160 
161  // exception #10 dimension 1 in the input point array must equal space dimension of the cell
162  FieldContainer<double> badPoints2(4, triBasis.getBaseCellTopology().getDimension() + 1);
163  INTREPID_TEST_COMMAND( triBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
164 
165  // exception #11 output values must be of rank-3 for OPERATOR_VALUE
166  FieldContainer<double> badVals1(4, 3);
167  INTREPID_TEST_COMMAND( triBasis.getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException );
168 
169  // exception #12 output values must be of rank-2 for OPERATOR_DIV
170  FieldContainer<double> badVals2(4, 3, 1);
171  INTREPID_TEST_COMMAND( triBasis.getValues(badVals2, triNodes, OPERATOR_VALUE), throwCounter, nException );
172 
173  // exception #13 incorrect 0th dimension of output array for OPERATOR_VALUE (must equal number of basis functions)
174  FieldContainer<double> badVals3(triBasis.getCardinality() + 1, triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension());
175  INTREPID_TEST_COMMAND( triBasis.getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException );
176 
177  // exception #14 incorrect 0th dimension of output array for OPERATOR_DIV (must equal number of basis functions)
178  FieldContainer<double> badVals4(triBasis.getCardinality() + 1, triNodes.dimension(0));
179  INTREPID_TEST_COMMAND( triBasis.getValues(badVals4, triNodes, OPERATOR_DIV), throwCounter, nException );
180 
181  // exception #15 incorrect 1st dimension of output array (must equal number of points)
182  FieldContainer<double> badVals5(triBasis.getCardinality(), triNodes.dimension(0) + 1, triBasis.getBaseCellTopology().getDimension());
183  INTREPID_TEST_COMMAND( triBasis.getValues(badVals5, triNodes, OPERATOR_VALUE), throwCounter, nException );
184 
185  // exception #16 incorrect 1st dimension of output array (must equal number of points)
186  FieldContainer<double> badVals6(triBasis.getCardinality(), triNodes.dimension(0) + 1);
187  INTREPID_TEST_COMMAND( triBasis.getValues(badVals6, triNodes, OPERATOR_DIV), throwCounter, nException );
188 
189  // exception #17: incorrect 2nd dimension of output array (must equal the space dimension)
190  FieldContainer<double> badVals7(triBasis.getCardinality(), triNodes.dimension(0), triBasis.getBaseCellTopology().getDimension() + 1);
191  INTREPID_TEST_COMMAND( triBasis.getValues(badVals7, triNodes, OPERATOR_VALUE), 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  if (throwCounter != nException) {
204  errorFlag++;
205  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
206  }
207 
208  *outStream \
209  << "\n"
210  << "===============================================================================\n"\
211  << "| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
212  << "===============================================================================\n";
213 
214  try{
215  std::vector<std::vector<int> > allTags = triBasis.getAllDofTags();
216 
217  // Loop over all tags, lookup the associated dof enumeration and then lookup the tag again
218  for (unsigned i = 0; i < allTags.size(); i++) {
219  int bfOrd = triBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
220 
221  std::vector<int> myTag = triBasis.getDofTag(bfOrd);
222  if( !( (myTag[0] == allTags[i][0]) &&
223  (myTag[1] == allTags[i][1]) &&
224  (myTag[2] == allTags[i][2]) &&
225  (myTag[3] == allTags[i][3]) ) ) {
226  errorFlag++;
227  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
228  *outStream << " getDofOrdinal( {"
229  << allTags[i][0] << ", "
230  << allTags[i][1] << ", "
231  << allTags[i][2] << ", "
232  << allTags[i][3] << "}) = " << bfOrd <<" but \n";
233  *outStream << " getDofTag(" << bfOrd << ") = { "
234  << myTag[0] << ", "
235  << myTag[1] << ", "
236  << myTag[2] << ", "
237  << myTag[3] << "}\n";
238  }
239  }
240 
241  // Now do the same but loop over basis functions
242  for( int bfOrd = 0; bfOrd < triBasis.getCardinality(); bfOrd++) {
243  std::vector<int> myTag = triBasis.getDofTag(bfOrd);
244  int myBfOrd = triBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
245  if( bfOrd != myBfOrd) {
246  errorFlag++;
247  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
248  *outStream << " getDofTag(" << bfOrd << ") = { "
249  << myTag[0] << ", "
250  << myTag[1] << ", "
251  << myTag[2] << ", "
252  << myTag[3] << "} but getDofOrdinal({"
253  << myTag[0] << ", "
254  << myTag[1] << ", "
255  << myTag[2] << ", "
256  << myTag[3] << "} ) = " << myBfOrd << "\n";
257  }
258  }
259  }
260  catch (const std::logic_error & err){
261  *outStream << err.what() << "\n\n";
262  errorFlag = -1000;
263  };
264 
265  *outStream \
266  << "\n"
267  << "===============================================================================\n"\
268  << "| TEST 3: correctness of basis function values |\n"\
269  << "===============================================================================\n";
270 
271  outStream -> precision(20);
272 
273  // VALUE: Correct basis values in (F,P,D) format: each group of two rows gives basis function
274  // values at vertices followed by midpoints. This is the same array format as the output from getValues.
275  double basisValues[] = {
276  // basis function 0 at 3 vertices followed by 3 midpoints
277  0.0,-1.0, 1.0,-1.0, 0.0, 0.0,
278  0.5,-1.0, 0.5,-0.5, 0.0,-0.5,
279  // basis function 1 at 3 vertices followed by 3 midpoints
280  0.0, 0.0, 1.0, 0.0, 0.0, 1.0,
281  0.5, 0.0, 0.5, 0.5, 0.0, 0.5,
282  // basis function 2 at 3 vertices followed by 3 midpoints
283  -1.0, 0.0, 0.0, 0.0, -1.0, 1.0,
284  -0.5, 0.0, -0.5, 0.5, -1.0, 0.5
285  };
286 
287  // DIV: each row gives the 3 correct values of the divergence of the 3 basis functions
288  double basisDivs[] = {
289  // 3 vertices
290  2.0, 2.0, 2.0,
291  2.0, 2.0, 2.0,
292  2.0, 2.0, 2.0,
293  // 3 edge centers
294  2.0, 2.0, 2.0,
295  2.0, 2.0, 2.0,
296  2.0, 2.0, 2.0,
297  };
298 
299  try{
300 
301  // Dimensions for the output arrays:
302  int numFields = triBasis.getCardinality();
303  int numPoints = triNodes.dimension(0);
304  int spaceDim = triBasis.getBaseCellTopology().getDimension();
305 
306  // Generic array for values and divs that will be properly sized before each call
308 
309  // Check VALUE of basis functions: resize vals to rank-3 container:
310  vals.resize(numFields, numPoints, spaceDim);
311  triBasis.getValues(vals, triNodes, OPERATOR_VALUE);
312  for (int i = 0; i < numFields; i++) {
313  for (int j = 0; j < numPoints; j++) {
314  for (int k = 0; k < spaceDim; k++) {
315  // basisValues are in (F,P,D) format and the multiindex is (i,j,k), here's the offset:
316  int l = k + j * spaceDim + i * spaceDim * numPoints;
317 
318  if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
319  errorFlag++;
320  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
321 
322  // Output the multi-index of the value where the error is:
323  *outStream << " address = "<< l <<"\n";
324  *outStream << " At multi-index { ";
325  *outStream << i << " ";*outStream << j << " ";*outStream << k << " ";
326  *outStream << "} computed value: " << vals(i,j,k)
327  << " but reference value: " << basisValues[l] << "\n";
328  }
329  }
330  }
331  }
332 
333 
334  // Check DIV of basis function: resize vals to rank-2 container
335  vals.resize(numFields, numPoints);
336  triBasis.getValues(vals, triNodes, OPERATOR_DIV);
337  for (int i = 0; i < numFields; i++) {
338  for (int j = 0; j < numPoints; j++) {
339  int l = i + j * numFields;
340  if (std::abs(vals(i,j) - basisDivs[l]) > INTREPID_TOL) {
341  errorFlag++;
342  *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n";
343 
344  // Output the multi-index of the value where the error is:
345  *outStream << " At (Field,Point,Dim) multi-index { ";
346  *outStream << i << " ";*outStream << j << " ";
347  *outStream << "} computed divergence component: " << vals(i,j)
348  << " but reference divergence component: " << basisDivs[l] << "\n";
349  }
350  }
351  }
352 
353  }
354 
355  // Catch unexpected errors
356  catch (const std::logic_error & err) {
357  *outStream << err.what() << "\n\n";
358  errorFlag = -1000;
359  };
360 
361  if (errorFlag != 0)
362  std::cout << "End Result: TEST FAILED\n";
363  else
364  std::cout << "End Result: TEST PASSED\n";
365 
366  // reset format state of std::cout
367  std::cout.copyfmt(oldFormatState);
368 
369  return errorFlag;
370 }
virtual const std::vector< int > & getDofTag(const int dofOrd)
DoF ordinal to DoF tag lookup.
virtual int getCardinality() const
Returns cardinality of the basis.
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
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(div)-compatible FEM basis of degree 1 on a Triangle cell...
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
Header file for the Intrepid::HDIV_TRI_I1_FEM class.
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...