This file is indexed.

/usr/include/dune/geometry/generalvertexorder.hh is in libdune-geometry-dev 2.3.1-1.

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
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:

#ifndef DUNE_GEOMETRY_GENERALVERTEXORDER_HH
#define DUNE_GEOMETRY_GENERALVERTEXORDER_HH

#include <algorithm>
#include <cassert>
#include <cstddef>
#include <functional>
#include <iterator>
#include <vector>

#include <dune/common/iteratorfacades.hh>

#include "type.hh"
#include <dune/geometry/referenceelements.hh>

namespace Dune {

  /**
   * \brief Algorithm to reduce vertex order information
   *
   * \code
   * #include <dune/geometry/generalvertexorder.hh>
   * \endcode
   *
   * \param inBegin Start of the range of ids to reduce.
   * \param inEnd   End of the range of ids to reduce.
   * \param outIt   Start of the sequence where to store the result.
   *
   * \c inBegin and \c inEnd must be ForwardIterators; their \c value_type may
   * constant.  \c outIt must be an OutputIterator and must allow \c
   * std::distance(inBegin,inEnd) increments.  Complexity is quadratic.
   *
   * \sa GeneralVertexOrder, VertexOrderByIdFactory
   */
  template<class InIterator, class OutIterator>
  void reduceOrder(const InIterator& inBegin, const InIterator& inEnd,
                   OutIterator outIt)
  {
    typedef std::less<
        typename std::iterator_traits<InIterator>::value_type
        > less_t;
    const less_t less = less_t();

    for(InIterator inIt = inBegin; inIt != inEnd; ++inIt, ++outIt)
      *outIt = std::count_if(inBegin, inEnd, std::bind2nd(less, *inIt));
  }

  //! Class providing information on the ordering of vertices
  /**
   * \tparam dim    Dimension of the entity this class provides ordering
   *                information for.
   * \tparam Index_ Type of the indices.  Must be integral, may be
   *                non-negative.
   *
   * This class provides ordering information for all codimensions, including
   * the element itself.
   *
   * \warning The Interface of the VertexOrder stuff is subject to change.  It
   *          is currently needed to use some global-valued finite elements
   *          from dune-localfunctions.
   *
   * \sa reduceOrder(), VertexOrderByIdFactory
   */
  template<std::size_t dim, class Index_ = std::size_t>
  class GeneralVertexOrder {
    typedef ReferenceElement<double, dim> RefElem;
    typedef ReferenceElements<double, dim> RefElems;

    const RefElem& refelem;
    GeometryType gt;
    std::vector<Index_> vertexOrder;

  public:
    //! Type of indices
    typedef Index_ Index;

    //! Iterate over the vertex indices of some sub-entity
    class iterator;

    //! export the dimension of the entity we provide information for
    static const std::size_t dimension = dim;
    //! get type of the entity's geometry
    const GeometryType &type() const { return gt; }

    //! construct a GeneralVertexOrder
    /**
     * \param gt_     Geometry type of the entity we provide information for.
     * \param inBegin Start of the range of vertex ids.
     * \param inEnd   End of the range of vertex ids.
     *
     * \c inBegin and \c inEnd denote the range of vertes ids to provide.
     * This class stores a reduced copy of the ids, converted to type Index.
     */
    template<class InIterator>
    GeneralVertexOrder(const GeometryType& gt_, const InIterator &inBegin,
                       const InIterator &inEnd) :
      refelem(RefElems::general(gt_)), gt(gt_),
      vertexOrder(refelem.size(dim))
    { reduceOrder(inBegin, inEnd, vertexOrder.begin()); }

    //! get begin iterator for the vertex indices of some sub-entity
    /**
     * \param codim     Codimension of the sub-entity.
     * \param subEntity Index of the sub-entity within that codimension.
     */
    iterator begin(std::size_t codim, std::size_t subEntity) const
    { return iterator(*this, codim, subEntity); }
    //! get end iterator for the vertex indices of some sub-entity
    /**
     * \param codim     Codimension of the sub-entity.
     * \param subEntity Index of the sub-entity within that codimension.
     */
    iterator end(std::size_t codim, std::size_t subEntity) const {
      return iterator(*this, codim, subEntity,
                      refelem.size(subEntity, codim, dim));
    }

    //! get a vector of reduced indices for some sub-entity
    /**
     * \param codim     Codimension of the sub-entity.
     * \param subEntity Index of the sub-entity within that codimension.

     * \param order     Where to store the result.  This function resizes the
     *                  vector to the suitable size.
     */
    void getReduced(std::size_t codim, std::size_t subEntity,
                    std::vector<Index>& order) const
    {
      order.resize(refelem.size(subEntity, codim, dim));
      reduceOrder(begin(codim, subEntity), end(codim, subEntity),
                  order.begin());
    }
  };

  //! Iterate over the vertex indices of some sub-entity
  /**
   * This is a random access iterator with constant \c value_type.
   */
  template<std::size_t dim, class Index_>
  class GeneralVertexOrder<dim, Index_>::iterator :
    public Dune::RandomAccessIteratorFacade<iterator, const Index_>
  {
    const GeneralVertexOrder *order;
    std::size_t codim;
    std::size_t subEntity;
    std::size_t vertex;

    iterator(const GeneralVertexOrder &order_, std::size_t codim_,
             std::size_t subEntity_, std::size_t vertex_ = 0) :
      order(&order_), codim(codim_), subEntity(subEntity_), vertex(vertex_)
    { }

  public:
    const Index &dereference() const {
      return order->vertexOrder[order->refelem.subEntity(subEntity, codim,
                                                         vertex, dim)];
    }
    const Index &elementAt(std::ptrdiff_t n) const {
      return order->vertexOrder[order->refelem.subEntity(subEntity, codim,
                                                         vertex+n, dim)];
    }
    bool equals(const iterator &other) const {
      return order == other.order && codim == other.codim &&
             subEntity == other.subEntity && vertex == other.vertex;
    }
    void increment() { ++vertex; }
    void decrement() { --vertex; }
    void advance(std::ptrdiff_t n) { vertex += n; }
    std::ptrdiff_t distanceTo(const iterator &other) const {
      // make sure we reference the same container
      assert(order == other.order && codim == other.codim &&
             subEntity == other.subEntity);
      if(vertex < other.vertex) return other.vertex - vertex;
      else return -static_cast<std::ptrdiff_t>(vertex - other.vertex);
    }

    friend class GeneralVertexOrder<dim, Index>;

    //! public default constructor
    /**
     * The contructed iterator object will have a singular value.  The only
     * valid operations will be assignment of a non-singular value and
     * destruction, all other operations will result in undefined behaviour.
     */
    iterator() { }
  };
} // namespace Dune

#endif // DUNE_GEOMETRY_GENERALVERTEXORDER_HH