$extrastylesheet
libMesh::QBase Class Reference

#include <quadrature.h>

Inheritance diagram for libMesh::QBase:

List of all members.

Public Member Functions

virtual ~QBase ()
virtual QuadratureType type () const =0
ElemType get_elem_type () const
unsigned int get_p_level () const
unsigned int n_points () const
unsigned int get_dim () const
const std::vector< Point > & get_points () const
std::vector< Point > & get_points ()
const std::vector< Real > & get_weights () const
std::vector< Real > & get_weights ()
Point qp (const unsigned int i) const
Real w (const unsigned int i) const
virtual void init (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
virtual void init (const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0)
Order get_order () const
void print_info (std::ostream &os=libMesh::out) const
void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
virtual bool shapes_need_reinit ()

Static Public Member Functions

static UniquePtr< QBasebuild (const std::string &name, const unsigned int _dim, const Order _order=INVALID_ORDER)
static UniquePtr< QBasebuild (const QuadratureType _qt, const unsigned int _dim, const Order _order=INVALID_ORDER)
static std::string get_info ()
static void print_info (std::ostream &out=libMesh::out)
static unsigned int n_objects ()
static void enable_print_counter_info ()
static void disable_print_counter_info ()

Public Attributes

bool allow_rules_with_negative_weights

Protected Types

typedef std::map< std::string,
std::pair< unsigned int,
unsigned int > > 
Counts

Protected Member Functions

 QBase (const unsigned int _dim, const Order _order=INVALID_ORDER)
virtual void init_0D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
virtual void init_1D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)=0
virtual void init_2D (const ElemType, unsigned int=0)
 libmesh_error_msg ("ERROR: Seems as if this quadrature rule \nis not implemented for 2D.")
virtual void init_3D (const ElemType, unsigned int=0)
void tensor_product_hex (const QBase &q1D)
void tensor_product_prism (const QBase &q1D, const QBase &q2D)
void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)

Protected Attributes

const unsigned int _dim
const Order _order
ElemType _type
unsigned int _p_level
std::vector< Point_points
std::vector< Real_weights

Static Protected Attributes

static Counts _counts
static Threads::atomic
< unsigned int > 
_n_objects
static Threads::spin_mutex _mutex
static bool _enable_print_counter = true

Friends

std::ostream & operator<< (std::ostream &os, const QBase &q)

Detailed Description

This is the QBase class. It provides the basic functionality from which various quadrature rules can be derived. The class contains dim dimensional points describing the quadrature locations (referenced to a master object) and associated weights.

Author:
Benjamin S. Kirk, 2002

Definition at line 56 of file quadrature.h.


Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited]

Data structure to log the information. The log is identified by the class name.

Definition at line 113 of file reference_counter.h.


Constructor & Destructor Documentation

libMesh::QBase::QBase ( const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
) [inline, protected]

Constructor. Protected to prevent instantiation of this base class.

Definition at line 358 of file quadrature.h.

virtual libMesh::QBase::~QBase ( ) [inline, virtual]

Destructor.

Definition at line 72 of file quadrature.h.

{}

Member Function Documentation

UniquePtr< QBase > libMesh::QBase::build ( const std::string &  name,
const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
) [static]

Builds a specific quadrature rule, identified through the name string. An UniquePtr<QBase> is returned to prevent a memory leak. This way the user need not remember to delete the object. Enables run-time decision of the quadrature rule. The input parameter name must be mappable through the Utility::string_to_enum<>() function.

Definition at line 42 of file quadrature_build.C.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().

{
  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
                       _dim,
                       _order);
}
UniquePtr< QBase > libMesh::QBase::build ( const QuadratureType  _qt,
const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
) [static]

Builds a specific quadrature rule, identified through the QuadratureType. An UniquePtr<QBase> is returned to prevent a memory leak. This way the user need not remember to delete the object. Enables run-time decision of the quadrature rule.

Definition at line 53 of file quadrature_build.C.

References libMesh::FIRST, libMesh::FORTYTHIRD, libmesh_error_msg(), libMesh::out, libMesh::QCLOUGH, libMesh::QCONICAL, libMesh::QGAUSS, libMesh::QGAUSS_LOBATTO, libMesh::QGRID, libMesh::QGRUNDMANN_MOLLER, libMesh::QJACOBI_1_0, libMesh::QJACOBI_2_0, libMesh::QMONOMIAL, libMesh::QSIMPSON, libMesh::QTRAP, libMesh::THIRD, and libMesh::TWENTYTHIRD.

