/usr/include/dolfin/geometry/SimplexQuadrature.h is in libdolfin-dev 1.4.0+dfsg-4.
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 | // Copyright (C) 2014 Anders Logg
//
// This file is part of DOLFIN.
//
// DOLFIN 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 3 of the License, or
// (at your option) any later version.
//
// DOLFIN 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 DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// First added: 2014-02-24
// Last changed: 2014-03-03
#ifndef __SIMPLEX_QUADRATURE_H
#define __SIMPLEX_QUADRATURE_H
#include <vector>
#include "Point.h"
namespace dolfin
{
class SimplexQuadrature
{
public:
/// Compute quadrature rule for simplex.
///
/// *Arguments*
/// coordinates (double *)
/// A flattened array of simplex coordinates of
/// dimension num_vertices x gdim = (tdim + 1)*gdim.
/// tdim (std::size_t)
/// The topological dimension of the simplex.
/// gdim (std::size_t)
/// The geometric dimension.
/// order (std::size_t)
/// The order of convergence of the quadrature rule.
///
/// *Returns*
/// std::pair<std::vector<double>, std::vector<double> >
/// An array of quadrature weights and a corresponding
/// flattened array of quadrature points.
static std::pair<std::vector<double>, std::vector<double> >
compute_quadrature_rule(const double* coordinates,
std::size_t tdim,
std::size_t gdim,
std::size_t order);
/// Compute quadrature rule for interval.
///
/// *Arguments*
/// coordinates (double *)
/// A flattened array of simplex coordinates of
/// dimension num_vertices x gdim = 2*gdim.
/// gdim (std::size_t)
/// The geometric dimension.
/// order (std::size_t)
/// The order of convergence of the quadrature rule.
///
/// *Returns*
/// std::pair<std::vector<double>, std::vector<double> >
/// An array of quadrature weights and a corresponding
/// flattened array of quadrature points.
static std::pair<std::vector<double>, std::vector<double> >
compute_quadrature_rule_interval(const double* coordinates,
std::size_t gdim,
std::size_t order);
/// Compute quadrature rule for triangle.
///
/// *Arguments*
/// coordinates (double *)
/// A flattened array of simplex coordinates of
/// dimension num_vertices x gdim = 3*gdim.
/// gdim (std::size_t)
/// The geometric dimension.
/// order (std::size_t)
/// The order of convergence of the quadrature rule.
///
/// *Returns*
/// std::pair<std::vector<double>, std::vector<double> >
/// An array of quadrature weights and a corresponding
/// flattened array of quadrature points.
static std::pair<std::vector<double>, std::vector<double> >
compute_quadrature_rule_triangle(const double* coordinates,
std::size_t gdim,
std::size_t order);
/// Compute quadrature rule for tetrahedron.
///
/// *Arguments*
/// coordinates (double *)
/// A flattened array of simplex coordinates of
/// dimension num_vertices x gdim = 4*gdim.
/// gdim (std::size_t)
/// The geometric dimension.
/// order (std::size_t)
/// The order of convergence of the quadrature rule.
///
/// *Returns*
/// std::pair<std::vector<double>, std::vector<double> >
/// An array of quadrature weights and a corresponding
/// flattened array of quadrature points.
static std::pair<std::vector<double>, std::vector<double> >
compute_quadrature_rule_tetrahedron(const double* coordinates,
std::size_t gdim,
std::size_t order);
};
}
#endif
|