Poster of Linux kernelThe best gift for a Linux geek
QConical

QConical

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

QConical -  

SYNOPSIS


#include <quadrature_conical.h>

Inherits QBase.  

Public Member Functions


QConical (const unsigned int _dim, const Order _order=INVALID_ORDER)

~QConical ()

QuadratureType type () const

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

void init (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)

Order get_order () const

void print_info (std::ostream &os=std::cout) const

void scale (std::pair< Real, Real > old_range, std::pair< Real, Real > new_range)
 

Static Public Member Functions


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

static AutoPtr< QBase > build (const QuadratureType _qt, const unsigned int _dim, const Order _order=INVALID_ORDER)

static void print_info ()

static std::string get_info ()

static unsigned int n_objects ()
 

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


virtual void init_0D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)

void increment_constructor_count (const std::string &name)

void increment_destructor_count (const std::string &name)
 

Protected Attributes


std::cerr<< 'ERROR: Seems as if this quadrature rule'<< std::endl<< ' is not implemented for 2D.'<< std::endl;libmesh_error();}#endif virtual void init_3D(const ElemType, unsigned int=0)#ifndef DEBUG{}#else{std::cerr<< 'ERROR: Seems as if this quadrature rule'<< std::endl<< ' is not implemented for 3D.'<< std::endl;libmesh_error();}#endif void tensor_product_quad(const QBase &q1D);void tensor_product_hex(const QBase &q1D);void tensor_product_prism(const QBase &q1D, const QBase &q2D);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
 

Private Member Functions


void init_1D (const ElemType, unsigned int=0)

void init_2D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)

void init_3D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)

void conical_product_tri (unsigned int p)

void conical_product_tet (unsigned int p)

void conical_product_pyramid (unsigned int p)
 

Friends


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

Detailed Description

This class implements the so-called conical product quadrature rules for Tri and Tet elements. These rules are generally non-optimal in the number of evaluation points, but have the nice property of having all positive weights and being well-defined to any order for which their underlying 1D Gauss and Jacobi quadrature rules are available.

The construction of these rules is given by e.g.

Stroud, A.H. 'Approximate Calculation of
 Multiple Integrals.', 1972 

Definition at line 40 of file quadrature_conical.h.  

Member Typedef Documentation

 

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.

Definition at line 105 of file reference_counter.h.  

Constructor & Destructor Documentation

 

QConical::QConical (const unsigned int_dim, const Order_order = INVALID_ORDER)Constructor. Declares the order of the quadrature rule.

Definition at line 34 of file quadrature_conical.C.

                                  : QBase(d,o)
{
}
 

QConical::~QConical ()Destructor.

Definition at line 42 of file quadrature_conical.C.

{
}
 

Member Function Documentation

 

AutoPtr< QBase > QBase::build (const std::string &name, const unsigned int_dim, const Order_order = INVALID_ORDER) [static, inherited]Builds a specific quadrature rule, identified through the name string. An AutoPtr<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 36 of file quadrature_build.C.

References Utility::string_to_enum< QuadratureType >().

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

{
  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
                       _dim,
                       _order);
}
 

AutoPtr< QBase > QBase::build (const QuadratureType_qt, const unsigned int_dim, const Order_order = INVALID_ORDER) [static, inherited]Builds a specific quadrature rule, identified through the QuadratureType. An AutoPtr<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 47 of file quadrature_build.C.

References libMeshEnums::FIRST, libMeshEnums::FORTYTHIRD, libMeshEnums::QCLOUGH, libMeshEnums::QGAUSS, libMeshEnums::QJACOBI_1_0, libMeshEnums::QJACOBI_2_0, libMeshEnums::QSIMPSON, libMeshEnums::QTRAP, libMeshEnums::THIRD, and libMeshEnums::TWENTYTHIRD.

