quadrature_build.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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_jacobi.h" 00027 #include "libmesh/quadrature_simpson.h" 00028 #include "libmesh/quadrature_trap.h" 00029 #include "libmesh/string_to_enum.h" 00030 00031 namespace libMesh 00032 { 00033 00034 00035 00036 //--------------------------------------------------------------- 00037 AutoPtr<QBase> QBase::build (const std::string &type, 00038 const unsigned int _dim, 00039 const Order _order) 00040 { 00041 return QBase::build (Utility::string_to_enum<QuadratureType> (type), 00042 _dim, 00043 _order); 00044 } 00045 00046 00047 00048 AutoPtr<QBase> QBase::build(const QuadratureType _qt, 00049 const unsigned int _dim, 00050 const Order _order) 00051 { 00052 switch (_qt) 00053 { 00054 00055 case QCLOUGH: 00056 { 00057 #ifdef DEBUG 00058 if (_order > TWENTYTHIRD) 00059 { 00060 libMesh::out << "WARNING: Clough quadrature implemented" << std::endl 00061 << " up to TWENTYTHIRD order." << std::endl; 00062 } 00063 #endif 00064 00065 AutoPtr<QBase> ap(new QClough(_dim, _order)); 00066 return ap; 00067 } 00068 00069 case QGAUSS: 00070 { 00071 00072 #ifdef DEBUG 00073 if (_order > FORTYTHIRD) 00074 { 00075 libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl 00076 << " up to FORTYTHIRD order." << std::endl; 00077 } 00078 #endif 00079 00080 AutoPtr<QBase> ap(new QGauss(_dim, _order)); 00081 return ap; 00082 } 00083 00084 case QJACOBI_1_0: 00085 { 00086 00087 #ifdef DEBUG 00088 if (_order > TWENTYTHIRD) 00089 { 00090 libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl 00091 << " up to TWENTYTHIRD order." << std::endl; 00092 } 00093 00094 if (_dim > 1) 00095 { 00096 libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl 00097 << " in 1D only." << std::endl; 00098 } 00099 #endif 00100 00101 AutoPtr<QBase> ap(new QJacobi(_dim, _order, 1, 0)); 00102 return ap; 00103 } 00104 00105 case QJACOBI_2_0: 00106 { 00107 00108 #ifdef DEBUG 00109 if (_order > TWENTYTHIRD) 00110 { 00111 libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl 00112 << " up to TWENTYTHIRD order." << std::endl; 00113 } 00114 00115 if (_dim > 1) 00116 { 00117 libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl 00118 << " in 1D only." << std::endl; 00119 } 00120 #endif 00121 00122 AutoPtr<QBase> ap(new QJacobi(_dim, _order, 2, 0)); 00123 return ap; 00124 } 00125 00126 case QSIMPSON: 00127 { 00128 00129 #ifdef DEBUG 00130 if (_order > THIRD) 00131 { 00132 libMesh::out << "WARNING: Simpson rule provides only" << std::endl 00133 << " THIRD order!" << std::endl; 00134 } 00135 #endif 00136 00137 AutoPtr<QBase> ap(new QSimpson(_dim)); 00138 return ap; 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 AutoPtr<QBase> ap(new QTrap(_dim)); 00153 return ap; 00154 } 00155 00156 00157 default: 00158 { 00159 libMesh::err << "ERROR: Bad qt=" << _qt << std::endl; 00160 libmesh_error(); 00161 } 00162 } 00163 00164 00165 libmesh_error(); 00166 AutoPtr<QBase> ap(NULL); 00167 return ap; 00168 } 00169 00170 } // namespace libMesh 00171
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: