Dune Core Modules (2.3.1)

generalvertexorder.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_GEOMETRY_GENERALVERTEXORDER_HH
5#define DUNE_GEOMETRY_GENERALVERTEXORDER_HH
6
7#include <algorithm>
8#include <cassert>
9#include <cstddef>
10#include <functional>
11#include <iterator>
12#include <vector>
13
15
16#include "type.hh"
17#include <dune/geometry/referenceelements.hh>
18
19namespace Dune {
20
38 template<class InIterator, class OutIterator>
39 void reduceOrder(const InIterator& inBegin, const InIterator& inEnd,
40 OutIterator outIt)
41 {
42 typedef std::less<
43 typename std::iterator_traits<InIterator>::value_type
44 > less_t;
45 const less_t less = less_t();
46
47 for(InIterator inIt = inBegin; inIt != inEnd; ++inIt, ++outIt)
48 *outIt = std::count_if(inBegin, inEnd, std::bind2nd(less, *inIt));
49 }
50
52
67 template<std::size_t dim, class Index_ = std::size_t>
71
72 const RefElem& refelem;
73 GeometryType gt;
74 std::vector<Index_> vertexOrder;
75
76 public:
78 typedef Index_ Index;
79
81 class iterator;
82
84 static const std::size_t dimension = dim;
86 const GeometryType &type() const { return gt; }
87
89
97 template<class InIterator>
98 GeneralVertexOrder(const GeometryType& gt_, const InIterator &inBegin,
99 const InIterator &inEnd) :
100 refelem(RefElems::general(gt_)), gt(gt_),
101 vertexOrder(refelem.size(dim))
102 { reduceOrder(inBegin, inEnd, vertexOrder.begin()); }
103
105
109 iterator begin(std::size_t codim, std::size_t subEntity) const
110 { return iterator(*this, codim, subEntity); }
112
116 iterator end(std::size_t codim, std::size_t subEntity) const {
117 return iterator(*this, codim, subEntity,
118 refelem.size(subEntity, codim, dim));
119 }
120
122
129 void getReduced(std::size_t codim, std::size_t subEntity,
130 std::vector<Index>& order) const
131 {
132 order.resize(refelem.size(subEntity, codim, dim));
133 reduceOrder(begin(codim, subEntity), end(codim, subEntity),
134 order.begin());
135 }
136 };
137
139
142 template<std::size_t dim, class Index_>
143 class GeneralVertexOrder<dim, Index_>::iterator :
144 public Dune::RandomAccessIteratorFacade<iterator, const Index_>
145 {
146 const GeneralVertexOrder *order;
147 std::size_t codim;
148 std::size_t subEntity;
149 std::size_t vertex;
150
151 iterator(const GeneralVertexOrder &order_, std::size_t codim_,
152 std::size_t subEntity_, std::size_t vertex_ = 0) :
153 order(&order_), codim(codim_), subEntity(subEntity_), vertex(vertex_)
154 { }
155
156 public:
157 const Index &dereference() const {
158 return order->vertexOrder[order->refelem.subEntity(subEntity, codim,
159 vertex, dim)];
160 }
161 const Index &elementAt(std::ptrdiff_t n) const {
162 return order->vertexOrder[order->refelem.subEntity(subEntity, codim,
163 vertex+n, dim)];
164 }
165 bool equals(const iterator &other) const {
166 return order == other.order && codim == other.codim &&
167 subEntity == other.subEntity && vertex == other.vertex;
168 }
169 void increment() { ++vertex; }
170 void decrement() { --vertex; }
171 void advance(std::ptrdiff_t n) { vertex += n; }
172 std::ptrdiff_t distanceTo(const iterator &other) const {
173 // make sure we reference the same container
174 assert(order == other.order && codim == other.codim &&
175 subEntity == other.subEntity);
176 if(vertex < other.vertex) return other.vertex - vertex;
177 else return -static_cast<std::ptrdiff_t>(vertex - other.vertex);
178 }
179
180 friend class GeneralVertexOrder<dim, Index>;
181
183
189 };
190} // namespace Dune
191
192#endif // DUNE_GEOMETRY_GENERALVERTEXORDER_HH
Iterate over the vertex indices of some sub-entity.
Definition: generalvertexorder.hh:145
iterator()
public default constructor
Definition: generalvertexorder.hh:188
Class providing information on the ordering of vertices.
Definition: generalvertexorder.hh:68
Index_ Index
Type of indices.
Definition: generalvertexorder.hh:78
const GeometryType & type() const
get type of the entity's geometry
Definition: generalvertexorder.hh:86
static const std::size_t dimension
export the dimension of the entity we provide information for
Definition: generalvertexorder.hh:84
void getReduced(std::size_t codim, std::size_t subEntity, std::vector< Index > &order) const
get a vector of reduced indices for some sub-entity
Definition: generalvertexorder.hh:129
iterator end(std::size_t codim, std::size_t subEntity) const
get end iterator for the vertex indices of some sub-entity
Definition: generalvertexorder.hh:116
GeneralVertexOrder(const GeometryType &gt_, const InIterator &inBegin, const InIterator &inEnd)
construct a GeneralVertexOrder
Definition: generalvertexorder.hh:98
iterator begin(std::size_t codim, std::size_t subEntity) const
get begin iterator for the vertex indices of some sub-entity
Definition: generalvertexorder.hh:109
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:431
This class provides access to geometric and topological properties of a reference element....
Definition: referenceelements.hh:58
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:86
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:132
This file implements iterator facade classes for writing stl conformant iterators.
Dune namespace.
Definition: alignment.hh:14
void reduceOrder(const InIterator &inBegin, const InIterator &inEnd, OutIterator outIt)
Algorithm to reduce vertex order information.
Definition: generalvertexorder.hh:39
Class providing access to the singletons of the reference elements. Special methods are available for...
Definition: referenceelements.hh:563
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)