|
Intrepid
|
00001 // @HEADER 00002 // ************************************************************************ 00003 // 00004 // Intrepid Package 00005 // Copyright (2007) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 00023 // USA 00024 // Questions? Contact Pavel Bochev (pbboche@sandia.gov) or 00025 // Denis Ridzal (dridzal@sandia.gov) or 00026 // Robert Kirby (robert.c.kirby@ttu.edu) 00027 // 00028 // ************************************************************************ 00029 // @HEADER 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 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided. 00056 int iprint = argc - 1; 00057 00058 Teuchos::RCP<std::ostream> outStream; 00059 Teuchos::oblackholestream bhs; // outputs nothing 00060 00061 if (iprint > 0) 00062 outStream = Teuchos::rcp(&std::cout, false); 00063 else 00064 outStream = Teuchos::rcp(&bhs, false); 00065 00066 // Save the format state of the original std::cout. 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 // First, get a reference quadrature rule 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 // Tabulate the basis functions at the cubature points 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 // Now let's compute the mass matrix 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 // reset format state of std::cout 00126 std::cout.copyfmt(oldFormatState); 00127 00128 return errorFlag; 00129 }
1.7.4