|
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_PointTools.hpp" 00042 #include "Intrepid_HGRAD_TET_Cn_FEM.hpp" 00043 #include "Shards_CellTopology.hpp" 00044 00045 #include <iostream> 00046 using namespace Intrepid; 00047 00048 00054 int main(int argc, char *argv[]) { 00055 00056 Teuchos::GlobalMPISession mpiSession(&argc, &argv); 00057 00058 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided. 00059 int iprint = argc - 1; 00060 00061 Teuchos::RCP<std::ostream> outStream; 00062 Teuchos::oblackholestream bhs; // outputs nothing 00063 00064 if (iprint > 0) 00065 outStream = Teuchos::rcp(&std::cout, false); 00066 else 00067 outStream = Teuchos::rcp(&bhs, false); 00068 00069 // Save the format state of the original std::cout. 00070 Teuchos::oblackholestream oldFormatState; 00071 oldFormatState.copyfmt(std::cout); 00072 00073 *outStream \ 00074 << "===============================================================================\n" \ 00075 << "| |\n" \ 00076 << "| Unit Test HGRAD_TET_Cn_FEM |\n" \ 00077 << "| |\n" \ 00078 << "| 1) Tests triangular Lagrange basis |\n" \ 00079 << "| |\n" \ 00080 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov) or |\n" \ 00081 << "| Denis Ridzal (dridzal@sandia.gov) or |\n" \ 00082 << "| Robert Kirby (robert.c.kirby@ttu.edu) |\n" \ 00083 << "| |\n" \ 00084 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00085 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00086 << "| |\n" \ 00087 << "===============================================================================\n"; 00088 00089 int errorFlag = 0; 00090 00091 // Let's instantiate a basis 00092 try { 00093 const int deg = 3; 00094 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_WARPBLEND ); 00095 00096 // Get a lattice 00097 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 ); 00098 const int nbf = myBasis.getCardinality(); 00099 FieldContainer<double> lattice( np_lattice , 3 ); 00100 PointTools::getLattice<double,FieldContainer<double> >( lattice , 00101 myBasis.getBaseCellTopology() , 00102 deg , 00103 0 , 00104 POINTTYPE_WARPBLEND ); 00105 FieldContainer<double> vals( nbf , np_lattice ); 00106 00107 myBasis.getValues( vals , lattice , OPERATOR_VALUE ); 00108 00109 // test for Kronecker property 00110 for (int i=0;i<nbf;i++) { 00111 for (int j=0;j<np_lattice;j++) { 00112 if ( i==j && std::abs( vals(i,j) - 1.0 ) > 100.0 * INTREPID_TOL ) { 00113 errorFlag++; 00114 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00115 *outStream << " Basis function " << i << " does not have unit value at its node\n"; 00116 } 00117 if ( i!=j && std::abs( vals(i,j) ) > 100.0 * INTREPID_TOL ) { 00118 errorFlag++; 00119 *outStream << std::setw(70) << "^^^^----FAILURE!" << "\n"; 00120 *outStream << " Basis function " << i << " does not vanish at node " << j << "\n"; 00121 *outStream << " Basis function value is " << vals(i,j) << "\n"; 00122 } 00123 } 00124 } 00125 } 00126 catch (std::exception err) { 00127 *outStream << err.what() << "\n\n"; 00128 errorFlag = -1000; 00129 } 00130 00131 try { 00132 const int deg = 4; 00133 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_WARPBLEND ); 00134 const std::vector<std::vector<std::vector<int> > >& dofdata = myBasis.getDofOrdinalData(); 00135 for (unsigned d=0;d<dofdata.size();d++) { 00136 std::cout << "Dimension " << d << "\n"; 00137 for (unsigned f=0;f<dofdata[d].size();f++) { 00138 std::cout << "\tFacet number " << f << "\n"; 00139 std::cout << "\t\tDOFS:\n"; 00140 for (unsigned n=0;n<dofdata[d][f].size();n++) { 00141 if ( dofdata[d][f][n] >= 0 ) { 00142 std::cout << "\t\t\t" << dofdata[d][f][n] << "\n"; 00143 } 00144 } 00145 } 00146 } 00147 } 00148 catch (std::exception err) { 00149 *outStream << err.what() << "\n\n"; 00150 errorFlag = -1000; 00151 } 00152 00153 try { 00154 const int deg = 3; 00155 Basis_HGRAD_TET_Cn_FEM<double,FieldContainer<double> > myBasis( deg , POINTTYPE_EQUISPACED ); 00156 00157 // Get a lattice 00158 const int np_lattice = PointTools::getLatticeSize( myBasis.getBaseCellTopology() , deg , 0 ); 00159 const int nbf = myBasis.getCardinality(); 00160 FieldContainer<double> lattice( np_lattice , 3 ); 00161 PointTools::getLattice<double,FieldContainer<double> >( lattice , 00162 myBasis.getBaseCellTopology() , 00163 deg , 00164 0 , 00165 POINTTYPE_EQUISPACED ); 00166 FieldContainer<double> vals( nbf , np_lattice , 3 ); 00167 00168 myBasis.getValues( vals , lattice , OPERATOR_GRAD ); 00169 00170 } 00171 catch (std::exception err) { 00172 *outStream << err.what() << "\n\n"; 00173 errorFlag = -1000; 00174 } 00175 00176 if (errorFlag != 0) 00177 std::cout << "End Result: TEST FAILED\n"; 00178 else 00179 std::cout << "End Result: TEST PASSED\n"; 00180 00181 // reset format state of std::cout 00182 std::cout.copyfmt(oldFormatState); 00183 00184 return errorFlag; 00185 }
1.7.4