|
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), 00025 // Robert C. Kirby (robert.c.kirby@ttu.edu) 00026 // Denis Ridzal (dridzal@sandia.gov), 00027 // Kara Peterson (kjpeter@sandia.gov). 00028 // 00029 // ************************************************************************ 00030 // @HEADER 00031 00054 // Intrepid includes 00055 #include "Intrepid_FunctionSpaceTools.hpp" 00056 #include "Intrepid_FieldContainer.hpp" 00057 #include "Intrepid_CellTools.hpp" 00058 #include "Intrepid_ArrayTools.hpp" 00059 #include "Intrepid_HGRAD_LINE_Cn_FEM.hpp" 00060 #include "Intrepid_HGRAD_QUAD_Cn_FEM.hpp" 00061 #include "Intrepid_HGRAD_HEX_Cn_FEM.hpp" 00062 #include "Intrepid_RealSpaceTools.hpp" 00063 #include "Intrepid_CubaturePolylib.hpp" 00064 #include "Intrepid_CubatureTensor.hpp" 00065 #include "Intrepid_Utils.hpp" 00066 00067 // Epetra includes 00068 #include "Epetra_Time.h" 00069 #include "Epetra_Map.h" 00070 #include "Epetra_FECrsMatrix.h" 00071 #include "Epetra_FEVector.h" 00072 #include "Epetra_SerialComm.h" 00073 00074 // Teuchos includes 00075 #include "Teuchos_oblackholestream.hpp" 00076 #include "Teuchos_RCP.hpp" 00077 #include "Teuchos_BLAS.hpp" 00078 00079 // Shards includes 00080 #include "Shards_CellTopology.hpp" 00081 00082 // EpetraExt includes 00083 #include "EpetraExt_RowMatrixOut.h" 00084 #include "EpetraExt_MultiVectorOut.h" 00085 00086 using namespace std; 00087 using namespace Intrepid; 00088 00089 int main(int argc, char *argv[]) { 00090 00091 //Check number of arguments 00092 if (argc < 2) { 00093 std::cout <<"\n>>> ERROR: Invalid number of arguments.\n\n"; 00094 std::cout <<"Usage:\n\n"; 00095 std::cout <<" ./Intrepid_example_Drivers_Example_05.exe min_deg max_deg verbose\n\n"; 00096 std::cout <<" where \n"; 00097 std::cout <<" int min_deg - beginning polynomial degree to check \n"; 00098 std::cout <<" int max_deg - final polynomial degree to check \n"; 00099 std::cout <<" verbose (optional) - any character, indicates verbose output \n\n"; 00100 exit(1); 00101 } 00102 00103 // This little trick lets us print to std::cout only if 00104 // a (dummy) command-line argument is provided. 00105 int iprint = argc - 1; 00106 Teuchos::RCP<std::ostream> outStream; 00107 Teuchos::oblackholestream bhs; // outputs nothing 00108 if (iprint > 1) 00109 outStream = Teuchos::rcp(&std::cout, false); 00110 else 00111 outStream = Teuchos::rcp(&bhs, false); 00112 00113 // Save the format state of the original std::cout. 00114 Teuchos::oblackholestream oldFormatState; 00115 oldFormatState.copyfmt(std::cout); 00116 00117 *outStream \ 00118 << "===============================================================================\n" \ 00119 << "| |\n" \ 00120 << "| Example: Check diagonalization of reference mass matrix |\n" \ 00121 << "| on line, quad, and hex |\n" \ 00122 << "| |\n" \ 00123 << "| Questions? Contact Pavel Bochev (pbboche@sandia.gov), |\n" \ 00124 << "| Denis Ridzal (dridzal@sandia.gov), |\n" \ 00125 << "| Kara Peterson (kjpeter@sandia.gov). |\n" \ 00126 << "| |\n" \ 00127 << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \ 00128 << "| Trilinos website: http://trilinos.sandia.gov |\n" \ 00129 << "| |\n" \ 00130 << "===============================================================================\n"; 00131 00132 00133 // ************************************ GET INPUTS ************************************** 00134 int min_degree = atoi(argv[1]); 00135 int max_degree = atoi(argv[2]); 00136 00137 00138 // *********************************** CELL TOPOLOGY ********************************** 00139 00140 // Get cell topology for base line 00141 typedef shards::CellTopology CellTopology; 00142 CellTopology line_2(shards::getCellTopologyData<shards::Line<2> >() ); 00143 CellTopology quad_4(shards::getCellTopologyData<shards::Quadrilateral<4> >() ); 00144 CellTopology hex_8(shards::getCellTopologyData<shards::Hexahedron<8> >() ); 00145 00146 std::vector<CellTopology> cts(3); 00147 00148 // loop over degrees 00149 for (int deg=min_degree;deg<=max_degree;deg++) { 00150 std::vector<Teuchos::RCP<Basis<double,FieldContainer<double> > > > bases; 00151 bases.push_back( Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<double, FieldContainer<double> >(deg,POINTTYPE_SPECTRAL)) ); 00152 bases.push_back( Teuchos::rcp(new Basis_HGRAD_QUAD_Cn_FEM<double, FieldContainer<double> >(deg,POINTTYPE_SPECTRAL)) ); 00153 bases.push_back( Teuchos::rcp(new Basis_HGRAD_HEX_Cn_FEM<double, FieldContainer<double> >(deg,POINTTYPE_SPECTRAL)) ); 00154 00155 // ************************************ CUBATURE ************************************** 00156 Teuchos::RCP<Cubature<double,FieldContainer<double>,FieldContainer<double> > > glcub 00157 = Teuchos::rcp(new CubaturePolylib<double,FieldContainer<double>,FieldContainer<double> >(2*deg-1,PL_GAUSS_LOBATTO) ); 00158 00159 std::vector<Teuchos::RCP<Cubature<double,FieldContainer<double>,FieldContainer<double> > > > 00160 cub_to_tensor; 00161 00162 // now we loop over spatial dimensions 00163 for (int sd=1;sd<=3;sd++) { 00164 // ************************************ CUBATURE ************************************** 00165 cub_to_tensor.push_back( glcub ); 00166 CubatureTensor<double,FieldContainer<double>,FieldContainer<double> > cubcur( cub_to_tensor ); 00167 00168 int cubDim = cubcur.getDimension(); 00169 int numCubPoints = cubcur.getNumPoints(); 00170 00171 FieldContainer<double> cubPoints(numCubPoints, cubDim); 00172 FieldContainer<double> cubWeights(numCubPoints); 00173 00174 cubcur.getCubature(cubPoints, cubWeights); 00175 00176 // ************************************ BASIS AT QP************************************** 00177 Teuchos::RCP<Basis<double,FieldContainer<double> > > basis_cur = bases[sd-1]; 00178 const int numFields = basis_cur->getCardinality(); 00179 FieldContainer<double> bf_at_cub_pts( numFields, numCubPoints ); 00180 FieldContainer<double> trans_bf_at_cub_pts( 1 , numFields,numCubPoints ); 00181 FieldContainer<double> w_trans_bf_at_cub_pts( 1, numFields , numCubPoints ); 00182 basis_cur->getValues( bf_at_cub_pts , cubPoints , OPERATOR_VALUE ); 00183 00184 // *********************************** FORM MASS MATRIX ***************************** 00185 FunctionSpaceTools::HGRADtransformVALUE<double>( trans_bf_at_cub_pts , 00186 bf_at_cub_pts ); 00187 cubWeights.resize(1,numCubPoints); 00188 FunctionSpaceTools::multiplyMeasure<double>( w_trans_bf_at_cub_pts , 00189 cubWeights , 00190 trans_bf_at_cub_pts ); 00191 cubWeights.resize(numCubPoints); 00192 00193 FieldContainer<double> mass_matrix( 1 , numFields, numFields ); 00194 FunctionSpaceTools::integrate<double>( mass_matrix , 00195 trans_bf_at_cub_pts , 00196 w_trans_bf_at_cub_pts , 00197 COMP_BLAS ); 00198 00199 // now we loop over the mass matrix and compare the nondiagonal entries to zero 00200 double max_offdiag = 0.0; 00201 for (int i=0;i<numFields;i++) { 00202 for (int j=0;j<numFields;j++) { 00203 if (i != j) { 00204 if ( abs(mass_matrix(0,i,j)) >= max_offdiag) { 00205 max_offdiag = abs(mass_matrix(0,i,j)); 00206 } 00207 } 00208 } 00209 } 00210 double min_diag = mass_matrix(0,0,0); 00211 for (int i=0;i<numFields;i++) { 00212 if ( mass_matrix(0,i,i) <= min_diag ) { 00213 min_diag = mass_matrix(0,i,i); 00214 } 00215 } 00216 *outStream << "Degree = " << deg << " and dimension = " << sd << "; Max offdiagonal" 00217 << " element is " << max_offdiag << " and min diagonal element is " << min_diag 00218 << ".\n"; 00219 00220 } 00221 00222 } 00223 00224 std::cout << "End Result: TEST PASSED\n"; 00225 00226 std::cout.copyfmt(oldFormatState); 00227 00228 return 0; 00229 } 00230
1.7.4