{
  switch (_qt)
    {

    case QCLOUGH:
      {
#ifdef DEBUG
        if (_order > TWENTYTHIRD)
          {
            libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
                         << " up to TWENTYTHIRD order." << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QClough(_dim, _order));
      }

    case QGAUSS:
      {

#ifdef DEBUG
        if (_order > FORTYTHIRD)
          {
            libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
                         << " up to FORTYTHIRD order." << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QGauss(_dim, _order));
      }

    case QJACOBI_1_0:
      {

#ifdef DEBUG
        if (_order > FORTYTHIRD)
          {
            libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
                         << " up to FORTYTHIRD order." << std::endl;
          }

        if (_dim > 1)
          {
            libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
                         << " in 1D only." << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QJacobi(_dim, _order, 1, 0));
      }

    case QJACOBI_2_0:
      {

#ifdef DEBUG
        if (_order > FORTYTHIRD)
          {
            libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
                         << " up to FORTYTHIRD order." << std::endl;
          }

        if (_dim > 1)
          {
            libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
                         << " in 1D only." << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QJacobi(_dim, _order, 2, 0));
      }

    case QSIMPSON:
      {

#ifdef DEBUG
        if (_order > THIRD)
          {
            libMesh::out << "WARNING: Simpson rule provides only" << std::endl
                         << " THIRD order!" << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QSimpson(_dim));
      }

    case QTRAP:
      {

#ifdef DEBUG
        if (_order > FIRST)
          {
            libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
                         << " FIRST order!" << std::endl;
          }
#endif

        return UniquePtr<QBase>(new QTrap(_dim));
      }

    case QGRID:
      return UniquePtr<QBase>(new QGrid(_dim, _order));

    case QGRUNDMANN_MOLLER:
      return UniquePtr<QBase>(new QGrundmann_Moller(_dim, _order));

    case QMONOMIAL:
      return UniquePtr<QBase>(new QMonomial(_dim, _order));

    case QGAUSS_LOBATTO:
      return UniquePtr<QBase>(new QGaussLobatto(_dim, _order));

    case QCONICAL:
      return UniquePtr<QBase>(new QConical(_dim, _order));

    default:
      libmesh_error_msg("ERROR: Bad qt=" << _qt);
    }


  libmesh_error_msg("We'll never get here!");
  return UniquePtr<QBase>();
}

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

{
  _enable_print_counter = true;
  return;
}
unsigned int libMesh::QBase::get_dim ( ) const [inline]
Returns:
the current element type we're set up for

Definition at line 106 of file quadrature.h.

References _type.

  { return _type; }
std::string libMesh::ReferenceCounter::get_info ( ) [static, inherited]

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)

  std::ostringstream oss;

  oss << '\n'
      << " ---------------------------------------------------------------------------- \n"
      << "| Reference count information                                                |\n"
      << " ---------------------------------------------------------------------------- \n";

  for (Counts::iterator it = _counts.begin();
       it != _counts.end(); ++it)
    {
      const std::string name(it->first);
      const unsigned int creations    = it->second.first;
      const unsigned int destructions = it->second.second;

      oss << "| " << name << " reference count information:\n"
          << "|  Creations:    " << creations    << '\n'
          << "|  Destructions: " << destructions << '\n';
    }

  oss << " ---------------------------------------------------------------------------- \n";

  return oss.str();

#else

  return "";

#endif
}
Order libMesh::QBase::get_order ( ) const [inline]
Returns:
the order of the quadrature rule.

Definition at line 183 of file quadrature.h.

References _order, and _p_level.

Referenced by libMesh::InfFE< Dim, T_radial, T_map >::attach_quadrature_rule().

{ return static_cast<Order>(_order + _p_level); }
unsigned int libMesh::QBase::get_p_level ( ) const [inline]
Returns:
the current p refinement level we're initialized with

Definition at line 112 of file quadrature.h.

References _p_level.

  { return _p_level; }
const std::vector<Point>& libMesh::QBase::get_points ( ) const [inline]
Returns:
a std::vector containing the quadrature point locations on a reference object.

Definition at line 131 of file quadrature.h.

References _points.

Referenced by libMesh::QClough::init_1D(), libMesh::QMonomial::init_1D(), libMesh::QClough::init_2D(), libMesh::QGauss::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QGauss::init_3D(), and libMesh::QMonomial::init_3D().

{ return _points;  }
std::vector<Point>& libMesh::QBase::get_points ( ) [inline]
Returns:
a std::vector containing the quadrature point locations on a reference object as a writeable reference.

Definition at line 137 of file quadrature.h.

References _points.

{ return _points;  }
std::vector<Real>& libMesh::QBase::get_weights ( ) [inline]
Returns:
a std::vector containing the quadrature weights.

Definition at line 147 of file quadrature.h.

References _weights.

{ return _weights; }
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 163 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.first++;
}
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 176 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.second++;
}
void libMesh::QBase::init ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
) [virtual]
void libMesh::QBase::init ( const Elem elem,
const std::vector< Real > &  vertex_distance_func,
unsigned int  p_level = 0 
) [virtual]

