This file is indexed.

/usr/include/libmesh/quadrature.h is in libmesh-dev 0.7.1-2ubuntu1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
// $Id: quadrature.h 4106 2010-11-16 16:35:07Z friedmud $

// The libMesh Finite Element Library.
// Copyright (C) 2002-2008 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
  
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
  
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
// Lesser General Public License for more details.
  
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA



#ifndef __quadrature_h__
#define __quadrature_h__

// C++ includes
#include <vector>
#include <string>
#include <utility>

// Local includes
#include "libmesh_common.h"
#include "reference_counted_object.h"
#include "point.h"
#include "enum_elem_type.h"
#include "enum_order.h"
#include "enum_quadrature_type.h"
#include "auto_ptr.h"

namespace libMesh
{



/**
 * This is the \p QBase class.  It provides the basic functionality
 * from which various quadrature rules can be derived.  The class contains
 * \p dim dimensional points describing the quadrature locations
 * (referenced to a master object) and associated weights.
 *
 * @author Benjamin S. Kirk, 2002
 */

// ------------------------------------------------------------
// QBase class definition

class QBase : public ReferenceCountedObject<QBase>
{
protected:

  /**
   * Constructor. Protected to prevent instantiation
   * of this base class.
   */
  QBase (const unsigned int _dim,
	 const Order _order=INVALID_ORDER);
  
public:
  
  /**
   * Destructor.
   */
  virtual ~QBase() {}

  /**
   * @returns the quadrature type in derived classes.
   */
  virtual QuadratureType type() const = 0;

  /**
   * Builds a specific quadrature rule, identified through the
   * \p name string.  An \p 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 \p name
   * must be mappable through the \p Utility::string_to_enum<>() 
   * function.
   */
  static AutoPtr<QBase> build (const std::string &name,
			       const unsigned int _dim,
			       const Order _order=INVALID_ORDER);

  /**
   * Builds a specific quadrature rule, identified through the
   * \p QuadratureType.  An \p 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.
   */
  static AutoPtr<QBase> build (const QuadratureType _qt,
			       const unsigned int _dim,
			       const Order _order=INVALID_ORDER);

  /**
   * @returns the current element type we're set up for
   */    
  ElemType get_elem_type() const
    { return _type; }

  /**
   * @returns the current p refinement level we're initialized with
   */    
  unsigned int get_p_level() const
    { return _p_level; }

  /**
   * @returns the number of points associated with the quadrature rule.
   */    
  unsigned int n_points() const
    { libmesh_assert (!_points.empty()); return _points.size(); }

  /**
   * @returns the dimension of the quadrature rule.
   */
  unsigned int get_dim() const { return _dim;  }

  /**
   * @returns a \p std::vector containing the quadrature point locations
   * on a reference object.
   */
  const std::vector<Point>& get_points() const { return _points;  }

  /**
   * @returns a \p std::vector containing the quadrature point locations
   * on a reference object as a writeable reference.
   */
  std::vector<Point>& get_points() { return _points;  }

  /**
   * @returns a \p std::vector containing the quadrature weights.
   */
  const std::vector<Real>& get_weights() const { return _weights; }

  /**
   * @returns a \p std::vector containing the quadrature weights.
   */
  std::vector<Real>& get_weights() { return _weights; }

  /**
   * @returns the \f$ i^{th} \f$ quadrature point on the reference object.
   */
  Point qp(const unsigned int i) const
    { libmesh_assert (i < _points.size()); return _points[i]; }

  /**
   * @returns the \f$ i^{th} \f$ quadrature weight.
   */
  Real w(const unsigned int i) const
    { libmesh_assert (i < _weights.size()); return _weights[i]; }
  
  /**
   * Initializes the data structures to contain a quadrature rule
   * for an object of type \p type.  
   */
  void init (const ElemType _type=INVALID_ELEM,
	     unsigned int p_level=0);

  /**
   * @returns the order of the quadrature rule.   
   */
  Order get_order() const { return static_cast<Order>(_order + _p_level); }
  
  /**
   * Prints information relevant to the quadrature rule, by default to
   * libMesh::out.
   */
  void print_info(std::ostream& os=libMesh::out) const;

  /**
   * 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.
   */
  void scale(std::pair<Real, Real> old_range,
	     std::pair<Real, Real> new_range);

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

  /**
   * 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.
   */
  virtual bool shapes_need_reinit() { return false; }

  /**
   * 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!
   */
  bool allow_rules_with_negative_weights;
  
protected:

  
  /**
   * 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.
   */
  virtual void init_0D (const ElemType _type=INVALID_ELEM,
			unsigned int p_level=0);

  /**
   * 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.
   */
  virtual void init_1D (const ElemType _type=INVALID_ELEM,
			unsigned int p_level=0) = 0;

  /**
   * 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 \p DEBUG defined) when called.
   */
  virtual void init_2D (const ElemType,
			unsigned int =0)
#ifndef DEBUG
  {}
#else
  {  
    libMesh::err << "ERROR: Seems as if this quadrature rule" << std::endl
	          << " is not implemented for 2D." << std::endl;
    libmesh_error();
  }
#endif

  /**
   * 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 \p DEBUG defined) when called.
   */
  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
  
  
  /**
   * Computes the tensor product of
   * two 1D rules and returns a 2D rule.
   * Used in the init_2D routines for
   * quadrilateral element types.
   */
  void tensor_product_quad (const QBase& q1D);

  /**
   * 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.
   */
  void tensor_product_hex (const QBase& q1D);
  
  /**
   * Computes the tensor product of
   * a 1D quadrature rule and a 2D
   * quadrature rule.
   * Used in the init_3D routines for
   * prismatic element types.
   */
  void tensor_product_prism (const QBase& q1D, const QBase& q2D);


  
  /**
   * The dimension
   */
  const unsigned int _dim;
  
  /**
   * The order of the quadrature rule.
   */ 
  const Order _order;

  /**
   * The type of element for which the current values have
   * been computed.
   */
  ElemType _type;

  /**
   * The p level of element for which the current values have
   * been computed.
   */
  unsigned int _p_level;

  /**
   * The reference element locations of the
   * quadrature points.
   */
  std::vector<Point> _points;

  /**
   * The value of the quadrature weights.
   */
  std::vector<Real> _weights;
};





// ------------------------------------------------------------
// QBase class members

inline
QBase::QBase(const unsigned int d,
	     const Order o) :
  allow_rules_with_negative_weights(true),
  _dim(d),
  _order(o),
  _type(INVALID_ELEM),
  _p_level(0)
{
}




inline
void QBase::print_info(std::ostream& os) const
{
  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 << ":\n"
	 << "  "
	 << _points[qp]
	 << " Weight:\n "
	 << "  w=" << _weights[qp] << "\n" << std::endl;
    }
}

} // namespace libMesh

#endif