$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 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 // forward declarations 00042 class Elem; 00043 00053 // ------------------------------------------------------------ 00054 // QBase class definition 00055 00056 class QBase : public ReferenceCountedObject<QBase> 00057 { 00058 protected: 00059 00064 QBase (const unsigned int _dim, 00065 const Order _order=INVALID_ORDER); 00066 00067 public: 00068 00072 virtual ~QBase() {} 00073 00077 virtual QuadratureType type() const = 0; 00078 00088 static UniquePtr<QBase> build (const std::string &name, 00089 const unsigned int _dim, 00090 const Order _order=INVALID_ORDER); 00091 00099 static UniquePtr<QBase> build (const QuadratureType _qt, 00100 const unsigned int _dim, 00101 const Order _order=INVALID_ORDER); 00102 00106 ElemType get_elem_type() const 00107 { return _type; } 00108 00112 unsigned int get_p_level() const 00113 { return _p_level; } 00114 00118 unsigned int n_points() const 00119 { libmesh_assert (!_points.empty()); 00120 return cast_int<unsigned int>(_points.size()); } 00121 00125 unsigned int get_dim() const { return _dim; } 00126 00131 const std::vector<Point>& get_points() const { return _points; } 00132 00137 std::vector<Point>& get_points() { return _points; } 00138 00142 const std::vector<Real>& get_weights() const { return _weights; } 00143 00147 std::vector<Real>& get_weights() { return _weights; } 00148 00152 Point qp(const unsigned int i) const 00153 { libmesh_assert_less (i, _points.size()); return _points[i]; } 00154 00158 Real w(const unsigned int i) const 00159 { libmesh_assert_less (i, _weights.size()); return _weights[i]; } 00160 00165 virtual void init (const ElemType type=INVALID_ELEM, 00166 unsigned int p_level=0); 00167 00176 virtual void init (const Elem &elem, 00177 const std::vector<Real> &vertex_distance_func, 00178 unsigned int p_level=0); 00179 00183 Order get_order() const { return static_cast<Order>(_order + _p_level); } 00184 00189 void print_info(std::ostream& os=libMesh::out) const; 00190 00197 void scale(std::pair<Real, Real> old_range, 00198 std::pair<Real, Real> new_range); 00199 00203 friend std::ostream& operator << (std::ostream& os, const QBase& q); 00204 00212 virtual bool shapes_need_reinit() { return false; } 00213 00229 bool allow_rules_with_negative_weights; 00230 00231 protected: 00232 00233 00239 virtual void init_0D (const ElemType type=INVALID_ELEM, 00240 unsigned int p_level=0); 00241 00249 virtual void init_1D (const ElemType type=INVALID_ELEM, 00250 unsigned int p_level=0) = 0; 00251 00260 virtual void init_2D (const ElemType, 00261 unsigned int =0) 00262 #ifndef DEBUG 00263 {} 00264 #else 00265 { 00266 libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 2D."); 00267 } 00268 #endif 00269 00278 virtual void init_3D (const ElemType, 00279 unsigned int =0) 00280 #ifndef DEBUG 00281 {} 00282 #else 00283 { 00284 libmesh_error_msg("ERROR: Seems as if this quadrature rule \nis not implemented for 3D."); 00285 } 00286 #endif 00287 00288 00295 void tensor_product_quad (const QBase& q1D); 00296 00303 void tensor_product_hex (const QBase& q1D); 00304 00312 void tensor_product_prism (const QBase& q1D, const QBase& q2D); 00313 00314 00315 00319 const unsigned int _dim; 00320 00324 const Order _order; 00325 00330 ElemType _type; 00331 00336 unsigned int _p_level; 00337 00342 std::vector<Point> _points; 00343 00347 std::vector<Real> _weights; 00348 }; 00349 00350 00351 00352 00353 00354 // ------------------------------------------------------------ 00355 // QBase class members 00356 00357 inline 00358 QBase::QBase(const unsigned int d, 00359 const Order o) : 00360 allow_rules_with_negative_weights(true), 00361 _dim(d), 00362 _order(o), 00363 _type(INVALID_ELEM), 00364 _p_level(0) 00365 { 00366 } 00367 00368 00369 00370 00371 inline 00372 void QBase::print_info(std::ostream& os) const 00373 { 00374 libmesh_assert(!_points.empty()); 00375 libmesh_assert(!_weights.empty()); 00376 00377 Real summed_weights=0; 00378 os << "N_Q_Points=" << this->n_points() << std::endl << std::endl; 00379 for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++) 00380 { 00381 os << " Point " << qpoint << ":\n" 00382 << " " 00383 << _points[qpoint] 00384 << "\n Weight:\n " 00385 << " w=" << _weights[qpoint] << "\n" << std::endl; 00386 00387 summed_weights += _weights[qpoint]; 00388 } 00389 os << "Summed Weights: " << summed_weights << std::endl; 00390 } 00391 00392 } // namespace libMesh 00393 00394 #endif // LIBMESH_QUADRATURE_H