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
00037 #include "Intrepid_FieldContainer.hpp"
00038 #include "Teuchos_oblackholestream.hpp"
00039 #include "Teuchos_RCP.hpp"
00040 #include "Teuchos_GlobalMPISession.hpp"
00041 #include "Intrepid_OrthogonalBases.hpp"
00042 #include <Intrepid_CubatureDirectTetDefault.hpp>
00043 #include <iostream>
00044 using namespace Intrepid;
00045
00050 int main(int argc, char *argv[]) {
00051
00052 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
00053
00054
00055 int iprint = argc - 1;
00056
00057 Teuchos::RCP<std::ostream> outStream;
00058 Teuchos::oblackholestream bhs;
00059
00060 if (iprint > 0)
00061 outStream = Teuchos::rcp(&std::cout, false);
00062 else
00063 outStream = Teuchos::rcp(&bhs, false);
00064
00065
00066 Teuchos::oblackholestream oldFormatState;
00067 oldFormatState.copyfmt(std::cout);
00068
00069 *outStream \
00070 << "===============================================================================\n" \
00071 << "| |\n" \
00072 << "| Unit Test OrthogonalBases |\n" \
00073 << "| |\n" \
00074 << "| 1) Tests orthogonality of tetrahedral orthogonal basis |\n" \
00075 << "| |\n" \
00076 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \
00077 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \
00078 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \
00079 << "| |\n" \
00080 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
00081 << "| Trilinos website: http://trilinos.sandia.gov |\n" \
00082 << "| |\n" \
00083 << "===============================================================================\n";
00084
00085 int errorFlag = 0;
00086
00087
00088
00089 CubatureDirectTetDefault<double,FieldContainer<double> > myCub(20);
00090 FieldContainer<double> cubPts( myCub.getNumPoints() , 3 );
00091 FieldContainer<double> cubWts( myCub.getNumPoints() );
00092
00093 myCub.getCubature( cubPts , cubWts );
00094
00095
00096 const int deg = 10;
00097 const int polydim = (deg+1)*(deg+2)*(deg+3)/6;
00098 FieldContainer<double> basisAtCubPts( polydim , myCub.getNumPoints() );
00099 OrthogonalBases::tabulateTetrahedron<double,FieldContainer<double>,FieldContainer<double> >( cubPts , deg , basisAtCubPts );
00100
00101
00102 for (int i=0;i<polydim;i++) {
00103 for (int j=0;j<polydim;j++) {
00104 double cur = 0;
00105 for (int k=0;k<myCub.getNumPoints();k++) {
00106 cur += cubWts(k) * basisAtCubPts( i , k ) * basisAtCubPts( j , k );
00107 }
00108 if (i != j && fabs( cur ) > 20.0 * INTREPID_TOL) {
00109 std::cout << INTREPID_TOL << std::endl;
00110 std::cout << i << " " << j << " " << cur << std::endl;
00111 errorFlag++;
00112 }
00113 else if (i == j && fabs( cur ) < 20.0 * INTREPID_TOL ) {
00114 std::cout << i << " " << j << " " << cur << std::endl;
00115 errorFlag++;
00116 }
00117
00118 }
00119 }
00120
00121
00122 if (errorFlag != 0)
00123 std::cout << "End Result: TEST FAILED\n";
00124 else
00125 std::cout << "End Result: TEST PASSED\n";
00126
00127
00128 std::cout.copyfmt(oldFormatState);
00129
00130 return errorFlag;
00131 }