{
  switch (_qt)
    {
      
    case QCLOUGH:
      {
#ifdef DEBUG
        if (_order > TWENTYTHIRD)
          {
            std::cout << 'WARNING: Clough quadrature implemented' << std::endl
                      << ' up to TWENTYTHIRD order.' << std::endl;
          }
#endif

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

    case QGAUSS:
      {

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

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

    case QJACOBI_1_0:
      {

#ifdef DEBUG
        if (_order > TWENTYTHIRD)
          {
            std::cout << 'WARNING: Jacobi(1,0) quadrature implemented' << std::endl
                      << ' up to TWENTYTHIRD order.' << std::endl;
          }

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

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

    case QJACOBI_2_0:
      {

#ifdef DEBUG
        if (_order > TWENTYTHIRD)
          {
            std::cout << 'WARNING: Jacobi(2,0) quadrature implemented' << std::endl
                      << ' up to TWENTYTHIRD order.' << std::endl;
          }

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

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

    case QSIMPSON:
      {

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

        AutoPtr<QBase> ap(new QSimpson(_dim));
        return ap;
      }

    case QTRAP:
      {

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

        AutoPtr<QBase> ap(new QTrap(_dim));
        return ap;
      }


    default:
      { 
        std::cerr << 'ERROR: Bad qt=' << _qt << std::endl;
        libmesh_error();
      }
    }


  libmesh_error();
  AutoPtr<QBase> ap(NULL);
  return ap;
}
 

void QConical::conical_product_pyramid (unsigned intp) [private]Implementation of conical product rule for a Pyramid in 3D of order = _order+2*p.

Definition at line 174 of file quadrature_conical.C.

References QBase::_points, QBase::_weights, QBase::get_dim(), QBase::n_points(), QBase::qp(), and QBase::w().

Referenced by init_3D().

{
  // Be sure the underlying rule object was built with the same dimension as the
  // rule we are about to construct.
  libmesh_assert (this->get_dim() == 3);
  
  QGauss  gauss1D(1,static_cast<Order>(_order+2*p));
  QJacobi jac1D(1,static_cast<Order>(_order+2*p),2,0);

  // These rules should have the same number of points
  libmesh_assert(gauss1D.n_points() == jac1D.n_points());
  
  // Save the number of points as a convenient variable
  const unsigned int n_points = gauss1D.n_points();

  // Resize the points and weights vectors
  _points.resize(n_points * n_points * n_points);
  _weights.resize(n_points * n_points * n_points);

  // Compute the conical product
  unsigned int qp = 0;
  for (unsigned int i=0; i<n_points; ++i)
    for (unsigned int j=0; j<n_points; ++j)
      for (unsigned int k=0; k<n_points; ++k, ++qp)
      {
        const Real xi=gauss1D.qp(i)(0);
        const Real yj=gauss1D.qp(j)(0);
        const Real zk=jac1D.qp(k)(0);
        
        _points[qp](0) = (1.-zk) * xi;
        _points[qp](1) = (1.-zk) * yj;
        _points[qp](2) = zk;
        _weights[qp]   = gauss1D.w(i) * gauss1D.w(j) * jac1D.w(k);
      }
  
  
}
 

void QConical::conical_product_tet (unsigned intp) [private]Implementation of conical product rule for a Tet in 3D of order = _order+2*p.

Definition at line 100 of file quadrature_conical.C.

References QBase::_points, QBase::_weights, QBase::get_dim(), QBase::n_points(), QBase::qp(), QBase::scale(), and QBase::w().

Referenced by init_3D().

{
  // Be sure the underlying rule object was built with the same dimension as the
  // rule we are about to construct.
  libmesh_assert (this->get_dim() == 3);
  
  QGauss  gauss1D(1,static_cast<Order>(_order+2*p));
  QJacobi jacA1D(1,static_cast<Order>(_order+2*p),1,0);
  QJacobi jacB1D(1,static_cast<Order>(_order+2*p),2,0);

  // The Gauss rule needs to be scaled to [0,1]
  std::pair<Real, Real> old_range(-1.0L, 1.0L);
  std::pair<Real, Real> new_range( 0.0L, 1.0L);
  gauss1D.scale(old_range,
                new_range);

  // Now construct the points and weights for the conical product rule.

  // All rules should have the same number of points
  libmesh_assert(gauss1D.n_points() == jacA1D.n_points());
  libmesh_assert(jacA1D.n_points()  == jacB1D.n_points());
  
  // Save the number of points as a convenient variable
  const unsigned int n_points = gauss1D.n_points();
  
  // All rules should be between x=0 and x=1
  libmesh_assert(gauss1D.qp(0)(0) >= 0.0); libmesh_assert(gauss1D.qp(n_points-1)(0) <= 1.0);
  libmesh_assert(jacA1D.qp(0)(0)  >= 0.0); libmesh_assert(jacA1D.qp(n_points-1)(0)  <= 1.0);
  libmesh_assert(jacB1D.qp(0)(0)  >= 0.0); libmesh_assert(jacB1D.qp(n_points-1)(0)  <= 1.0);

  // Resize the points and weights vectors
  _points.resize(n_points * n_points * n_points);
  _weights.resize(n_points * n_points * n_points);

  // Compute the conical product
  unsigned int gp = 0;
  for (unsigned int i=0; i<n_points; i++)
    for (unsigned int j=0; j<n_points; j++)
      for (unsigned int k=0; k<n_points; k++)
      {
        _points[gp](0) = jacB1D.qp(k)(0);                                                  //t[k];
        _points[gp](1) = jacA1D.qp(j)(0)  * (1.-jacB1D.qp(k)(0));                         //s[j]*(1.-t[k]);
        _points[gp](2) = gauss1D.qp(i)(0) * (1.-jacA1D.qp(j)(0)) * (1.-jacB1D.qp(k)(0)); //r[i]*(1.-s[j])*(1.-t[k]);
        _weights[gp]   = gauss1D.w(i)     * jacA1D.w(j)          * jacB1D.w(k);          //A[i]*B[j]*C[k];
        gp++;
      }
}
 

void QConical::conical_product_tri (unsigned intp) [private]Implementation of conical product rule for a Tri in 2D of order = _order+2*p.

Definition at line 51 of file quadrature_conical.C.

References QBase::_points, QBase::_weights, QBase::get_dim(), QBase::n_points(), QBase::qp(), QBase::scale(), and QBase::w().

Referenced by init_2D().

{
  // Be sure the underlying rule object was built with the same dimension as the
  // rule we are about to construct.
  libmesh_assert (this->get_dim() == 2);
  
  QGauss  gauss1D(1,static_cast<Order>(_order+2*p));
  QJacobi jac1D(1,static_cast<Order>(_order+2*p),1,0);
              
  // The Gauss rule needs to be scaled to [0,1]
  std::pair<Real, Real> old_range(-1.0L, 1.0L);
  std::pair<Real, Real> new_range( 0.0L, 1.0L);
  gauss1D.scale(old_range,
                new_range);

  // Now construct the points and weights for the conical product rule.

  // Both rules should have the same number of points.
  libmesh_assert(gauss1D.n_points() == jac1D.n_points());

  // Save the number of points as a convenient variable
  const unsigned int n_points = gauss1D.n_points();
  
  // Both rules should be between x=0 and x=1
  libmesh_assert(gauss1D.qp(0)(0) >= 0.0); libmesh_assert(gauss1D.qp(n_points-1)(0) <= 1.0);
  libmesh_assert(jac1D.qp(0)(0)   >= 0.0); libmesh_assert(jac1D.qp(n_points-1)(0) <= 1.0);

  // Resize the points and weights vectors
  _points.resize(n_points * n_points);
  _weights.resize(n_points * n_points);

  // Compute the conical product
  unsigned int gp = 0;
  for (unsigned int i=0; i<n_points; i++)
    for (unsigned int j=0; j<n_points; j++)
      {
        _points[gp](0) = jac1D.qp(j)(0);                          //s[j];
        _points[gp](1) = gauss1D.qp(i)(0) * (1.-jac1D.qp(j)(0)); //r[i]*(1.-s[j]);
        _weights[gp]   = gauss1D.w(i) * jac1D.w(j);              //A[i]*B[j];
        gp++;
      }
}
 

unsigned int QBase::get_dim () const [inline, inherited]Returns:

the dimension of the quadrature rule.

Definition at line 121 of file quadrature.h.

Referenced by InfFE< Dim, T_radial, T_map >::attach_quadrature_rule(), conical_product_pyramid(), conical_product_tet(), and conical_product_tri().

{ return _dim;  }
 

ElemType QBase::get_elem_type () const [inline, inherited]Returns:

the current element type we're set up for

Definition at line 103 of file quadrature.h.

    { return _type; }
 

std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.

Definition at line 45 of file reference_counter.C.

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

Referenced by ReferenceCounter::print_info().

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

  std::ostringstream out;
  
  out << '
      << ' ---------------------------------------------------------------------------- 
      << '| Reference count information                                                |
      << ' ---------------------------------------------------------------------------- ;
  
  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;

      out << '| ' << name << ' reference count information:
          << '|  Creations:    ' << creations    << '
          << '|  Destructions: ' << destructions << ';
    }
  
  out << ' ---------------------------------------------------------------------------- ;

  return out.str();

#else

  return '';
  
#endif
}
 

Order QBase::get_order () const [inline, inherited]Returns:

the order of the quadrature rule.

Definition at line 167 of file quadrature.h.

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

{ return static_cast<Order>(_order + _p_level); }
 

unsigned int QBase::get_p_level () const [inline, inherited]Returns:

the current p refinement level we're initialized with

Definition at line 109 of file quadrature.h.

    { return _p_level; }
 

const std::vector<Point>& QBase::get_points () const [inline, inherited]Returns:

a std::vector containing the quadrature point locations on a reference object.

Definition at line 127 of file quadrature.h.

References QBase::_points.

Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QMonomial::init_2D(), QGauss::init_2D(), QClough::init_2D(), QMonomial::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), InfFE< Dim, T_radial, T_map >::reinit(), and REINIT_ERROR().

{ return _points;  }
 

std::vector<Point>& QBase::get_points () [inline, inherited]Returns:

a std::vector containing the quadrature point locations on a reference object as a writeable reference.

Definition at line 133 of file quadrature.h.

References QBase::_points.

{ return _points;  }
 

const std::vector<Real>& QBase::get_weights () const [inline, inherited]Returns:

a std::vector containing the quadrature weights.

Definition at line 138 of file quadrature.h.

References QBase::_weights.

Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QMonomial::init_2D(), QGauss::init_2D(), QClough::init_2D(), QMonomial::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), and REINIT_ERROR().

