50 #include "Teuchos_oblackholestream.hpp"
51 #include "Teuchos_RCP.hpp"
52 #include "Teuchos_GlobalMPISession.hpp"
55 using namespace Intrepid;
57 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
63 catch (const std::logic_error & err) { \
65 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
66 *outStream << err.what() << '\n'; \
67 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
71 int main(
int argc,
char *argv[]) {
73 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
77 int iprint = argc - 1;
78 Teuchos::RCP<std::ostream> outStream;
79 Teuchos::oblackholestream bhs;
81 outStream = Teuchos::rcp(&std::cout,
false);
83 outStream = Teuchos::rcp(&bhs,
false);
86 Teuchos::oblackholestream oldFormatState;
87 oldFormatState.copyfmt(std::cout);
90 <<
"===============================================================================\n" \
92 <<
"| Unit Test (Basis_HDIV_TRI_I1_FEM) |\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" \
98 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
99 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
100 <<
"| Kara Peterson (kjpeter@sandia.gov). |\n" \
102 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
103 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
105 <<
"===============================================================================\n"\
106 <<
"| TEST 1: Basis creation, exception testing |\n"\
107 <<
"===============================================================================\n";
115 int throwCounter = 0;
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;
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;
135 INTREPID_TEST_COMMAND( triBasis.
getValues(vals, triNodes, OPERATOR_GRAD), throwCounter, nException );
138 INTREPID_TEST_COMMAND( triBasis.
getValues(vals, triNodes, OPERATOR_CURL), throwCounter, nException );
143 INTREPID_TEST_COMMAND( triBasis.
getDofOrdinal(3,0,0), throwCounter, nException );
145 INTREPID_TEST_COMMAND( triBasis.
getDofOrdinal(1,1,1), throwCounter, nException );
147 INTREPID_TEST_COMMAND( triBasis.
getDofOrdinal(0,2,1), throwCounter, nException );
149 INTREPID_TEST_COMMAND( triBasis.
getDofTag(3), throwCounter, nException );
151 INTREPID_TEST_COMMAND( triBasis.
getDofTag(-1), throwCounter, nException );
153 INTREPID_TEST_COMMAND( triBasis.
getDofOrdinal(2,0,0), throwCounter, nException );
155 #ifdef HAVE_INTREPID_DEBUG
159 INTREPID_TEST_COMMAND( triBasis.
getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
163 INTREPID_TEST_COMMAND( triBasis.
getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
167 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals1, triNodes, OPERATOR_VALUE), throwCounter, nException );
171 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals2, triNodes, OPERATOR_VALUE), throwCounter, nException );
175 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals3, triNodes, OPERATOR_VALUE), throwCounter, nException );
179 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals4, triNodes, OPERATOR_DIV), throwCounter, nException );
183 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals5, triNodes, OPERATOR_VALUE), throwCounter, nException );
187 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals6, triNodes, OPERATOR_DIV), throwCounter, nException );
191 INTREPID_TEST_COMMAND( triBasis.
getValues(badVals7, triNodes, OPERATOR_VALUE), throwCounter, nException );
195 catch (
const std::logic_error & err) {
196 *outStream <<
"UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
197 *outStream << err.what() <<
'\n';
198 *outStream <<
"-------------------------------------------------------------------------------" <<
"\n\n";
203 if (throwCounter != nException) {
205 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
210 <<
"===============================================================================\n"\
211 <<
"| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
212 <<
"===============================================================================\n";
215 std::vector<std::vector<int> > allTags = triBasis.
getAllDofTags();
218 for (
unsigned i = 0; i < allTags.size(); i++) {
219 int bfOrd = triBasis.
getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
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]) ) ) {
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 <<
") = { "
237 << myTag[3] <<
"}\n";
243 std::vector<int> myTag = triBasis.
getDofTag(bfOrd);
244 int myBfOrd = triBasis.
getDofOrdinal(myTag[0], myTag[1], myTag[2]);
245 if( bfOrd != myBfOrd) {
247 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
248 *outStream <<
" getDofTag(" << bfOrd <<
") = { "
252 << myTag[3] <<
"} but getDofOrdinal({"
256 << myTag[3] <<
"} ) = " << myBfOrd <<
"\n";
260 catch (
const std::logic_error & err){
261 *outStream << err.what() <<
"\n\n";
267 <<
"===============================================================================\n"\
268 <<
"| TEST 3: correctness of basis function values |\n"\
269 <<
"===============================================================================\n";
271 outStream -> precision(20);
275 double basisValues[] = {
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,
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,
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
288 double basisDivs[] = {
303 int numPoints = triNodes.dimension(0);
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++) {
316 int l = k + j * spaceDim + i * spaceDim * numPoints;
318 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
320 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
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) {
342 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
356 catch (
const std::logic_error & err) {
357 *outStream << err.what() <<
"\n\n";
362 std::cout <<
"End Result: TEST FAILED\n";
364 std::cout <<
"End Result: TEST PASSED\n";
367 std::cout.copyfmt(oldFormatState);
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...