Initializes the data structures for a specific, potentially cut element. The array vertex_distance_func contains vertex values of a signed distance function that cuts the element. This interface is indended to be extended by derived classes that can cut the element into subelements, for example, and constuct a composite quadrature rule for the cut element.

Definition at line 72 of file quadrature.C.

References init(), and libMesh::Elem::type().

{
  // dispatch generic implementation
  this->init(elem.type(), p_level);
}
void libMesh::QBase::init_0D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
) [protected, virtual]

Initializes the 0D quadrature rule by filling the points and weights vectors with the appropriate values. Generally this is just one point with weight 1.

Definition at line 82 of file quadrature.C.

References _points, and _weights.

Referenced by init().

{
  _points.resize(1);
  _weights.resize(1);
  _points[0] = Point(0.);
  _weights[0] = 1.0;
}
virtual void libMesh::QBase::init_1D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
) [protected, pure virtual]

Initializes the 1D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. It is assumed that derived quadrature rules will at least define the init_1D function, therefore it is pure virtual.

Implemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QJacobi, libMesh::QGrid, libMesh::QSimpson, libMesh::QConical, libMesh::QGauss, libMesh::QClough, libMesh::QTrap, and libMesh::QGaussLobatto.

Referenced by init().

virtual void libMesh::QBase::init_2D ( const ElemType  ,
unsigned int  = 0 
) [inline, protected, virtual]

Initializes the 2D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. Should not be pure virtual since a derived quadrature rule may only be defined in 1D. If not redefined, gives an error (when DEBUG defined) when called.

Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QGrid, libMesh::QConical, libMesh::QSimpson, libMesh::QGauss, libMesh::QClough, libMesh::QTrap, and libMesh::QGaussLobatto.

Definition at line 260 of file quadrature.h.

Referenced by init().

  {}
virtual void libMesh::QBase::init_3D ( const ElemType  ,
unsigned int  = 0 
) [inline, protected, virtual]

Initializes the 3D quadrature rule by filling the points and weights vectors with the appropriate values. The order of the rule will be defined by the implementing class. Should not be pure virtual since a derived quadrature rule may only be defined in 1D. If not redefined, gives an error (when DEBUG defined) when called.

Reimplemented in libMesh::QGrundmann_Moller, libMesh::QMonomial, libMesh::QConical, libMesh::QGrid, libMesh::QSimpson, libMesh::QGauss, libMesh::QClough, libMesh::QTrap, and libMesh::QGaussLobatto.

Definition at line 278 of file quadrature.h.

Referenced by init().

  {}
static unsigned int libMesh::ReferenceCounter::n_objects ( ) [inline, static, inherited]

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

Referenced by libMesh::LibMeshInit::~LibMeshInit().

  { return _n_objects; }
unsigned int libMesh::QBase::n_points ( ) const [inline]
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out) [static, inherited]

Prints the reference information, by default to libMesh::out.

Definition at line 88 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

Referenced by libMesh::LibMeshInit::~LibMeshInit().

void libMesh::QBase::print_info ( std::ostream &  os = libMesh::out) const [inline]

Prints information relevant to the quadrature rule, by default to libMesh::out.

Definition at line 372 of file quadrature.h.

References _points, _weights, libMesh::libmesh_assert(), n_points(), and libMesh::Real.

Referenced by libMesh::operator<<().

{
  libmesh_assert(!_points.empty());
  libmesh_assert(!_weights.empty());

  Real summed_weights=0;
  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
  for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
    {
      os << " Point " << qpoint << ":\n"
         << "  "
         << _points[qpoint]
         << "\n Weight:\n "
         << "  w=" << _weights[qpoint] << "\n" << std::endl;

      summed_weights += _weights[qpoint];
    }
  os << "Summed Weights: " << summed_weights << std::endl;
}
Point libMesh::QBase::qp ( const unsigned int  i) const [inline]
Returns:
the $ i^{th} $ quadrature point on the reference object.

Definition at line 152 of file quadrature.h.

References _points.

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), tensor_product_hex(), and tensor_product_prism().

  { libmesh_assert_less (i, _points.size()); return _points[i]; }
void libMesh::QBase::scale ( std::pair< Real, Real old_range,
std::pair< Real, Real new_range 
)

Maps the points of a 1D interval quadrature rule (typically [-1,1]) to any other 1D interval (typically [0,1]) and scales the weights accordingly. The quadrature rule will be mapped from the entries of old_range to the entries of new_range.

Definition at line 93 of file quadrature.C.

References _dim, _points, _weights, and libMesh::Real.

Referenced by libMesh::QConical::conical_product_tet(), and libMesh::QConical::conical_product_tri().