{ return _weights; }
 

std::vector<Real>& QBase::get_weights () [inline, inherited]Returns:

a std::vector containing the quadrature weights.

Definition at line 143 of file quadrature.h.

References QBase::_weights.

{ return _weights; }
 

void 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 149 of file reference_counter.h.

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

Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().

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

  p.first++;
}
 

void 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 167 of file reference_counter.h.

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

Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().

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

  p.second++;
}
 

void QBase::init (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [inherited]Initializes the data structures to contain a quadrature rule for an object of type type.

Definition at line 26 of file quadrature.C.

References QBase::init_0D(), QBase::init_1D(), and QBase::init_2D().

Referenced by FE< Dim, T >::edge_reinit(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), QGauss::QGauss(), QJacobi::QJacobi(), QSimpson::QSimpson(), QTrap::QTrap(), InfFE< Dim, T_radial, T_map >::reinit(), and REINIT_ERROR().

{
  // check to see if we have already
  // done the work for this quadrature rule
  if (t == _type && p == _p_level)
    return;
  else
    {
      _type = t;
      _p_level = p;
    }
    
  
  
  switch(_dim)
    {
    case 0:
      this->init_0D(_type,_p_level);

      return;
      
    case 1:
      this->init_1D(_type,_p_level);

      return;
      
    case 2:
      this->init_2D(_type,_p_level);

      return;

    case 3:
      this->init_3D(_type,_p_level);

      return;

    default:
      libmesh_error();
    }
}
 

void QBase::init_0D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [protected, virtual, inherited]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 70 of file quadrature.C.

References QBase::_points, and QBase::_weights.

Referenced by QBase::init().

{
  _points.resize(1);
  _weights.resize(1);
  _points[0] = Point(0.);
  _weights[0] = 1.0;
}
 

void QConical::init_1D (const ElemType_type, unsignedp_level = 0) [inline, private, 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.

Implements QBase.

Definition at line 62 of file quadrature_conical.h.

  {
    // See about making this non-pure virtual in the base class
    libmesh_error();
  }
 

void QConical::init_2D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [private, virtual]The conical product rules are defined in 2D only for Tris.

Reimplemented from QBase.

Definition at line 26 of file quadrature_conical_2D.C.

References conical_product_tri(), libMeshEnums::TRI3, and libMeshEnums::TRI6.

{
  switch (_type)
    {
    case TRI3:
    case TRI6:
      {
        this->conical_product_tri(p);
        return;
        
      } // end case TRI3, TRI6


      
      //---------------------------------------------
      // Unsupported element type
    default:
      {
        std::cerr << 'ERROR: Unsupported element type: ' << _type << std::endl;
        libmesh_error();
      }
    } // end switch (_type)

  // We must have returned or errored-out by this point.  If not,
  // throw an error now.
  libmesh_error();
  return;
}
 

void QConical::init_3D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [private]The conical product rules are defined in 3D only for Tets.

Definition at line 26 of file quadrature_conical_3D.C.

References conical_product_pyramid(), conical_product_tet(), libMeshEnums::PYRAMID5, libMeshEnums::TET10, and libMeshEnums::TET4.

{
  switch (_type)
    {
    case TET4:
    case TET10:
      {
        this->conical_product_tet(p);
        return;
        
      } // end case TET4, TET10

    case PYRAMID5:
      {
        this->conical_product_pyramid(p);
        return;
        
      } // end case PYRAMID5

      
      //---------------------------------------------
      // Unsupported element type
    default:
      {
        std::cerr << 'ERROR: Unsupported element type: ' << _type << std::endl;
        libmesh_error();
      }
    } // end switch (_type)

  // We must have returned or errored-out by this point.  If not,
  // throw an error now.
  libmesh_error();
  return;
}
 

static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 76 of file reference_counter.h.

References ReferenceCounter::_n_objects.

Referenced by System::read_serialized_blocked_dof_objects(), and System::write_serialized_blocked_dof_objects().

  { return _n_objects; }
 

unsigned int QBase::n_points () const [inline, inherited]Returns:

the number of points associated with the quadrature rule.

Definition at line 115 of file quadrature.h.

References QBase::_points.

Referenced by FEBase::coarsened_dof_values(), conical_product_pyramid(), conical_product_tet(), conical_product_tri(), FEMSystem::eulerian_residual(), InfFE< Dim, T_radial, T_map >::init_face_shape_functions(), FEMSystem::mass_residual(), and QBase::print_info().

    { libmesh_assert (!_points.empty()); return _points.size(); }
 

void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.

Definition at line 83 of file reference_counter.C.

References ReferenceCounter::get_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
  
  std::cout << ReferenceCounter::get_info();
  
#endif
}
 

void QBase::print_info (std::ostream &os = std::cout) const [inline, inherited]Prints information relevant to the quadrature rule.

Definition at line 350 of file quadrature.h.

References QBase::_points, QBase::_weights, QBase::n_points(), and QBase::qp().

Referenced by operator<<().

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

  os << 'N_Q_Points=' << this->n_points() << std::endl << std::endl;
  for (unsigned int qp=0; qp<this->n_points(); qp++)
    {
      os << ' Point ' << qp << ':
         << '  '
         << _points[qp]
         << ' Weight:'
         << '  w=' << _weights[qp] << ' << std::endl;
    }
}
 

