51 #include "Teuchos_oblackholestream.hpp"
52 #include "Teuchos_RCP.hpp"
53 #include "Teuchos_GlobalMPISession.hpp"
56 #include "Shards_CellTopology.hpp"
59 using namespace Intrepid;
67 int main(
int argc,
char *argv[]) {
69 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
72 int iprint = argc - 1;
74 Teuchos::RCP<std::ostream> outStream;
75 Teuchos::oblackholestream bhs;
78 outStream = Teuchos::rcp(&std::cout,
false);
80 outStream = Teuchos::rcp(&bhs,
false);
83 Teuchos::oblackholestream oldFormatState;
84 oldFormatState.copyfmt(std::cout);
87 <<
"===============================================================================\n" \
89 <<
"| Unit Test HGRAD_TET_Cn_FEM |\n" \
91 <<
"| 1) Tests triangular Lagrange basis |\n" \
93 <<
"| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
94 <<
"| Denis Ridzal (dridzal@sandia.gov) or |\n" \
95 <<
"| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
97 <<
"| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
98 <<
"| Trilinos website: http://trilinos.sandia.gov |\n" \
100 <<
"===============================================================================\n";
113 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
117 POINTTYPE_WARPBLEND );
120 myBasis.
getValues( vals , lattice , OPERATOR_VALUE );
123 for (
int i=0;i<nbf;i++) {
124 for (
int j=0;j<np_lattice;j++) {
125 if ( i==j && std::abs( vals(i,j) - 1.0 ) > 100.0 * INTREPID_TOL ) {
127 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
128 *outStream <<
" Basis function " << i <<
" does not have unit value at its node\n";
130 if ( i!=j && std::abs( vals(i,j) ) > 100.0 * INTREPID_TOL ) {
132 *outStream << std::setw(70) <<
"^^^^----FAILURE!" <<
"\n";
133 *outStream <<
" Basis function " << i <<
" does not vanish at node " << j <<
"\n";
134 *outStream <<
" Basis function value is " << vals(i,j) <<
"\n";
139 catch (
const std::exception & err) {
140 *outStream << err.what() <<
"\n\n";
147 const std::vector<std::vector<std::vector<int> > >& dofdata = myBasis.
getDofOrdinalData();
148 for (
unsigned d=0;d<dofdata.size();d++) {
149 std::cout <<
"Dimension " << d <<
"\n";
150 for (
unsigned f=0;f<dofdata[d].size();f++) {
151 std::cout <<
"\tFacet number " << f <<
"\n";
152 std::cout <<
"\t\tDOFS:\n";
153 for (
unsigned n=0;n<dofdata[d][f].size();n++) {
154 if ( dofdata[d][f][n] >= 0 ) {
155 std::cout <<
"\t\t\t" << dofdata[d][f][n] <<
"\n";
161 catch (
const std::exception & err) {
162 *outStream << err.what() <<
"\n\n";
174 PointTools::getLattice<double,FieldContainer<double> >( lattice ,
178 POINTTYPE_EQUISPACED );
181 myBasis.
getValues( vals , lattice , OPERATOR_GRAD );
184 catch (
const std::exception & err) {
185 *outStream << err.what() <<
"\n\n";
190 std::cout <<
"End Result: TEST FAILED\n";
192 std::cout <<
"End Result: TEST PASSED\n";
195 std::cout.copyfmt(oldFormatState);
void getValues(ArrayScalar &outputValues, const ArrayScalar &inputPoints, const EOperator operatorType) const
Evaluation of a FEM basis on a reference Triangle cell.
virtual int getCardinality() const
Returns cardinality of the basis.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Tetrahedron ce...
Header file for utility class to provide multidimensional containers.
virtual const shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation http://trilin...
Header file for the Intrepid::HGRAD_TET_Cn_FEM class.
virtual const std::vector< std::vector< std::vector< int > > > & getDofOrdinalData()
DoF tag to ordinal data structure.