51#include "Teuchos_oblackholestream.hpp"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_GlobalMPISession.hpp"
57#include "Shards_CellTopology.hpp"
60using namespace Intrepid;
66int main(
int argc,
char *argv[]) {
68 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
71 int iprint = argc - 1;
73 Teuchos::RCP<std::ostream> outStream;
74 Teuchos::oblackholestream bhs;
77 outStream = Teuchos::rcp(&std::cout,
false);
79 outStream = Teuchos::rcp(&bhs,
false);
82 Teuchos::oblackholestream oldFormatState;
83 oldFormatState.copyfmt(std::cout);
86 <<
"===============================================================================\n" \
88 <<
"| Unit Test HCURL_TRI_In_FEM |\n" \
90 <<
"| 1) Tests triangular Nedelec-Thomas basis |\n" \
92 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
93 <<
"| Denis Ridzal (dridzal@sandia.gov) or |\n" \
94 <<
"| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
96 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
97 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
99 <<
"===============================================================================\n";
116 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
120 POINTTYPE_EQUISPACED );
122 myBasisDIV.
getValues( myBasisValuesDIV , lattice , OPERATOR_VALUE );
123 myBasisCURL.
getValues( myBasisValuesCURL, lattice , OPERATOR_VALUE );
125 for (
int i=0;i<myBasisValuesDIV.
dimension(0);i++) {
126 for (
int j=0;j<myBasisValuesDIV.
dimension(1);j++) {
127 if (std::abs( myBasisValuesDIV(i,j,1) + myBasisValuesCURL(i,j,0) ) > INTREPID_TOL ) {
129 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
131 *outStream <<
" At multi-index { ";
132 *outStream << i <<
" " << j <<
" and component 0";
133 *outStream <<
"} computed value: " << myBasisValuesCURL(i,j,0)
134 <<
" but correct value: " << -myBasisValuesDIV(i,j,1) <<
"\n";
135 *outStream <<
"Difference: " << std::abs( myBasisValuesCURL(i,j,0) + myBasisValuesDIV(i,j,1) ) <<
"\n";
137 if (std::abs( myBasisValuesDIV(i,j,0) - myBasisValuesCURL(i,j,1) ) > INTREPID_TOL ) {
139 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
141 *outStream <<
" At multi-index { ";
142 *outStream << i <<
" " << j <<
" and component 1";
143 *outStream <<
"} computed value: " << myBasisValuesCURL(i,j,1)
144 <<
" but correct value: " << myBasisValuesDIV(i,j,0) <<
"\n";
145 *outStream <<
"Difference: " << std::abs( myBasisValuesCURL(i,j,1) - myBasisValuesDIV(i,j,1) ) <<
"\n";
150 catch (
const std::exception & err) {
151 *outStream << err.what() <<
"\n\n";
161 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
165 POINTTYPE_EQUISPACED );
169 myBasis.
getValues( myBasisValues, lattice , OPERATOR_VALUE );
171 const double fiat_vals[] = {
1722.000000000000000e+00, 0.000000000000000e+00,
1735.000000000000000e-01, 2.500000000000000e-01,
174-1.000000000000000e+00, -1.000000000000000e+00,
1752.500000000000000e-01, 0.000000000000000e+00,
176-5.000000000000000e-01, -5.000000000000000e-01,
1770.000000000000000e+00, 0.000000000000000e+00,
178-1.000000000000000e+00, 0.000000000000000e+00,
1795.000000000000000e-01, 2.500000000000000e-01,
1802.000000000000000e+00, 2.000000000000000e+00,
181-5.000000000000000e-01, 0.000000000000000e+00,
1822.500000000000000e-01, 2.500000000000000e-01,
1830.000000000000000e+00, 0.000000000000000e+00,
1840.000000000000000e+00, 0.000000000000000e+00,
1850.000000000000000e+00, 2.500000000000000e-01,
1860.000000000000000e+00, 2.000000000000000e+00,
1875.000000000000000e-01, 0.000000000000000e+00,
188-2.500000000000000e-01, 2.500000000000000e-01,
1891.000000000000000e+00, 0.000000000000000e+00,
1900.000000000000000e+00, 0.000000000000000e+00,
1910.000000000000000e+00, -5.000000000000000e-01,
1920.000000000000000e+00, -1.000000000000000e+00,
193-2.500000000000000e-01, 0.000000000000000e+00,
194-2.500000000000000e-01, 2.500000000000000e-01,
195-2.000000000000000e+00, 0.000000000000000e+00,
1960.000000000000000e+00, 1.000000000000000e+00,
1970.000000000000000e+00, 5.000000000000000e-01,
1980.000000000000000e+00, 0.000000000000000e+00,
199-2.500000000000000e-01, -5.000000000000000e-01,
200-2.500000000000000e-01, -2.500000000000000e-01,
201-2.000000000000000e+00, -2.000000000000000e+00,
2020.000000000000000e+00, -2.000000000000000e+00,
2030.000000000000000e+00, -2.500000000000000e-01,
2040.000000000000000e+00, 0.000000000000000e+00,
205-2.500000000000000e-01, -5.000000000000000e-01,
2065.000000000000000e-01, 5.000000000000000e-01,
2071.000000000000000e+00, 1.000000000000000e+00,
2080.000000000000000e+00, 0.000000000000000e+00,
2090.000000000000000e+00, -7.500000000000000e-01,
2100.000000000000000e+00, 0.000000000000000e+00,
2111.500000000000000e+00, 0.000000000000000e+00,
2127.500000000000000e-01, 7.500000000000000e-01,
2130.000000000000000e+00, 0.000000000000000e+00,
2140.000000000000000e+00, 0.000000000000000e+00,
2150.000000000000000e+00, 1.500000000000000e+00,
2160.000000000000000e+00, 0.000000000000000e+00,
217-7.500000000000000e-01, 0.000000000000000e+00,
2187.500000000000000e-01, 7.500000000000000e-01,
2190.000000000000000e+00, 0.000000000000000e+00
223 for (
int i=0;i<myBasisValues.dimension(0);i++) {
224 for (
int j=0;j<myBasisValues.dimension(1);j++) {
225 for (
int k=0;k<myBasisValues.dimension(2);k++) {
226 if (std::abs( myBasisValues(i,j,k) - fiat_vals[cur] ) > INTREPID_TOL ) {
228 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
231 *outStream <<
" At multi-index { ";
232 *outStream << i <<
" " << j <<
" " << k;
233 *outStream <<
"} computed value: " << myBasisValues(i,j,k)
234 <<
" but correct value: " << fiat_vals[cur] <<
"\n";
235 *outStream <<
"Difference: " << std::abs( myBasisValues(i,j,k) - fiat_vals[cur] ) <<
"\n";
242 catch (
const std::exception & err) {
243 *outStream << err.what() <<
"\n\n";
252 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
256 POINTTYPE_EQUISPACED );
260 myBasis.
getValues( myBasisValues, lattice , OPERATOR_CURL );
262 const double fiat_curls[] = {
2637.000000000000000e+00,
2642.500000000000000e+00,
265-2.000000000000000e+00,
2662.500000000000000e+00,
267-2.000000000000000e+00,
268-2.000000000000000e+00,
269-2.000000000000000e+00,
2702.500000000000000e+00,
2717.000000000000000e+00,
272-2.000000000000000e+00,
2732.500000000000000e+00,
274-2.000000000000000e+00,
275-2.000000000000000e+00,
2762.500000000000000e+00,
2777.000000000000000e+00,
278-2.000000000000000e+00,
2792.500000000000000e+00,
280-2.000000000000000e+00,
281-2.000000000000000e+00,
282-2.000000000000000e+00,
283-2.000000000000000e+00,
2842.500000000000000e+00,
2852.500000000000000e+00,
2867.000000000000000e+00,
287-2.000000000000000e+00,
288-2.000000000000000e+00,
289-2.000000000000000e+00,
2902.500000000000000e+00,
2912.500000000000000e+00,
2927.000000000000000e+00,
2937.000000000000000e+00,
2942.500000000000000e+00,
295-2.000000000000000e+00,
2962.500000000000000e+00,
297-2.000000000000000e+00,
298-2.000000000000000e+00,
299-9.000000000000000e+00,
300-4.500000000000000e+00,
3010.000000000000000e+00,
3020.000000000000000e+00,
3034.500000000000000e+00,
3049.000000000000000e+00,
3059.000000000000000e+00,
3060.000000000000000e+00,
307-9.000000000000000e+00,
3084.500000000000000e+00,
309-4.500000000000000e+00,
314 for (
int i=0;i<myBasisValues.dimension(0);i++) {
315 for (
int j=0;j<myBasisValues.dimension(1);j++) {
316 if (std::abs( myBasisValues(i,j) - fiat_curls[cur] ) > INTREPID_TOL ) {
318 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
321 *outStream <<
" At multi-index { ";
322 *outStream << i <<
" " << j;
323 *outStream <<
"} computed value: " << myBasisValues(i,j)
324 <<
" but correct value: " << fiat_curls[cur] <<
"\n";
325 *outStream <<
"Difference: " << std::abs( myBasisValues(i,j) - fiat_curls[cur] ) <<
"\n";
331 catch (
const std::exception & err) {
332 *outStream << err.what() <<
"\n\n";
337 std::cout <<
"End Result: TEST FAILED\n";
339 std::cout <<
"End Result: TEST PASSED\n";
342 std::cout.copyfmt(oldFormatState);
int main(int argc, char *argv[])
Performs a code-code comparison to FIAT for Nedelec bases on triangles (values and curls)
Header file for utility class to provide multidimensional containers.
Header file for the Intrepid::HCURL_TRI_In_FEM class.
Header file for the Intrepid::HDIV_TRI_In_FEM class.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
Implementation of the default H(div)-compatible Raviart-Thomas basis of arbitrary degree on Triangle ...
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
virtual int getCardinality() const
Returns cardinality of the basis.
Implementation of a templated lexicographical container for a multi-indexed scalar quantity....
int dimension(const int whichDim) const
Returns the specified dimension.