Point QBase::qp (const unsigned inti) const [inline, inherited]Returns:

the $ i^{th} $ quadrature point on the reference object.

Definition at line 148 of file quadrature.h.

References QBase::_points.

Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), and QBase::print_info().

    { libmesh_assert (i < _points.size()); return _points[i]; }
 

void QBase::scale (std::pair< Real, Real >old_range, std::pair< Real, Real >new_range) [inherited]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 81 of file quadrature.C.

References QBase::_points, and QBase::_weights.

Referenced by conical_product_tet(), and conical_product_tri().

{
  // Make sure we are in 1D
  libmesh_assert(_dim == 1);
  
  // Make sure that we have sane ranges 
  libmesh_assert(new_range.second > new_range.first);
  libmesh_assert(old_range.second > old_range.first);

  // Make sure there are some points
  libmesh_assert(_points.size() > 0);

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

  // Compute the scale factor and scale the weights
  const Real scfact = (new_range.second - new_range.first) /
                      (old_range.second - old_range.first);

  for (unsigned int i=0; i<_points.size(); i++)
    _weights[i] *= scfact;
}
 

QuadratureType QConical::type () const [inline, virtual]Returns:

the QuadratureType for this class

Implements QBase.

Definition at line 58 of file quadrature_conical.h.

References libMeshEnums::QCONICAL.

