/usr/include/libmesh/quadrature_monomial.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 | // $Id: quadrature_monomial.h 3874 2010-07-02 21:57:26Z roystgnr $
// 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_monomial_h__
#define __quadrature_monomial_h__
// Local includes
#include "quadrature.h"
namespace libMesh
{
/**
* This class defines alternate quadrature rules on
* "tensor-product" elements (QUADs and HEXes) which can be
* useful when integrating monomial finite element bases.
*
* While tensor product rules are ideal for integrating
* bi/tri-linear, bi/tri-quadratic, etc. (i.e. tensor product)
* bases (which consist of incomplete polynomials up to degree=
* dim*p) they are not optimal for the MONOMIAL or FEXYZ bases,
* which consist of complete polynomials of degree=p.
*
* This class is implemented to provide quadrature rules which are
* more efficient than tensor product rules when they are available,
* and fall back on Gaussian quadrature rules when necessary.
*
* A number of these rules have been helpfully collected in electronic form by:
*
* Prof. Ronald Cools
* Katholieke Universiteit Leuven,
* Dept. Computerwetenschappen
* http://www.cs.kuleuven.ac.be/~nines/research/ecf/ecf.html
*
* (A username and password to access the tables is available by request.)
*
* We also provide the original reference for each rule, as available,
* in the source code file.
*
* @author John W. Peterson, 2008
*/
class QMonomial : public QBase
{
public:
/**
* Constructor. Declares the order of the quadrature rule.
*/
QMonomial (const unsigned int _dim,
const Order _order=INVALID_ORDER);
/**
* Destructor.
*/
~QMonomial();
/**
* @returns \p QMONOMIAL
*/
QuadratureType type() const { return QMONOMIAL; }
private:
void init_1D (const ElemType,
unsigned int =0)
{
// See about making this non-pure virtual in the base class?
libmesh_error();
}
/**
* More efficient rules for QUADs
*/
void init_2D (const ElemType _type=INVALID_ELEM,
unsigned int p_level=0);
/**
* More efficient rules for HEXes
*/
void init_3D (const ElemType _type=INVALID_ELEM,
unsigned int p_level=0);
/**
* Wissmann published three interesting "partially symmetric" rules
* for integrating degree 4, 6, and 8 polynomials exactly on QUADs.
* These rules have all positive weights, all points inside the
* reference element, and have fewer points than tensor-product
* rules of equivalent order, making them superior to those rules
* for monomial bases.
*
* J. W. Wissman and T. Becker, Partially symmetric cubature
* formulas for even degrees of exactness, SIAM J. Numer. Anal. 23
* (1986), 676--685.
*/
void wissmann_rule(const Real rule_data[][3],
const unsigned int n_pts);
/**
* Stroud's rules for QUADs and HEXes can have one of several
* different types of symmetry. The rule_symmetry array describes
* how the different lines of the rule_data array are to be
* applied. The different rule_symmetry possibilities are:
* 0) Origin or single-point: (x,y)
* Fully-symmetric, 3 cases:
* 1) (x,y) -> (x,y), (-x,y), (x,-y), (-x,-y)
* (y,x), (-y,x), (y,-x), (-y,-x)
* 2) (x,x) -> (x,x), (-x,x), (x,-x), (-x,-x)
* 3) (x,0) -> (x,0), (-x,0), (0, x), ( 0,-x)
* 4) Rotational Invariant, (x,y) -> (x,y), (-x,-y), (-y, x), (y,-x)
* 5) Partial Symmetry, (x,y) -> (x,y), (-x, y) [x!=0]
* 6) Rectangular Symmetry, (x,y) -> (x,y), (-x, y), (-x,-y), (x,-y)
* 7) Central Symmetry, (0,y) -> (0,y), ( 0,-y)
*
* Not all rules with these symmetries are due to Stroud, however,
* his book is probably the most frequently-cited compendium of
* quadrature rules and later authors certainly built upon his work.
*/
void stroud_rule(const Real rule_data[][3],
const unsigned int* rule_symmetry,
const unsigned int n_pts);
/**
* Rules from Kim and Song, Comm. Korean Math. Soc vol. 13, no. 4, 1998, pp. 913-931.
* The rules are obtained by considering the group G^{rot} of rotations of the reference
* hex, and the invariant polynomials of this group.
*
* In Kim and Song's rules, quadrauture points are described by the following points
* and their unique permutations under the G^{rot} group:
*
* 0.) (0,0,0) ( 1 perm ) -> [0, 0, 0]
* 1.) (x,0,0) ( 6 perms) -> [x, 0, 0], [0, -x, 0], [-x, 0, 0], [0, x, 0], [0, 0, -x], [0, 0, x]
* 2.) (x,x,0) (12 perms) -> [x, x, 0], [x, -x, 0], [-x, -x, 0], [-x, x, 0], [x, 0, -x], [x, 0, x],
* [0, x, -x], [0, x, x], [0, -x, -x], [-x, 0, -x], [0, -x, x], [-x, 0, x]
* 3.) (x,y,0) (24 perms) -> [x, y, 0], [y, -x, 0], [-x, -y, 0], [-y, x, 0], [x, 0, -y], [x, -y, 0],
* [x, 0, y], [0, y, -x], [-x, y, 0],
* [0, y, x], [y, 0, -x], [0, -y, -x], [-y, 0, -x], [y, x, 0], [-y, -x, 0],
* [y, 0, x], [0, -y, x], [-y, 0, x],
* [-x, 0, y], [0, -x, -y], [0, -x, y], [-x, 0, -y], [0, x, y], [0, x, -y]
* 4.) (x,x,x) ( 8 perms) -> [x, x, x], [x, -x, x], [-x, -x, x], [-x, x, x], [x, x, -x], [x, -x, -x],
* [-x, x, -x], [-x, -x, -x]
* 5.) (x,x,z) (24 perms) -> [x, x, z], [x, -x, z], [-x, -x, z], [-x, x, z], [x, z, -x], [x, -x, -z],
* [x, -z, x], [z, x, -x], [-x, x, -z],
* [-z, x, x], [x, -z, -x], [-z, -x, -x], [-x, z, -x], [x, x, -z],
* [-x, -x, -z], [x, z, x], [z, -x, x],
* [-x, -z, x], [-x, z, x], [z, -x, -x], [-z, -x, x], [-x, -z, -x],
* [z, x, x], [-z, x, -x]
* 6.) (x,y,z) (24 perms) -> [x, y, z], [y, -x, z], [-x, -y, z], [-y, x, z], [x, z, -y], [x, -y, -z],
* [x, -z, y], [z, y, -x], [-x, y, -z],
* [-z, y, x], [y, -z, -x], [-z, -y, -x], [-y, z, -x], [y, x, -z], [-y, -x, -z],
* [y, z, x], [z, -y, x], [-y, -z, x], [-x, z, y], [z, -x, -y],
* [-z, -x, y], [-x, -z, -y], [z, x, y], [-z, x, -y]
*
* Only two of Kim and Song's rules are particularly useful for FEM calculations: the degree 7,
* 38-point rule and their degree 8, 47-point rule. The others either contain negative weights or points
* outside the reference interval. The points and weights, to 32 digits, were obtained from:
* Ronald Cools' website (http://www.cs.kuleuven.ac.be/~nines/research/ecf/ecf.html) and
* the unique permutations of G^{rot} were computed by me [JWP] using Maple.
*/
void kim_rule(const Real rule_data[][4],
const unsigned int* rule_id,
const unsigned int n_pts);
};
} // namespace libMesh
#endif
|