libMesh::QGrundmann_Moller Class Reference

#include <quadrature_gm.h>

Inheritance diagram for libMesh::QGrundmann_Moller:

Public Member Functions

 QGrundmann_Moller (const unsigned int _dim, const Order _order=INVALID_ORDER)
 
 ~QGrundmann_Moller ()
 
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=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 AutoPtr< QBasebuild (const std::string &name, const unsigned int _dim, const Order _order=INVALID_ORDER)
 
static AutoPtr< QBasebuild (const QuadratureType _qt, const unsigned int _dim, const Order _order=INVALID_ORDER)
 
static void print_info (std::ostream &out=libMesh::out)
 
static std::string get_info ()
 
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

virtual void init_0D (const ElemType type=INVALID_ELEM, unsigned int p_level=0)
 
virtual void init_2D (const ElemType, unsigned int=0)
 
void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

libMesh::err<< "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{libMesh::err<< "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
 
static bool _enable_print_counter = true
 

Private Member Functions

void init_1D (const ElemType, unsigned int=0)
 
void init_3D (const ElemType _type=INVALID_ELEM, unsigned int p_level=0)
 
void gm_rule (unsigned int s)
 
void compose_all (unsigned int s, unsigned int p, std::vector< std::vector< unsigned int > > &result)
 

Detailed Description

This class implements the Grundmann-Moller quadrature rules for tetrahedra. The GM rules are well-defined for simplices of arbitrary dimension and to any order, but the rules by Dunavant for two-dimensional simplices are in general superior. This is primarily due to the fact that the GM rules contain a significant proportion of negative weights, making them susceptible to round-off error at high-order.

The GM rules are interesting in 3D because they overlap with the conical product rules at higher order while having significantly fewer evaluation points, making them potentially much more efficient. The table below gives a comparison between the number of points in a conical product (CP) rule and the GM rule of equivalent order. The GM rules are defined to be exact for polynomials of degree d=2*s+1, s=0,1,2,3,... The table also gives the percentage of each GM rule's weights which are negative. Although the percentage of negative weights does not grow particularly quickly, the amplification factor (a measure of the effect of round-off) defined as

amp. factor = $ \frac{1}{V} \sum \|w_i\|, $

where V is the volume of the reference element, does grow quickly. (A rule with all positive has has an amplification factor of 1.0 by definition.)

  s  | d     | N. CP        | N. GM   | % neg wts | amp. factor
-----------------------------------------------------------------
  0  | 1     |              | 1       |           |
  1  | 2-3   |              | 5       |           |
  2  | 4-5   |              | 15      |           |
  3  | 6-7   |              | 35      | 31.43     |   11.94
  4  | 8-9   |  5^3=125     | 70      | 34.29     |   25.35
  5  | 10-11 |  6^3=216     | 126     | 36.51     |   54.14
  6  | 12-13 |  7^3=343     | 210     | 38.10     |  116.30
  7  | 14-15 |  8^3=512     | 330     | 39.39     |  251.10
  8  | 16-17 |  9^3=729     | 495     | 40.40     |  544.68
  9  | 18-19 | 10^3=1,000   | 715     | 41.26     | 1186.16
 10  | 20-21 | 11^3=1,331   | 1,001   | 41.96     | 2591.97
 11  | 22-23 | 12^3=1,728   | 1,365   | 42.56     | 5680.75
 ...
 16  | 32-33 | 17^3=4,913   | 4,845   |
 17  | 34-35 | 18^3=5,832   | 5,985   | <= Cross-over point, CP has fewer points for d >= 34
 18  | 36-37 | 19^3=6,859   | 7,315   |
 ...
 21  | 42-43 | 22^3=10,648  | 12,650  |

Reference: Axel Grundmann and Michael M"{o}ller, "Invariant Integration Formulas for the N-Simplex by Combinatorial Methods," SIAM Journal on Numerical Analysis, Volume 15, Number 2, April 1978, pages 282-290.

Reference LGPL Fortran90 code by John Burkardt can be found here: http://people.scs.fsu.edu/~burkardt/f_src/gm_rules/gm_rules.html

