This file is indexed.

/usr/include/dolfin/mesh/MeshQuality.h is in libdolfin-dev 2016.2.0-2.

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
// Copyright (C) 2013 Garth N. Wells
//
// 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:  2013-10-07
// Last changed:

#ifndef __MESH_QUALITY_H
#define __MESH_QUALITY_H

#include <string>
#include <utility>
#include <vector>
#include <boost/multi_array.hpp>
#include <memory>
#include "Cell.h"

namespace dolfin
{

    class Mesh;

  /// The class provides functions to quantify mesh quality

  class MeshQuality
  {
  public:

    /// Compute the radius ratio for all cells.
    ///
    /// *Returns*
    ///     CellFunction<double>
    ///         The cell radius ratio radius ratio geometric_dimension *
    ///         * inradius / circumradius (geometric_dimension
    ///         is normalization factor). It has range zero to one.
    ///         Zero indicates a degenerate element.
    ///
    /// *Example*
    ///     .. note::
    ///
    ///         std::shared_ptr<Mesh> mesh(new UnitCubeMesh(4, 4, 4));
    ///         CellFunction<double> = MeshQuality::radius_ratio(mesh);
    static CellFunction<double>
      radius_ratios(std::shared_ptr<const Mesh> mesh);


    /// Compute the minimum and maximum radius ratio of cells
    /// (across all processes)
    ///
    /// *Returns*
    ///     std::pair<double, double>
    ///         The [minimum, maximum] cell radii ratio (geometric_dimension *
    ///         * inradius / circumradius, geometric_dimension
    ///         is normalization factor). It has range zero to one.
    ///         Zero indicates a degenerate element.
    ///
    /// *Example*
    ///     .. note::
    ///
    ///         Mesh  UnitCubeMesh(4, 4, 4);
    ///         std::pair<double, double> ratios
    ///            = MeshQuality::radius_ratio_min_max(mesh);
    ///         double min_ratio = ratios.first;
    ///         double max_ratio = ratios.second;
    static std::pair<double, double> radius_ratio_min_max(const Mesh& mesh);


    /// Create (ratio, number of cells) data for creating a histogram
    /// of cell quality
    static std::pair<std::vector<double>, std::vector<double> >
      radius_ratio_histogram_data(const Mesh& mesh,
                                  std::size_t num_bins = 50);

    /// Create Matplotlib string to plot cell quality histogram
    static std::string
      radius_ratio_matplotlib_histogram(const Mesh& mesh,
					std::size_t num_intervals = 50);

    /// Get internal dihedral angles of a tetrahedral cell
    static void dihedral_angles(const Cell& cell, std::vector<double>& dh_angle); 

    /// Get internal minimum and maximum dihedral angles of a 3D mesh
    static std::pair<double, double> 
        dihedral_angles_min_max(const Mesh& mesh);

    /// Create (dihedral angles, number of cells) data for creating a histogram
    /// of dihedral
    static std::pair<std::vector<double>, std::vector<double> >
      dihedral_angles_histogram_data(const Mesh& mesh,
                                  std::size_t num_bins = 100);

    /// Create Matplotlib string to plot dihedral angles quality histogram
    static std::string
      dihedral_angles_matplotlib_histogram(const Mesh& mesh,
                                std::size_t num_intervals = 100);
  };

}

#endif