00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00038 #include "Intrepid_FieldContainer.hpp"
00039 #include "Teuchos_oblackholestream.hpp"
00040 #include "Teuchos_RCP.hpp"
00041 #include "Teuchos_GlobalMPISession.hpp"
00042 #include "Intrepid_OrthogonalBases.hpp"
00043 #include <Intrepid_CubatureDirectTriDefault.hpp>
00044 #include <iostream>
00045 using namespace Intrepid;
00046
00051 int main(int argc, char *argv[]) {
00052
00053 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00054
00055
00056 int iprint = argc - 1;
00057
00058 Teuchos::RCP<std::ostream> outStream;
00059 Teuchos::oblackholestream bhs;
00060
00061 if (iprint > 0)
00062 outStream = Teuchos::rcp(&std::cout, false);
00063 else
00064 outStream = Teuchos::rcp(&bhs, false);
00065
00066
00067 Teuchos::oblackholestream oldFormatState;
00068 oldFormatState.copyfmt(std::cout);
00069
00070 *outStream \
00071 << "===============================================================================\n" \
00072 << "| |\n" \
00073 << "| Unit Test OrthogonalBases |\n" \
00074 << "| |\n" \
00075 << "| 1) Tests orthogonality of triangular orthogonal basis (Dubiner) |\n" \
00076 << "| |\n" \
00077 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
00078 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
00079 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
00080 << "| |\n" \
00081 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00082 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00083 << "| |\n" \
00084 << "===============================================================================\n";
00085
00086 int errorFlag = 0;
00087
00088
00089
00090 CubatureDirectTriDefault<double,FieldContainer<double> > myCub(20);
00091 FieldContainer<double> cubPts( myCub.getNumPoints() , 2 );
00092 FieldContainer<double> cubWts( myCub.getNumPoints() );
00093
00094 myCub.getCubature( cubPts , cubWts );
00095
00096
00097 const int deg = 10;
00098 const int polydim = (deg+1)*(deg+2)/2;
00099 FieldContainer<double> basisAtCubPts( polydim , myCub.getNumPoints() );
00100 OrthogonalBases::tabulateTriangle<double,FieldContainer<double>,FieldContainer<double> >( cubPts , deg , basisAtCubPts );
00101
00102
00103 for (int i=0;i<polydim;i++) {
00104 for (int j=0;j<polydim;j++) {
00105 double cur = 0;
00106 for (int k=0;k<myCub.getNumPoints();k++) {
00107 cur += cubWts(k) * basisAtCubPts( i , k ) * basisAtCubPts( j , k );
00108 }
00109 if (i != j && fabs( cur ) > INTREPID_TOL) {
00110 errorFlag++;
00111 }
00112 else if (i == j && fabs( cur ) < INTREPID_TOL ) {
00113 errorFlag++;
00114 }
00115
00116 }
00117 }
00118
00119
00120 if (errorFlag != 0)
00121 std::cout << "End Result: TEST FAILED\n";
00122 else
00123 std::cout << "End Result: TEST PASSED\n";
00124
00125
00126 std::cout.copyfmt(oldFormatState);
00127
00128 return errorFlag;
00129 }