Author
John W. Peterson, 2008

Definition at line 100 of file quadrature_gm.h.

Member Typedef Documentation

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

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::QGrundmann_Moller::QGrundmann_Moller ( const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
)

Constructor. Declares the order of the quadrature rule.

Definition at line 32 of file quadrature_gm.C.

33  : QBase(d,o)
34 {
35 }
libMesh::QGrundmann_Moller::~QGrundmann_Moller ( )

Destructor.

Definition at line 39 of file quadrature_gm.C.

40 {
41 }

Member Function Documentation

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

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 37 of file quadrature_build.C.

Referenced by libMesh::InfFE< friend_Dim, friend_T_radial, friend_T_map >::attach_quadrature_rule().

40 {
41  return QBase::build (Utility::string_to_enum<QuadratureType> (type),
42  _dim,
43  _order);
44 }
AutoPtr< QBase > libMesh::QBase::build ( const QuadratureType  _qt,
const unsigned int  _dim,
const Order  _order = INVALID_ORDER 
)
staticinherited

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 48 of file quadrature_build.C.

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

51 {
52  switch (_qt)
53  {
54 
55  case QCLOUGH:
56  {
57 #ifdef DEBUG
58  if (_order > TWENTYTHIRD)
59  {
60  libMesh::out << "WARNING: Clough quadrature implemented" << std::endl
61  << " up to TWENTYTHIRD order." << std::endl;
62  }
63 #endif
64 
65  AutoPtr<QBase> ap(new QClough(_dim, _order));
66  return ap;
67  }
68 
69  case QGAUSS:
70  {
71 
72 #ifdef DEBUG
73  if (_order > FORTYTHIRD)
74  {
75  libMesh::out << "WARNING: Gauss quadrature implemented" << std::endl
76  << " up to FORTYTHIRD order." << std::endl;
77  }
78 #endif
79 
80  AutoPtr<QBase> ap(new QGauss(_dim, _order));
81  return ap;
82  }
83 
84  case QJACOBI_1_0:
85  {
86 
87 #ifdef DEBUG
88  if (_order > TWENTYTHIRD)
89  {
90  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
91  << " up to TWENTYTHIRD order." << std::endl;
92  }
93 
94  if (_dim > 1)
95  {
96  libMesh::out << "WARNING: Jacobi(1,0) quadrature implemented" << std::endl
97  << " in 1D only." << std::endl;
98  }
99 #endif
100 
101  AutoPtr<QBase> ap(new QJacobi(_dim, _order, 1, 0));
102  return ap;
103  }
104 
105  case QJACOBI_2_0:
106  {
107 
108 #ifdef DEBUG
109  if (_order > TWENTYTHIRD)
110  {
111  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
112  << " up to TWENTYTHIRD order." << std::endl;
113  }
114 
115  if (_dim > 1)
116  {
117  libMesh::out << "WARNING: Jacobi(2,0) quadrature implemented" << std::endl
118  << " in 1D only." << std::endl;
119  }
120 #endif
121 
122  AutoPtr<QBase> ap(new QJacobi(_dim, _order, 2, 0));
123  return ap;
124  }
125 
126  case QSIMPSON:
127  {
128 
129 #ifdef DEBUG
130  if (_order > THIRD)
131  {
132  libMesh::out << "WARNING: Simpson rule provides only" << std::endl
133  << " THIRD order!" << std::endl;
134  }
135 #endif
136 
137  AutoPtr<QBase> ap(new QSimpson(_dim));
138  return ap;
139  }
140 
141  case QTRAP:
142  {
143 
144 #ifdef DEBUG
145  if (_order > FIRST)
146  {
147  libMesh::out << "WARNING: Trapezoidal rule provides only" << std::endl
148  << " FIRST order!" << std::endl;
149  }
150 #endif
151 
152  AutoPtr<QBase> ap(new QTrap(_dim));
153  return ap;
154  }
155 
156 
157  default:
158  {
159  libMesh::err << "ERROR: Bad qt=" << _qt << std::endl;
160  libmesh_error();
161  }
162  }
163 
164 
165  libmesh_error();
166  AutoPtr<QBase> ap(NULL);
167  return ap;
168 }
void libMesh::QGrundmann_Moller::compose_all ( unsigned int  s,
unsigned int  p,
std::vector< std::vector< unsigned int > > &  result 
)
private

