Dune Core Modules (unstable)

boundaryiterators.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 
6 #ifndef DUNE_GRID_IO_FILE_VTK_BOUNDARYITERATORS_HH
7 #define DUNE_GRID_IO_FILE_VTK_BOUNDARYITERATORS_HH
8 
9 #include <iterator>
10 #include <memory>
11 
13 
14 #include <dune/grid/io/file/vtk/corner.hh>
15 #include <dune/grid/io/file/vtk/corneriterator.hh>
16 #include <dune/grid/io/file/vtk/functionwriter.hh>
17 
18 namespace Dune {
19 
22 
28  namespace VTK {
29 
31 
35  template<typename GV>
37  : public ForwardIteratorFacade
38  < BoundaryIterator<GV>,
39  const typename GV::Intersection,
40  const typename GV::Intersection&,
41  typename std::iterator_traits<typename GV::template Codim<0>::
42  Iterator>::difference_type>
43  {
44  public:
45  // reiterator the facades typedefs here
47  typedef const typename GV::Intersection Value;
48  typedef Value& Reference;
49  typedef typename GV::template Codim<0>::Iterator ElementIterator;
50  typedef typename GV::IntersectionIterator IntersectionIterator;
51  typedef typename std::iterator_traits<ElementIterator>::difference_type
52  DifferenceType;
53 
54  private:
55  typedef ForwardIteratorFacade<DerivedType, Value, Reference,
56  DifferenceType> Facade;
57 
58  const GV* gv;
59  ElementIterator eit;
60  std::shared_ptr<IntersectionIterator> iit;
61 
62  bool valid() const {
63  // we're valid if we're passed-the-end
64  if(eit == gv->template end<0>()) return true;
65  // or if we're on a boundary
66  if((*iit)->boundary() && !(*iit)->neighbor()) return true;
67  // otherwise we're invalid
68  return false;
69  }
70 
71  void basic_increment() {
72  ++*iit;
73  if(*iit == gv->iend(*eit)) {
74  iit.reset();
75  ++eit;
76  if(eit != gv->template end<0>())
77  iit.reset(new IntersectionIterator(gv->ibegin(*eit)));
78  }
79  }
80 
81  public:
82  Reference dereference() const {
83  return **iit;
84  }
85  bool equals(const DerivedType& other) const {
86  if(eit != other.eit) return false;
87 
88  // this is a bit tricky, since we may not compare iit if we are
89  // passed-the-end
90  bool mePassedTheEnd = eit == gv->template end<0>();
91  bool otherPassedTheEnd = other.eit == other.gv->template end<0>();
92 
93  // both passed-the-end => consider them equal
94  if(mePassedTheEnd && otherPassedTheEnd) return true;
95 
96  // one passed the end => not equal
97  if(mePassedTheEnd || otherPassedTheEnd) return false;
98 
99  // none passed-the-end => do their iit iterators match?
100  return *iit == *other.iit;
101  }
102 
103  void increment() {
104  basic_increment();
105  while(!valid()) basic_increment();
106  }
107 
109 
113  BoundaryIterator(const GV& gv_, const ElementIterator& eit_,
114  const IntersectionIterator& iit_)
115  : gv(&gv_), eit(eit_), iit(new IntersectionIterator(iit_))
116  {
117  while(!valid()) basic_increment();
118  }
120 
125  BoundaryIterator(const GV& gv_, const ElementIterator& eit_)
126  : gv(&gv_), eit(eit_)
127  {
128  if(eit != gv->template end<0>())
129  iit.reset(new IntersectionIterator(gv->ibegin(*eit)));
130 
131  while(!valid()) basic_increment();
132  }
134 
138  BoundaryIterator(const GV& gv_, bool end = false)
139  : gv(&gv_), eit(end ? gv->template end<0>() : gv->template begin<0>())
140  {
141  if(eit != gv->template end<0>())
142  iit.reset(new IntersectionIterator(gv->ibegin(*eit)));
143 
144  while(!valid()) basic_increment();
145  }
146  };
147 
148  template<typename ElementIndexSet>
149  class IntersectionIndexSet {
150  const ElementIndexSet& eis;
151 
152  public:
153  IntersectionIndexSet(const ElementIndexSet& eis_)
154  : eis(eis_)
155  { }
156  };
157 
158  template<typename GV>
159  class NonConformingBoundaryIteratorFactory {
160  const GV& gv;
161 
162  public:
163  static const unsigned dimCell = GV::dimension-1;
164 
165  typedef typename GV::Intersection Cell;
166  typedef BoundaryIterator<GV> CellIterator;
167 
168  typedef VTK::Corner<Cell> Corner;
169  typedef VTK::CornerIterator<CellIterator> CornerIterator;
170 
171  typedef Corner Point;
172  typedef CornerIterator PointIterator;
173 
174  typedef NonConformingConnectivityWriter<Cell> ConnectivityWriter;
175  typedef typename GV::Communication Communication;
176 
177  explicit NonConformingBoundaryIteratorFactory(const GV& gv_)
178  : gv(gv_)
179  { }
180 
181  CellIterator beginCells() const {
182  return CellIterator(gv);
183  }
184  CellIterator endCells() const {
185  return CellIterator(gv, true);
186  }
187 
188  CornerIterator beginCorners() const {
189  return CornerIterator(beginCells(), endCells());
190  }
191  CornerIterator endCorners() const {
192  return CornerIterator(endCells());
193  }
194 
195  PointIterator beginPoints() const { return beginCorners(); }
196  PointIterator endPoints() const { return endCorners(); }
197 
198  ConnectivityWriter makeConnectivity() const {
199  return ConnectivityWriter();
200  }
201  const Communication& comm() const {
202  return gv.comm();
203  }
204  };
205 
206  } // namespace VTK
207 
209 
210 } // namespace Dune
211 
212 #endif // DUNE_GRID_IO_FILE_VTK_BOUNDARYITERATORS_HH
Base class for stl conformant forward iterators.
Definition: iteratorfacades.hh:141
iterate over the GridViews boundary intersections
Definition: boundaryiterators.hh:43
BoundaryIterator(const GV &gv_, const ElementIterator &eit_)
construct a BoundaryIterator
Definition: boundaryiterators.hh:125
BoundaryIterator(const GV &gv_, bool end=false)
construct a BoundaryIterator
Definition: boundaryiterators.hh:138
BoundaryIterator(const GV &gv_, const ElementIterator &eit_, const IntersectionIterator &iit_)
construct a BoundaryIterator
Definition: boundaryiterators.hh:113
concept Intersection
Model of an intersection.
Definition: intersection.hh:23
concept IntersectionIterator
Model of an intersection iterator.
Definition: intersectioniterator.hh:21
This file implements iterator facade classes for writing stl conformant iterators.
Dune namespace.
Definition: alignedallocator.hh:13
Static tag representing a codimension.
Definition: dimension.hh:24
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 15, 22:30, 2024)