51 #include "Teuchos_oblackholestream.hpp"
52 #include "Teuchos_RCP.hpp"
53 #include "Teuchos_GlobalMPISession.hpp"
56 using namespace Intrepid;
58 #define INTREPID_TEST_COMMAND( S , throwCounter, nException ) \
64 catch (const std::logic_error & err) { \
66 *outStream << "Expected Error " << nException << " -------------------------------------------------------------\n"; \
67 *outStream << err.what() << '\n'; \
68 *outStream << "-------------------------------------------------------------------------------" << "\n\n"; \
72 int main(
int argc,
char *argv[]) {
74 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
78 int iprint = argc - 1;
79 Teuchos::RCP<std::ostream> outStream;
80 Teuchos::oblackholestream bhs;
82 outStream = Teuchos::rcp(&std::cout,
false);
84 outStream = Teuchos::rcp(&bhs,
false);
87 Teuchos::oblackholestream oldFormatState;
88 oldFormatState.copyfmt(std::cout);
91 <<
"===============================================================================\n" \
93 <<
"| Unit Test (Basis_HCURL_QUAD_In_FEM) |\n" \
95 <<
"| 1) Conversion of Dof tags into Dof ordinals and back |\n" \
96 <<
"| 2) Basis values for VALUE and CURL 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";
111 shards::CellTopology line(shards::getCellTopologyData< shards::Line<> >());
118 int throwCounter = 0;
122 quadNodes(0,0) = -1.0; quadNodes(0,1) = -1.0;
123 quadNodes(1,0) = 0.0; quadNodes(1,1) = -1.0;
124 quadNodes(2,0) = 1.0; quadNodes(2,1) = -1.0;
125 quadNodes(3,0) = -1.0; quadNodes(3,1) = 0.0;
126 quadNodes(4,0) = 0.0; quadNodes(4,1) = 0.0;
127 quadNodes(5,0) = 1.0; quadNodes(5,1) = 0.0;
128 quadNodes(6,0) = -1.0; quadNodes(6,1) = 1.0;
129 quadNodes(7,0) = 0.0; quadNodes(7,1) = 1.0;
130 quadNodes(8,0) = 1.0; quadNodes(8,1) = 1.0;
139 vals.
resize(quadBasis.getCardinality(), quadNodes.dimension(0), 4 );
140 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_GRAD), throwCounter, nException );
144 vals.
resize(quadBasis.getCardinality(), quadNodes.dimension(0) );
145 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_DIV), throwCounter, nException );
150 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(3,0,0), throwCounter, nException );
152 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(1,1,1), throwCounter, nException );
154 INTREPID_TEST_COMMAND( quadBasis.getDofOrdinal(0,4,1), throwCounter, nException );
156 INTREPID_TEST_COMMAND( quadBasis.getDofTag(12), throwCounter, nException );
158 INTREPID_TEST_COMMAND( quadBasis.getDofTag(-1), throwCounter, nException );
160 #ifdef HAVE_INTREPID_DEBUG
164 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
168 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
172 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals1, quadNodes, OPERATOR_VALUE), throwCounter, nException );
176 INTREPID_TEST_COMMAND( quadBasis.getValues(badCurls1, quadNodes, OPERATOR_CURL), throwCounter, nException );
179 FieldContainer<double> badVals2(quadBasis.getCardinality() + 1, quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension());
180 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals2, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
183 FieldContainer<double> badVals3(quadBasis.getCardinality(), quadNodes.dimension(0) + 1, quadBasis.getBaseCellTopology().getDimension() );
184 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals3, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
187 FieldContainer<double> badVals4(quadBasis.getCardinality(), quadNodes.dimension(0), quadBasis.getBaseCellTopology().getDimension() - 1);
188 INTREPID_TEST_COMMAND( quadBasis.getValues(badVals4, quadNodes, OPERATOR_VALUE), throwCounter, nException ) ;
192 vals.
resize(quadBasis.getCardinality(),
193 quadNodes.dimension(0),
194 Intrepid::getDkCardinality(OPERATOR_D2, quadBasis.getBaseCellTopology().getDimension()));
195 INTREPID_TEST_COMMAND( quadBasis.getValues(vals, quadNodes, OPERATOR_D2), throwCounter, nException );
199 catch (
const std::logic_error & err) {
200 *outStream <<
"UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
201 *outStream << err.what() <<
'\n';
202 *outStream <<
"-------------------------------------------------------------------------------" <<
"\n\n";
208 if (throwCounter != nException) {
210 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
216 <<
"===============================================================================\n"\
217 <<
"| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
218 <<
"===============================================================================\n";
221 std::vector<std::vector<int> > allTags = quadBasis.getAllDofTags();
224 for (
unsigned i = 0; i < allTags.size(); i++) {
225 int bfOrd = quadBasis.getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
227 std::vector<int> myTag = quadBasis.getDofTag(bfOrd);
228 if( !( (myTag[0] == allTags[i][0]) &&
229 (myTag[1] == allTags[i][1]) &&
230 (myTag[2] == allTags[i][2]) &&
231 (myTag[3] == allTags[i][3]) ) ) {
233 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
234 *outStream <<
" getDofOrdinal( {"
235 << allTags[i][0] <<
", "
236 << allTags[i][1] <<
", "
237 << allTags[i][2] <<
", "
238 << allTags[i][3] <<
"}) = " << bfOrd <<
" but \n";
239 *outStream <<
" getDofTag(" << bfOrd <<
") = { "
243 << myTag[3] <<
"}\n";
248 for(
int bfOrd = 0; bfOrd < quadBasis.getCardinality(); bfOrd++) {
249 std::vector<int> myTag = quadBasis.getDofTag(bfOrd);
250 int myBfOrd = quadBasis.getDofOrdinal(myTag[0], myTag[1], myTag[2]);
251 if( bfOrd != myBfOrd) {
253 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
254 *outStream <<
" getDofTag(" << bfOrd <<
") = { "
258 << myTag[3] <<
"} but getDofOrdinal({"
262 << myTag[3] <<
"} ) = " << myBfOrd <<
"\n";
266 catch (
const std::logic_error & err){
267 *outStream << err.what() <<
"\n\n";
273 <<
"===============================================================================\n"\
274 <<
"| TEST 3: correctness of basis function values |\n"\
275 <<
"===============================================================================\n";
277 outStream -> precision(20);
280 double basisValues[] = {
308 double basisCurls[] = {
309 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5,
310 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
311 -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5, -0.5,
312 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5
318 int numFields = quadBasis.getCardinality();
319 int numPoints = quadNodes.dimension(0);
320 int spaceDim = quadBasis.getBaseCellTopology().getDimension();
326 vals.
resize(numFields, numPoints, spaceDim);
327 quadBasis.getValues(vals, quadNodes, OPERATOR_VALUE);
328 for (
int i = 0; i < numFields; i++) {
329 for (
int j = 0; j < numPoints; j++) {
330 for (
int k = 0; k < spaceDim; k++) {
333 int l = k + i * spaceDim * numPoints + j * spaceDim;
334 if (std::abs(vals(i,j,k) - basisValues[l]) > INTREPID_TOL) {
336 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
339 *outStream <<
" At multi-index { ";
340 *outStream << i <<
" ";*outStream << j <<
" ";*outStream << k <<
" ";
341 *outStream <<
"} computed value: " << vals(i,j,k)
342 <<
" but reference value: " << basisValues[l] <<
"\n";
349 vals.
resize(numFields, numPoints);
350 quadBasis.getValues(vals, quadNodes, OPERATOR_CURL);
351 for (
int i = 0; i < numFields; i++) {
352 for (
int j = 0; j < numPoints; j++) {
353 int l = j + i * numPoints;
354 if (std::abs(vals(i,j) - basisCurls[l]) > INTREPID_TOL) {
356 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
359 *outStream <<
" At multi-index { ";
360 *outStream << i <<
" ";*outStream << j <<
" ";
361 *outStream <<
"} computed curl component: " << vals(i,j)
362 <<
" but reference curl component: " << basisCurls[l] <<
"\n";
370 catch (
const std::logic_error & err) {
371 *outStream << err.what() <<
"\n\n";
376 std::cout <<
"End Result: TEST FAILED\n";
378 std::cout <<
"End Result: TEST PASSED\n";
381 std::cout.copyfmt(oldFormatState);
Header file for utility class to provide multidimensional containers.
void resize(const int dim0)
Resizes FieldContainer to a rank-1 container with the specified dimension, initialized by 0...
Header file for the Intrepid::HCURL_QUAD_In_FEM class.
Implementation of the default H(div)-compatible FEM basis of degree 1 on Quadrilateral cell...