Routine which generates p-compositions of a given order, s, as well as permutations thereof. This routine is called internally by the gm_rule() routine, you should not call this yourself!

Definition at line 144 of file quadrature_gm.C.

References libMesh::libmesh_assert_greater().

Referenced by gm_rule().

147 {
148  // Clear out results remaining from previous calls
149  result.clear();
150 
151  // Allocate storage for a workspace. The workspace will periodically
152  // be copied into the result container.
153  std::vector<unsigned int> workspace(p);
154 
155  // The first result is always (s,0,...,0)
156  workspace[0] = s;
157  result.push_back(workspace);
158 
159  // the value of the first non-zero entry
160  unsigned int head_value=s;
161 
162  // When head_index=-1, it refers to "off the front" of the array. Therefore,
163  // this needs to be a regular int rather than unsigned. I initially tried to
164  // do this with head_index unsigned and an else statement below, but then there
165  // is the special case: (1,0,...,0) which does not work correctly.
166  int head_index = -1;
167 
168  // At the end, all the entries will be in the final slot of workspace
169  while (workspace.back() != s)
170  {
171  // Uncomment for debugging
172  //libMesh::out << "previous head_value=" << head_value << " -> ";
173 
174  // If the previous head value is still larger than 1, reset the index
175  // to "off the front" of the array
176  if (head_value > 1)
177  head_index = -1;
178 
179  // Either move the index onto the front of the array or on to
180  // the next value.
181  head_index++;
182 
183  // Get current value of the head entry
184  head_value = workspace[head_index];
185 
186  // Uncomment for debugging
187  //std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(libMesh::out, " "));
188  //libMesh::out << ", head_index=" << head_index;
189  //libMesh::out << ", head_value=" << head_value << " -> ";
190 
191  // Put a zero into the head_index of the array. If head_index==0,
192  // this will be overwritten in the next line with head_value-1.
193  workspace[head_index] = 0;
194 
195  // The initial entry gets the current head value, minus 1.
196  // If head_value > 1, the next loop iteration will start back
197  // at workspace[0] again.
198  libmesh_assert_greater (head_value, 0);
199  workspace[0] = head_value - 1;
200 
201  // Increment the head+1 value
202  workspace[head_index+1] += 1;
203 
204  // Save this composition in the results
205  result.push_back(workspace);
206 
207  // Uncomment for debugging
208  //std::copy(workspace.begin(), workspace.end(), std::ostream_iterator<int>(libMesh::out, " "));
209  //libMesh::out<<"\n";
210  }
211 }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

107 {
108  _enable_print_counter = false;
109  return;
110 }
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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.

101 {
102  _enable_print_counter = true;
103  return;
104 }
unsigned int libMesh::QBase::get_dim ( ) const
inlineinherited
ElemType libMesh::QBase::get_elem_type ( ) const
inlineinherited
Returns
the current element type we're set up for

Definition at line 104 of file quadrature.h.

105  { return _type; }
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

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().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
Order libMesh::QBase::get_order ( ) const
inlineinherited
Returns
the order of the quadrature rule.

Definition at line 169 of file quadrature.h.

Referenced by libMesh::InfFE< friend_Dim, friend_T_radial, friend_T_map >::attach_quadrature_rule().

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

Definition at line 110 of file quadrature.h.

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

Definition at line 129 of file quadrature.h.

References libMesh::QBase::_points.

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

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

Definition at line 135 of file quadrature.h.

References libMesh::QBase::_points.

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

Definition at line 140 of file quadrature.h.

References libMesh::QBase::_weights.

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

140 { return _weights; }
std::vector<Real>& libMesh::QBase::get_weights ( )
inlineinherited
Returns
a std::vector containing the quadrature weights.