{
  // Make sure we are in 1D
  libmesh_assert_equal_to (_dim, 1);

  Real
    h_new = new_range.second - new_range.first,
    h_old = old_range.second - old_range.first;

  // Make sure that we have sane ranges
  libmesh_assert_greater (h_new, 0.);
  libmesh_assert_greater (h_old, 0.);

  // Make sure there are some points
  libmesh_assert_greater (_points.size(), 0);

  // Compute the scale factor
  Real scfact = h_new/h_old;

  // We're mapping from old_range -> new_range
  for (unsigned int i=0; i<_points.size(); i++)
    {
      _points[i](0) = new_range.first +
        (_points[i](0) - old_range.first) * scfact;

      // Scale the weights
      _weights[i] *= scfact;
    }
}
virtual bool libMesh::QBase::shapes_need_reinit ( ) [inline, virtual]

Returns true if the shape functions need to be recalculated.

This can happen if the number of points or their positions change.

By default this will return false.

Definition at line 212 of file quadrature.h.

{ return false; }
void libMesh::QBase::tensor_product_hex ( const QBase q1D) [protected]

Computes the tensor product quadrature rule [q1D x q1D x q1D] from the 1D rule q1D. Used in the init_3D routines for hexahedral element types.

Definition at line 154 of file quadrature.C.

References _points, _weights, n_points(), qp(), and w().

Referenced by libMesh::QGaussLobatto::init_3D(), libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGrid::init_3D().

{
  const unsigned int np = q1D.n_points();

  _points.resize(np * np * np);

  _weights.resize(np * np * np);

  unsigned int q=0;

  for (unsigned int k=0; k<np; k++)
    for (unsigned int j=0; j<np; j++)
      for (unsigned int i=0; i<np; i++)
        {
          _points[q](0) = q1D.qp(i)(0);
          _points[q](1) = q1D.qp(j)(0);
          _points[q](2) = q1D.qp(k)(0);

          _weights[q] = q1D.w(i) * q1D.w(j) * q1D.w(k);

          q++;
        }
}
void libMesh::QBase::tensor_product_prism ( const QBase q1D,
const QBase q2D 
) [protected]

Computes the tensor product of a 1D quadrature rule and a 2D quadrature rule. Used in the init_3D routines for prismatic element types.

Definition at line 181 of file quadrature.C.

References _points, _weights, n_points(), qp(), and w().

Referenced by libMesh::QTrap::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), and libMesh::QGrid::init_3D().

{
  const unsigned int n_points1D = q1D.n_points();
  const unsigned int n_points2D = q2D.n_points();

  _points.resize  (n_points1D * n_points2D);
  _weights.resize (n_points1D * n_points2D);

  unsigned int q=0;

  for (unsigned int j=0; j<n_points1D; j++)
    for (unsigned int i=0; i<n_points2D; i++)
      {
        _points[q](0) = q2D.qp(i)(0);
        _points[q](1) = q2D.qp(i)(1);
        _points[q](2) = q1D.qp(j)(0);

        _weights[q] = q2D.w(i) * q1D.w(j);

        q++;
      }

}
Real libMesh::QBase::w ( const unsigned int  i) const [inline]
Returns:
the $ i^{th} $ quadrature weight.

Definition at line 158 of file quadrature.h.

References _weights.

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), tensor_product_hex(), and tensor_product_prism().

  { libmesh_assert_less (i, _weights.size()); return _weights[i]; }

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const QBase q 
) [friend]

Same as above, but allows you to use the stream syntax.

Definition at line 208 of file quadrature.C.

{
  q.print_info(os);
  return os;
}

Member Data Documentation

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

Flag to control whether reference count information is printed when print_info is called.

Definition at line 137 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited]

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 126 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

unsigned int libMesh::QBase::_p_level [protected]

The p level of element for which the current values have been computed.

Definition at line 336 of file quadrature.h.

Referenced by get_order(), get_p_level(), and init().

The type of element for which the current values have been computed.

Definition at line 330 of file quadrature.h.

Referenced by get_elem_type(), and init().

Flag (default true) controlling the use of quadrature rules with negative weights. Set this to false to ONLY use (potentially) safer but more expensive rules with all positive weights.

Negative weights typically appear in Gaussian quadrature rules over three-dimensional elements. Rules with negative weights can be unsuitable for some problems. For example, it is possible for a rule with negative weights to obtain a negative result when integrating a positive function.

A particular example: if rules with negative weights are not allowed, a request for TET,THIRD (5 points) will return the TET,FIFTH (14 points) rule instead, nearly tripling the computational effort required!

Definition at line 229 of file quadrature.h.

Referenced by libMesh::QGrundmann_Moller::init_2D(), libMesh::QGauss::init_3D(), libMesh::QMonomial::init_3D(), and libMesh::QGrundmann_Moller::init_3D().


The documentation for this class was generated from the following files: