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 (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[]) {
73 Teuchos::oblackholestream oldFormatState;
74 oldFormatState.copyfmt(std::cout);
76 #ifndef KOKKOS_ENABLE_CUDA
77 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
81 int iprint = argc - 1;
82 Teuchos::RCP<std::ostream> outStream;
83 Teuchos::oblackholestream bhs;
85 outStream = Teuchos::rcp(&std::cout,
false);
87 outStream = Teuchos::rcp(&bhs,
false);
92 <<
"===============================================================================\n" \
94 <<
"| Unit Test (Basis_HGRAD_LINE_C1_FEM) |\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" \
99 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \
100 <<
"| Denis Ridzal (dridzal@sandia.gov), |\n" \
101 <<
"| Kara Peterson (dridzal@sandia.gov). |\n" \
103 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
104 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
106 <<
"===============================================================================\n"\
107 <<
"| TEST 1: Basis creation, exception testing |\n"\
108 <<
"===============================================================================\n";
114 #ifndef KOKKOS_ENABLE_CUDA
117 int throwCounter = 0;
121 lineNodes(0,0) = -1.0;
122 lineNodes(1,0) = 1.0;
123 lineNodes(2,0) = 0.0;
124 lineNodes(3,0) = 0.5;
136 INTREPID_TEST_COMMAND( lineBasis.
getDofOrdinal(2,0,0), throwCounter, nException );
138 INTREPID_TEST_COMMAND( lineBasis.
getDofOrdinal(1,1,1), throwCounter, nException );
140 INTREPID_TEST_COMMAND( lineBasis.
getDofOrdinal(0,4,0), throwCounter, nException );
142 INTREPID_TEST_COMMAND( lineBasis.
getDofTag(5), throwCounter, nException );
144 INTREPID_TEST_COMMAND( lineBasis.
getDofTag(-1), throwCounter, nException );
146 #ifdef HAVE_INTREPID_DEBUG
150 INTREPID_TEST_COMMAND( lineBasis.
getValues(vals, badPoints1, OPERATOR_VALUE), throwCounter, nException );
154 INTREPID_TEST_COMMAND( lineBasis.
getValues(vals, badPoints2, OPERATOR_VALUE), throwCounter, nException );
158 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals1, lineNodes, OPERATOR_VALUE), throwCounter, nException );
162 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals2, lineNodes, OPERATOR_GRAD), throwCounter, nException );
165 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals2, lineNodes, OPERATOR_DIV), throwCounter, nException );
168 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals2, lineNodes, OPERATOR_CURL), throwCounter, nException );
171 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals2, lineNodes, OPERATOR_D2), throwCounter, nException );
175 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals3, lineNodes, OPERATOR_VALUE), throwCounter, nException );
179 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals4, lineNodes, OPERATOR_VALUE), throwCounter, nException );
183 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals5, lineNodes, OPERATOR_GRAD), throwCounter, nException );
187 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals6, lineNodes, OPERATOR_D2), throwCounter, nException );
190 INTREPID_TEST_COMMAND( lineBasis.
getValues(badVals6, lineNodes, OPERATOR_D3), throwCounter, nException );
194 catch (std::logic_error err) {
195 *outStream <<
"UNEXPECTED ERROR !!! ----------------------------------------------------------\n";
196 *outStream << err.what() <<
'\n';
197 *outStream <<
"-------------------------------------------------------------------------------" <<
"\n\n";
202 if (throwCounter != nException) {
204 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
209 <<
"===============================================================================\n"\
210 <<
"| TEST 2: correctness of tag to enum and enum to tag lookups |\n"\
211 <<
"===============================================================================\n";
214 std::vector<std::vector<int> > allTags = lineBasis.
getAllDofTags();
217 for (
unsigned i = 0; i < allTags.size(); i++) {
218 int bfOrd = lineBasis.
getDofOrdinal(allTags[i][0], allTags[i][1], allTags[i][2]);
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]) ) ) {
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 <<
") = { "
236 << myTag[3] <<
"}\n";
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) {
246 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
247 *outStream <<
" getDofTag(" << bfOrd <<
") = { "
251 << myTag[3] <<
"} but getDofOrdinal({"
255 << myTag[3] <<
"} ) = " << myBfOrd <<
"\n";
259 catch (std::logic_error err){
260 *outStream << err.what() <<
"\n\n";
266 <<
"===============================================================================\n"\
267 <<
"| TEST 3: correctness of basis function values |\n"\
268 <<
"===============================================================================\n";
270 outStream -> precision(20);
273 double basisValues[] = {
281 double basisDerivs[] = {
292 int numPoints = lineNodes.dimension(0);
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) {
306 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
309 *outStream <<
" At multi-index { ";
310 *outStream << i <<
" ";*outStream << j <<
" ";
311 *outStream <<
"} computed value: " << vals(i,j)
312 <<
" but reference value: " << basisValues[l] <<
"\n";
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) {
326 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
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) {
346 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
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) {
368 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
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) {
389 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
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";
403 for(EOperator op = OPERATOR_D2; op < OPERATOR_MAX; op++) {
407 vals.
resize(numFields, numPoints, DkCardin);
409 lineBasis.
getValues(vals, lineNodes, op);
410 for (
int i = 0; i < vals.
size(); i++) {
411 if (std::abs(vals[i]) > INTREPID_TOL) {
413 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
416 std::vector<int> myIndex;
419 *outStream <<
" At multi-index { ";
420 for(
int j = 0; j < vals.
rank(); j++) {
421 *outStream << myIndex[j] <<
" ";
423 *outStream <<
"} computed D"<< ord <<
" component: " << vals[i]
424 <<
" but reference D" << ord <<
" component: 0 \n";
431 catch (std::logic_error err) {
432 *outStream << err.what() <<
"\n\n";
438 std::cout <<
"End Result: TEST FAILED\n";
440 std::cout <<
"End Result: TEST PASSED\n";
443 std::cout.copyfmt(oldFormatState);
int main(int argc, char *argv[])
Performs a code-code comparison to FIAT for Nedelec bases on tets (values and curls)
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.