Definition at line 145 of file quadrature.h.

References libMesh::QBase::_weights.

145 { return _weights; }
void libMesh::QGrundmann_Moller::gm_rule ( unsigned int  s)
private

This routine is called from the different cases of init_3D(). It actually fills the _points and _weights vectors for a given rule index, s.

Definition at line 47 of file quadrature_gm.C.

References libMesh::QBase::_points, libMesh::QBase::_weights, compose_all(), libMesh::dim, std::max(), libMesh::Real, and libMesh::MeshTools::weight().

Referenced by init_3D().

48 {
49  // A GM rule of index s can integrate polynomials of degree 2*s+1 exactly
50  const unsigned int degree = 2*s+1;
51 
52  // Here we are considering only tetrahedra rules, so dim==3
53  const unsigned int dim = 3;
54 
55  // The number of points for rule of index s is
56  // (dim+1+s)! / (dim+1)! / s!
57  // In 3D, this is = 1/24 * P_{i=1}^4 (s+i)
58  const unsigned int n_pts = (s+4)*(s+3)*(s+2)*(s+1) / 24;
59  //libMesh::out << "n_pts=" << n_pts << std::endl;
60 
61  // Allocate space for points and weights
62  _points.resize(n_pts);
63  _weights.resize(n_pts);
64 
65  // (-1)^i -> This one flips sign at each iteration of the i-loop below.
66  int one_pm=1;
67 
68  // Where we store all the integer point compositions/permutations
69  std::vector<std::vector<unsigned int> > permutations;
70 
71  // Index into the vector where we should start adding the next round of points/weights
72  std::size_t offset=0;
73 
74  // Implement the GM formula 4.1 on page 286 of the paper
75  for (unsigned int i=0; i<=s; ++i)
76  {
77  // Get all the ordered compositions (and their permutations)
78  // of |beta| = s-i into dim+1=4 parts
79  compose_all(s-i, dim+1, permutations);
80  //libMesh::out << "n. permutations=" << permutations.size() << std::endl;
81 
82  for (unsigned int p=0; p<permutations.size(); ++p)
83  {
84  // We use the first dim=3 entries of each permutation to
85  // construct an integration point.
86  for (unsigned int j=0; j<3; ++j)
87  _points[offset+p](j) =
88  static_cast<Real>(2.*permutations[p][j] + 1.) /
89  static_cast<Real>( degree + dim - 2.*i );
90  }
91 
92  // Compute the weight for this i, being careful to avoid overflow.
93  // This technique is borrowed from Burkardt's code as well.
94  // Use once for each of the points obtained from the permutations array.
95  Real weight = one_pm;
96 
97  // This for loop needs to run for dim, degree, or dim+degree-i iterations,
98  // whichever is largest.
99  const unsigned int weight_loop_index =
100  std::max(dim, std::max(degree, degree+dim-i));
101 
102  for (unsigned int j=1; j<=weight_loop_index; ++j)
103  {
104  if (j <= degree) // Accumulate (d+n-2i)^d term
105  weight *= static_cast<Real>(degree+dim-2*i);
106 
107  if (j <= 2*s) // Accumulate 2^{-2s}
108  weight *= 0.5;
109 
110  if (j <= i) // Accumulate (i!)^{-1}
111  weight /= static_cast<Real>(j);
112 
113  if (j <= degree+dim-i) // Accumulate ( (d+n-i)! )^{-1}
114  weight /= static_cast<Real>(j);
115  }
116 
117  // This is the weight for each of the points computed previously
118  for (unsigned int j=0; j<permutations.size(); ++j)
119  _weights[offset+j] = weight;
120 
121  // Change sign for next iteration
122  one_pm = -one_pm;
123 
124  // Update offset for the next set of points
125  offset += permutations.size();
126  }
127 }
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

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().

164 {
165  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
166  std::pair<unsigned int, unsigned int>& p = _counts[name];
167 
168  p.first++;
169 }
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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().