{ return QCONICAL; }
 

Real QBase::w (const unsigned inti) const [inline, inherited]Returns:

the $ i^{th} $ quadrature weight.

Definition at line 154 of file quadrature.h.

References QBase::_weights.

Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), and QMonomial::stroud_rule().

    { libmesh_assert (i < _weights.size()); return _weights[i]; }
 

Friends And Related Function Documentation

 

std::ostream& operator<< (std::ostream &os, const QBase &q) [friend, inherited]Same as above, but allows you to use the stream syntax.

Definition at line 196 of file quadrature.C.

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

Member Data Documentation

 

ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.

Definition at line 110 of file reference_counter.h.

Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().  

Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.

Definition at line 123 of file reference_counter.h.  

Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.

Definition at line 118 of file reference_counter.h.

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

std::cerr<< 'ERROR: Seems as if this quadrature rule' << std::endl << ' is not implemented for 2D.' << std::endl; libmesh_error(); }#endif virtual void init_3D (const ElemType, unsigned int =0)#ifndef DEBUG {}#else { std::cerr << 'ERROR: Seems as if this quadrature rule' << std::endl << ' is not implemented for 3D.' << std::endl; libmesh_error(); }#endif void tensor_product_quad (const QBase& q1D); void tensor_product_hex (const QBase& q1D); void tensor_product_prism (const QBase& q1D, const QBase& q2D); const unsigned int _dim; const Order _order; ElemType _type; unsigned int _p_level; std::vector<Point> QBase::_points [protected, inherited]

