$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 // C++ includes 00020 00021 00022 00023 // Local includes 00024 #include "libmesh/quadrature_clough.h" 00025 #include "libmesh/quadrature_gauss.h" 00026 #include "libmesh/quadrature_gm.h" 00027 #include "libmesh/quadrature_grid.h" 00028 #include "libmesh/quadrature_jacobi.h" 00029 #include "libmesh/quadrature_monomial.h" 00030 #include "libmesh/quadrature_simpson.h" 00031 #include "libmesh/quadrature_trap.h" 00032 #include "libmesh/quadrature_gauss_lobatto.h" 00033 #include "libmesh/quadrature_conical.h" 00034 #include "libmesh/string_to_enum.h" 00035 00036 namespace libMesh 00037 { 00038 00039 00040 00041 //--------------------------------------------------------------- 00042 UniquePtr<QBase> QBase::build (const std::string &type, 00043 const unsigned int _dim, 00044 const Order _order) 00045 { 00046 return QBase::build (Utility::string_to_enum<QuadratureType> (type), 00047 _dim, 00048 _order); 00049 } 00050 00051 00052 00053 UniquePtr<QBase> QBase::build(const QuadratureType _qt, 00054 const unsigned int _dim, 00055 const Order _order) 00056 { 00057 switch (_qt) 00058 { 00059 00060 case QCLOUGH: 00061 { 00062 #ifdef DEBUG 00063 if (_order > TWENTYTHIRD) 00064 { 00065 libMesh::out << "WARNING: Clough quadrature implemented" << std::endl 00066 << " up to TWENTYTHIRD order." << std::endl; 00067 } 00068 #endif 00069 00070 return UniquePtr<QBase>(new QClough(_dim, _order)); 00071 } 00072 00073 case QGAUSS: 00074 { 00075 00076 #ifdef DEBUG 00077 if (_order > FORTYTHIRD) 00078 { 00079 libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl 00080 << " up to FORTYTHIRD order." << std::endl; 00081 } 00082 #endif 00083 00084 return UniquePtr<QBase>(new QGauss(_dim, _order)); 00085 } 00086 00087 case QJACOBI_1_0: 00088 { 00089 00090 #ifdef DEBUG 00091 if (_order > FORTYTHIRD) 00092 { 00093 libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl 00094 << " up to FORTYTHIRD order." << std::endl; 00095 } 00096 00097 if (_dim > 1) 00098 { 00099 libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl 00100 << " in 1D only." << std::endl; 00101 } 00102 #endif 00103 00104 return UniquePtr<QBase>(new QJacobi(_dim, _order, 1, 0)); 00105 } 00106 00107 case QJACOBI_2_0: 00108 { 00109 00110 #ifdef DEBUG 00111 if (_order > FORTYTHIRD) 00112 { 00113 libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl 00114 << " up to FORTYTHIRD order." << std::endl; 00115 } 00116 00117 if (_dim > 1) 00118 { 00119 libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl 00120 << " in 1D only." << std::endl; 00121 } 00122 #endif 00123 00124 return UniquePtr<QBase>(new QJacobi(_dim, _order, 2, 0)); 00125 } 00126 00127 case QSIMPSON: 00128 { 00129 00130 #ifdef DEBUG 00131 if (_order > THIRD) 00132 { 00133 libMesh::out << "WARNING: Simpson rule provides only" << std::endl 00134 << " THIRD order!" << std::endl; 00135 } 00136 #endif 00137 00138 return UniquePtr<QBase>(new QSimpson(_dim)); 00139 } 00140 00141 case QTRAP: 00142 { 00143 00144 #ifdef DEBUG 00145 if (_order > FIRST) 00146 { 00147 libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl 00148 << " FIRST order!" << std::endl; 00149 } 00150 #endif 00151 00152 return UniquePtr<QBase>(new QTrap(_dim)); 00153 } 00154 00155 case QGRID: 00156 return UniquePtr<QBase>(new QGrid(_dim, _order)); 00157 00158 case QGRUNDMANN_MOLLER: 00159 return UniquePtr<QBase>(new QGrundmann_Moller(_dim, _order)); 00160 00161 case QMONOMIAL: 00162 return UniquePtr<QBase>(new QMonomial(_dim, _order)); 00163 00164 case QGAUSS_LOBATTO: 00165 return UniquePtr<QBase>(new QGaussLobatto(_dim, _order)); 00166 00167 case QCONICAL: 00168 return UniquePtr<QBase>(new QConical(_dim, _order)); 00169 00170 default: 00171 libmesh_error_msg("ERROR: Bad qt=" << _qt); 00172 } 00173 00174 00175 libmesh_error_msg("We'll never get here!"); 00176 return UniquePtr<QBase>(); 00177 } 00178 00179 } // namespace libMesh