177 {
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  std::pair<unsigned int, unsigned int>& p = _counts[name];
180 
181  p.second++;
182 }
void libMesh::QBase::init ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
inherited

Initializes the data structures to contain a quadrature rule for an object of type type.

Definition at line 27 of file quadrature.C.

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

Referenced by libMesh::QClough::init_1D(), libMesh::QTrap::init_2D(), libMesh::QClough::init_2D(), libMesh::QGrid::init_2D(), libMesh::QGauss::init_2D(), libMesh::QSimpson::init_2D(), libMesh::QMonomial::init_2D(), libMesh::QTrap::init_3D(), libMesh::QGrid::init_3D(), libMesh::QGauss::init_3D(), libMesh::QSimpson::init_3D(), libMesh::QMonomial::init_3D(), libMesh::QGauss::QGauss(), libMesh::QJacobi::QJacobi(), libMesh::QSimpson::QSimpson(), and libMesh::QTrap::QTrap().

29 {
30  // check to see if we have already
31  // done the work for this quadrature rule
32  if (t == _type && p == _p_level)
33  return;
34  else
35  {
36  _type = t;
37  _p_level = p;
38  }
39 
40 
41 
42  switch(_dim)
43  {
44  case 0:
45  this->init_0D(_type,_p_level);
46 
47  return;
48 
49  case 1:
50  this->init_1D(_type,_p_level);
51 
52  return;
53 
54  case 2:
55  this->init_2D(_type,_p_level);
56 
57  return;
58 
59  case 3:
60  this->init_3D(_type,_p_level);
61 
62  return;
63 
64  default:
65  libmesh_error();
66  }
67 }
void libMesh::QBase::init_0D ( const ElemType  type = INVALID_ELEM,
unsigned int  p_level = 0 
)
protectedvirtualinherited

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 71 of file quadrature.C.

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

Referenced by libMesh::QBase::init().

73 {
74  _points.resize(1);
75  _weights.resize(1);
76  _points[0] = Point(0.);
77  _weights[0] = 1.0;
78 }
void libMesh::QGrundmann_Moller::init_1D ( const ElemType  type,
unsigned  p_level = 0 
)
inlineprivatevirtual

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 libMesh::QBase.

Definition at line 123 of file quadrature_gm.h.

125  {
126  // See about making this non-pure virtual in the base class
127  libmesh_error();
128  }
virtual void libMesh::QBase::init_2D ( const ElemType  ,
unsigned int  = 0 
)
inlineprotectedvirtualinherited

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::QMonomial, libMesh::QConical, libMesh::QSimpson, libMesh::QGauss, libMesh::QClough, libMesh::QGrid, and libMesh::QTrap.

Definition at line 246 of file quadrature.h.

Referenced by libMesh::QBase::init().

249  {}
void libMesh::QGrundmann_Moller::init_3D ( const ElemType  _type = INVALID_ELEM,
unsigned int  p_level = 0 
)
private

The GM rules are only defined for 3D since better 2D rules for simplexes are available.

Definition at line 28 of file quadrature_gm_3D.C.

References libMesh::QBase::allow_rules_with_negative_weights, libMesh::err, gm_rule(), libMeshEnums::TET10, and libMeshEnums::TET4.

30 {
31  // Nearly all GM rules contain negative weights, so if you are not
32  // allowing rules with negative weights, we cannot continue!
34  {
35  libMesh::err << "You requested a Grundmann-Moller rule but\n"
36  << "are not allowing rules with negative weights!\n"
37  << "Either select a different quadrature class or\n"
38  << "set allow_rules_with_negative_weights==true."
39  << std::endl;
40 
41  libmesh_error();
42  }
43 
44  switch (type_in)
45  {
46  case TET4:
47  case TET10:
48  {
49  // Untested above _order=23 but should work...
50  gm_rule( (_order + 2*p)/2 );
51  return;
52 
53  } // end case TET4, TET10
54 
55 
56 
57  //---------------------------------------------
58  // Unsupported element type
59  default:
60  {
61  libMesh::err << "ERROR: Unsupported element type: " << type_in << std::endl;
62  libmesh_error();
63  }
64  } // end switch (type_in)
65 
66  // We must have returned or errored-out by this point. If not,
67  // throw an error now.
68  libmesh_error();
69  return;
70 }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

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.

