quadrature.h
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 00020 #ifndef LIBMESH_QUADRATURE_H 00021 #define LIBMESH_QUADRATURE_H 00022 00023 // Local includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/reference_counted_object.h" 00026 #include "libmesh/point.h" 00027 #include "libmesh/enum_elem_type.h" 00028 #include "libmesh/enum_order.h" 00029 #include "libmesh/enum_quadrature_type.h" 00030 #include "libmesh/auto_ptr.h" 00031 00032 // C++ includes 00033 #include <vector> 00034 #include <string> 00035 #include <utility> 00036 00037 namespace libMesh 00038 { 00039 00040 00041 00051 // ------------------------------------------------------------ 00052 // QBase class definition 00053 00054 class QBase : public ReferenceCountedObject<QBase> 00055 { 00056 protected: 00057 00062 QBase (const unsigned int _dim, 00063 const Order _order=INVALID_ORDER); 00064 00065 public: 00066 00070 virtual ~QBase() {} 00071 00075 virtual QuadratureType type() const = 0; 00076 00086 static AutoPtr<QBase> build (const std::string &name, 00087 const unsigned int _dim, 00088 const Order _order=INVALID_ORDER); 00089 00097 static AutoPtr<QBase> build (const QuadratureType _qt, 00098 const unsigned int _dim, 00099 const Order _order=INVALID_ORDER); 00100 00104 ElemType get_elem_type() const 00105 { return _type; } 00106 00110 unsigned int get_p_level() const 00111 { return _p_level; } 00112 00116 unsigned int n_points() const 00117 { libmesh_assert (!_points.empty()); 00118 return libmesh_cast_int<unsigned int>(_points.size()); } 00119 00123 unsigned int get_dim() const { return _dim; } 00124 00129 const std::vector<Point>& get_points() const { return _points; } 00130 00135 std::vector<Point>& get_points() { return _points; } 00136 00140 const std::vector<Real>& get_weights() const { return _weights; } 00141 00145 std::vector<Real>& get_weights() { return _weights; } 00146 00150 Point qp(const unsigned int i) const 00151 { libmesh_assert_less (i, _points.size()); return _points[i]; } 00152 00156 Real w(const unsigned int i) const 00157 { libmesh_assert_less (i, _weights.size()); return _weights[i]; } 00158 00163 void init (const ElemType type=INVALID_ELEM, 00164 unsigned int p_level=0); 00165 00169 Order get_order() const { return static_cast<Order>(_order + _p_level); } 00170 00175 void print_info(std::ostream& os=libMesh::out) const; 00176 00183 void scale(std::pair<Real, Real> old_range, 00184 std::pair<Real, Real> new_range); 00185 00189 friend std::ostream& operator << (std::ostream& os, const QBase& q); 00190 00198 virtual bool shapes_need_reinit() { return false; } 00199 00215 bool allow_rules_with_negative_weights; 00216 00217 protected: 00218 00219 00225 virtual void init_0D (const ElemType type=INVALID_ELEM, 00226 unsigned int p_level=0); 00227 00235 virtual void init_1D (const ElemType type=INVALID_ELEM, 00236 unsigned int p_level=0) = 0; 00237 00246 virtual void init_2D (const ElemType, 00247 unsigned int =0) 00248 #ifndef DEBUG 00249 {} 00250 #else 00251 { 00252 libMesh::err << "ERROR: Seems as if this quadrature rule" << std::endl 00253 << " is not implemented for 2D." << std::endl; 00254 libmesh_error(); 00255 } 00256 #endif 00257 00266 virtual void init_3D (const ElemType, 00267 unsigned int =0) 00268 #ifndef DEBUG 00269 {} 00270 #else 00271 { 00272 libMesh::err << "ERROR: Seems as if this quadrature rule" << std::endl 00273 << " is not implemented for 3D." << std::endl; 00274 libmesh_error(); 00275 } 00276 #endif 00277 00278 00285 void tensor_product_quad (const QBase& q1D); 00286 00293 void tensor_product_hex (const QBase& q1D); 00294 00302 void tensor_product_prism (const QBase& q1D, const QBase& q2D); 00303 00304 00305 00309 const unsigned int _dim; 00310 00314 const Order _order; 00315 00320 ElemType _type; 00321 00326 unsigned int _p_level; 00327 00332 std::vector<Point> _points; 00333 00337 std::vector<Real> _weights; 00338 }; 00339 00340 00341 00342 00343 00344 // ------------------------------------------------------------ 00345 // QBase class members 00346 00347 inline 00348 QBase::QBase(const unsigned int d, 00349 const Order o) : 00350 allow_rules_with_negative_weights(true), 00351 _dim(d), 00352 _order(o), 00353 _type(INVALID_ELEM), 00354 _p_level(0) 00355 { 00356 } 00357 00358 00359 00360 00361 inline 00362 void QBase::print_info(std::ostream& os) const 00363 { 00364 libmesh_assert(!_points.empty()); 00365 libmesh_assert(!_weights.empty()); 00366 00367 os << "N_Q_Points=" << this->n_points() << std::endl << std::endl; 00368 for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++) 00369 { 00370 os << " Point " << qpoint << ":\n" 00371 << " " 00372 << _points[qpoint] 00373 << " Weight:\n " 00374 << " w=" << _weights[qpoint] << "\n" << std::endl; 00375 } 00376 } 00377 00378 } // namespace libMesh 00379 00380 #endif // LIBMESH_QUADRATURE_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:48 UTC
Hosted By: