|
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 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 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided. 00055 int iprint = argc - 1; 00056 00057 Teuchos::RCP<std::ostream> outStream; 00058 Teuchos::oblackholestream bhs; // outputs nothing 00059 00060 if (iprint > 0) 00061 outStream = Teuchos::rcp(&std::cout, false); 00062 else 00063 outStream = Teuchos::rcp(&bhs, false); 00064 00065 // Save the format state of the original std::cout. 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 // First, get a reference quadrature rule 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 // Tabulate the basis functions at the cubature points 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 // Now let's compute the mass matrix 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 // reset format state of std::cout 00128 std::cout.copyfmt(oldFormatState); 00129 00130 return errorFlag; 00131 }
1.7.4