80  { return _n_objects; }
unsigned int libMesh::QBase::n_points ( ) const
inlineinherited
Returns
the number of points associated with the quadrature rule.

Definition at line 116 of file quadrature.h.

References libMesh::QBase::_points, and libMesh::libmesh_assert().

Referenced by libMesh::QConical::conical_product_pyramid(), libMesh::QConical::conical_product_tet(), libMesh::QConical::conical_product_tri(), libMesh::ProjectFEMSolution::operator()(), and libMesh::QBase::print_info().

117  { libmesh_assert (!_points.empty());
118  return libmesh_cast_int<unsigned int>(_points.size()); }
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

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().

89 {
91 }
void libMesh::QBase::print_info ( std::ostream &  os = libMesh::out) const
inlineinherited

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

Definition at line 362 of file quadrature.h.

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

Referenced by libMesh::operator<<().

363 {
364  libmesh_assert(!_points.empty());
365  libmesh_assert(!_weights.empty());
366 
367  os << "N_Q_Points=" << this->n_points() << std::endl << std::endl;
368  for (unsigned int qpoint=0; qpoint<this->n_points(); qpoint++)
369  {
370  os << " Point " << qpoint << ":\n"
371  << " "
372  << _points[qpoint]
373  << " Weight:\n "
374  << " w=" << _weights[qpoint] << "\n" << std::endl;
375  }
376 }
Point libMesh::QBase::qp ( const unsigned int  i) const
inlineinherited
Returns
the $ i^{th} $ quadrature point on the reference object.

Definition at line 150 of file quadrature.h.

References libMesh::QBase::_points.

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

151  { 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 
)
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 82 of file quadrature.C.

References libMesh::QBase::_points, libMesh::QBase::_weights, libMesh::libmesh_assert_greater(), and libMesh::Real.

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

84 {
85  // Make sure we are in 1D
86  libmesh_assert_equal_to (_dim, 1);
87 
88  // Make sure that we have sane ranges
89  libmesh_assert_greater (new_range.second, new_range.first);
90  libmesh_assert_greater (old_range.second, old_range.first);
91 
92  // Make sure there are some points
93  libmesh_assert_greater (_points.size(), 0);
94 
95  // We're mapping from old_range -> new_range
96  for (unsigned int i=0; i<_points.size(); i++)
97  {
98  _points[i](0) =
99  (_points[i](0) - old_range.first) *
100  (new_range.second - new_range.first) /
101  (old_range.second - old_range.first) +
102  new_range.first;
103  }
104 
105  // Compute the scale factor and scale the weights
106  const Real scfact = (new_range.second - new_range.first) /
107  (old_range.second - old_range.first);
108 
109  for (unsigned int i=0; i<_points.size(); i++)
110  _weights[i] *= scfact;
111 }
virtual bool libMesh::QBase::shapes_need_reinit ( )
inlinevirtualinherited

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 198 of file quadrature.h.

198 { return false; }
QuadratureType libMesh::QGrundmann_Moller::type ( ) const
inlinevirtual
Returns
QGRUNDMANN_MOLLER

Implements libMesh::QBase.

Definition at line 118 of file quadrature_gm.h.

References libMeshEnums::QGRUNDMANN_MOLLER.

118 { return QGRUNDMANN_MOLLER; }
Real libMesh::QBase::w ( const unsigned int  i) const
inlineinherited
Returns
the $ i^{th} $ quadrature weight.

Definition at line 156 of file quadrature.h.

References libMesh::QBase::_weights.

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

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

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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().

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

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
staticprotectedinherited

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().

libMesh::err<< "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 { libMesh::err << "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> libMesh::QBase::_points
protectedinherited
bool libMesh::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 215 of file quadrature.h.

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


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:58:01 UTC

Hosted By:
SourceForge.net Logo