Definition at line 320 of file quadrature.h.

Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), QGauss::dunavant_rule(), QGauss::dunavant_rule2(), QBase::get_points(), QGrundmann_Moller::gm_rule(), QBase::init_0D(), QTrap::init_1D(), QSimpson::init_1D(), QJacobi::init_1D(), QGrid::init_1D(), QGauss::init_1D(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), QBase::n_points(), QBase::print_info(), QBase::qp(), QBase::scale(), QMonomial::stroud_rule(), and QMonomial::wissmann_rule().  

std::vector<Real> QBase::_weights [protected, inherited]The value of the quadrature weights.

Definition at line 325 of file quadrature.h.

Referenced by conical_product_pyramid(), conical_product_tet(), conical_product_tri(), QGauss::dunavant_rule(), QGauss::dunavant_rule2(), QBase::get_weights(), QGrundmann_Moller::gm_rule(), QBase::init_0D(), QTrap::init_1D(), QSimpson::init_1D(), QJacobi::init_1D(), QGrid::init_1D(), QGauss::init_1D(), QClough::init_1D(), QTrap::init_2D(), QSimpson::init_2D(), QMonomial::init_2D(), QGrid::init_2D(), QGauss::init_2D(), QClough::init_2D(), QTrap::init_3D(), QSimpson::init_3D(), QMonomial::init_3D(), QGrid::init_3D(), QGauss::init_3D(), QGauss::keast_rule(), QMonomial::kim_rule(), QBase::print_info(), QBase::scale(), QMonomial::stroud_rule(), QBase::w(), and QMonomial::wissmann_rule().  

bool QBase::allow_rules_with_negative_weights [inherited]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 203 of file quadrature.h.

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

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Public Member Functions
Static Public Member Functions
Public Attributes
Protected Types
Protected Member Functions
Protected Attributes
Static Protected Attributes
Private Member Functions
Friends
Detailed Description
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.
Constructor & Destructor Documentation
QConical::QConical (const unsigned int_dim, const Order_order = INVALID_ORDER)Constructor. Declares the order of the quadrature rule.
QConical::~QConical ()Destructor.
Member Function Documentation
AutoPtr< QBase > QBase::build (const std::string &name, const unsigned int_dim, const Order_order = INVALID_ORDER) [static, inherited]Builds a specific quadrature rule, identified through the name string. An AutoPtr<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.
AutoPtr< QBase > QBase::build (const QuadratureType_qt, const unsigned int_dim, const Order_order = INVALID_ORDER) [static, inherited]Builds a specific quadrature rule, identified through the QuadratureType. An AutoPtr<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.
void QConical::conical_product_pyramid (unsigned intp) [private]Implementation of conical product rule for a Pyramid in 3D of order = _order+2*p.
void QConical::conical_product_tet (unsigned intp) [private]Implementation of conical product rule for a Tet in 3D of order = _order+2*p.
void QConical::conical_product_tri (unsigned intp) [private]Implementation of conical product rule for a Tri in 2D of order = _order+2*p.
unsigned int QBase::get_dim () const [inline, inherited]Returns:
ElemType QBase::get_elem_type () const [inline, inherited]Returns:
std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.
Order QBase::get_order () const [inline, inherited]Returns:
unsigned int QBase::get_p_level () const [inline, inherited]Returns:
const std::vector<Point>& QBase::get_points () const [inline, inherited]Returns:
std::vector<Point>& QBase::get_points () [inline, inherited]Returns:
const std::vector<Real>& QBase::get_weights () const [inline, inherited]Returns:
std::vector<Real>& QBase::get_weights () [inline, inherited]Returns:
void 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.
void 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.
void QBase::init (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [inherited]Initializes the data structures to contain a quadrature rule for an object of type type.
void QBase::init_0D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [protected, virtual, inherited]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.
void QConical::init_1D (const ElemType_type, unsignedp_level = 0) [inline, private, 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.
void QConical::init_2D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [private, virtual]The conical product rules are defined in 2D only for Tris.
void QConical::init_3D (const ElemType_type = INVALID_ELEM, unsigned intp_level = 0) [private]The conical product rules are defined in 3D only for Tets.
static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.
unsigned int QBase::n_points () const [inline, inherited]Returns:
void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.
void QBase::print_info (std::ostream &os = std::cout) const [inline, inherited]Prints information relevant to the quadrature rule.
Point QBase::qp (const unsigned inti) const [inline, inherited]Returns:
void QBase::scale (std::pair< Real, Real >old_range, std::pair< Real, Real >new_range) [inherited]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.
QuadratureType QConical::type () const [inline, virtual]Returns:
Real QBase::w (const unsigned inti) const [inline, inherited]Returns:
Friends And Related Function Documentation
std::ostream& operator<< (std::ostream &os, const QBase &q) [friend, inherited]Same as above, but allows you to use the stream syntax.
Member Data Documentation
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.
std::cerr<< 'ERROR: Seems as if this quadrature rule' << std::endl << ' is not implemented for 2D.' << std::endl; libmesh_error(); }#endif virtual void init_3D (const ElemType, unsigned int =0)#ifndef DEBUG {}#else { std::cerr << 'ERROR: Seems as if this quadrature rule' << std::endl << ' is not implemented for 3D.' << std::endl; libmesh_error(); }#endif void tensor_product_quad (const QBase& q1D); void tensor_product_hex (const QBase& q1D); void tensor_product_prism (const QBase& q1D, const QBase& q2D); const unsigned int _dim; const Order _order; ElemType _type; unsigned int _p_level; std::vector<Point> QBase::_points [protected, inherited]
std::vector<Real> QBase::_weights [protected, inherited]The value of the quadrature weights.
bool QBase::allow_rules_with_negative_weights [inherited]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.
Author

This document was created by man2html, using the manual pages.
Time: 21:53:29 GMT, April 16, 2011