DUNE PDELab (2.7)

treeiterator.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_ALBERTA_TREEITERATOR_HH
5#define DUNE_ALBERTA_TREEITERATOR_HH
6
7#include <dune/common/hybridutilities.hh>
8#include <dune/common/std/utility.hh>
10
13
14#if HAVE_ALBERTA
15
16namespace Dune
17{
18
19 // AlbertaMarkerVector
20 // -------------------
21
30 template< int dim, int dimworld >
32 {
34
36
37 //friend class AlbertaGrid< dim, dimworld >;
38
39 static const int dimension = Grid::dimension;
40
41 typedef Alberta::HierarchyDofNumbering< dimension > DofNumbering;
42 typedef Alberta::ElementInfo< dimension > ElementInfo;
43
44 template< bool >
45 struct NoMarkSubEntities;
46 template< bool >
47 struct MarkSubEntities;
48
49 public:
51 explicit AlbertaMarkerVector ( const DofNumbering &dofNumbering )
52 : dofNumbering_( dofNumbering )
53 {
54 for( int codim = 0; codim <= dimension; ++codim )
55 marker_[ codim ] = 0;
56 }
57
58 AlbertaMarkerVector ( const This &other )
59 : dofNumbering_( other.dofNumbering_ )
60 {
61 for( int codim = 0; codim <= dimension; ++codim )
62 marker_[ codim ] = 0;
63 }
64
65 ~AlbertaMarkerVector ()
66 {
67 clear();
68 }
69
70 private:
71 This &operator= ( const This & );
72
73 public:
75 template< int codim >
76 bool subEntityOnElement ( const ElementInfo &elementInfo, int subEntity ) const;
77
78 template< int firstCodim, class Iterator >
79 void markSubEntities ( const Iterator &begin, const Iterator &end );
80
81 void clear ()
82 {
83 for( int codim = 0; codim <= dimension; ++codim )
84 {
85 if( marker_[ codim ] != 0 )
86 delete[] marker_[ codim ];
87 marker_[ codim ] = 0;
88 }
89 }
90
92 bool up2Date () const
93 {
94 return (marker_[ dimension ] != 0);
95 }
96
98 void print ( std::ostream &out = std::cout ) const;
99
100 private:
101 const DofNumbering &dofNumbering_;
102 int *marker_[ dimension+1 ];
103 };
104
105
106
107 // AlbertaMarkerVector::NoMarkSubEntities
108 // --------------------------------------
109
110 template< int dim, int dimworld >
111 template< bool >
112 struct AlbertaMarkerVector< dim, dimworld >::NoMarkSubEntities
113 {
114 template< int firstCodim, class Iterator >
115 static void mark ( const DofNumbering &dofNumbering, int *(&marker)[ dimension + 1 ],
116 const Iterator &begin, const Iterator &end )
117 {}
118 };
119
120
121
122 // AlbertaMarkerVector::MarkSubEntities
123 // ------------------------------------
124
125 template< int dim, int dimworld >
126 template< bool >
127 struct AlbertaMarkerVector< dim, dimworld >::MarkSubEntities
128 {
129 template< int codim >
130 struct Codim
131 {
132 static const int numSubEntities = Alberta::NumSubEntities< dimension, codim >::value;
133
134 typedef Alberta::ElementInfo< dimension > ElementInfo;
135
136 static void apply ( const DofNumbering &dofNumbering,
137 int *(&marker)[ dimension + 1 ],
138 const ElementInfo &elementInfo )
139 {
140 int *array = marker[ codim ];
141
142 const int index = dofNumbering( elementInfo, 0, 0 );
143 for( int i = 0; i < numSubEntities; ++i )
144 {
145 int &mark = array[ dofNumbering( elementInfo, codim, i ) ];
146 mark = std::max( index, mark );
147 }
148 }
149 };
150
151 template< int firstCodim, class Iterator >
152 static void mark ( const DofNumbering &dofNumbering, int *(&marker)[ dimension + 1 ],
153 const Iterator &begin, const Iterator &end )
154 {
155 for( int codim = firstCodim; codim <= dimension; ++codim )
156 {
157 const int size = dofNumbering.size( codim );
158 marker[ codim ] = new int[ size ];
159
160 int *array = marker[ codim ];
161 for( int i = 0; i < size; ++i )
162 array[ i ] = -1;
163 }
164
165 for( Iterator it = begin; it != end; ++it )
166 {
167 const ElementInfo &elementInfo = it->impl().elementInfo();
168 Hybrid::forEach( Std::make_index_sequence< dimension+1-firstCodim >{},
169 [ & ]( auto i ){ Codim< i+firstCodim >::apply( dofNumbering, marker, elementInfo ); } );
170 }
171 }
172 };
173
174
175
176 // AlbertaGridTreeIterator
177 // -----------------------
178
182 template< int codim, class GridImp, bool leafIterator >
184 {
186
187 public:
188 static const int dimension = GridImp::dimension;
189 static const int codimension = codim;
190 static const int dimensionworld = GridImp::dimensionworld;
191
192 private:
193 friend class AlbertaGrid< dimension, dimensionworld >;
194
195 static const int numSubEntities
196 = Alberta::NumSubEntities< dimension, codimension >::value;
197
198 public:
199 typedef Alberta::MeshPointer< dimension > MeshPointer;
200 typedef typename MeshPointer::MacroIterator MacroIterator;
201
202 typedef typename GridImp::template Codim< codim >::Entity Entity;
204 typedef typename EntityObject::ImplementationType EntityImp;
205 typedef typename EntityImp::ElementInfo ElementInfo;
206
208
210
213
215 This &operator= ( const This &other );
216
218 AlbertaGridTreeIterator ( const GridImp &grid, int travLevel );
219
222 const MarkerVector *marker,
223 int travLevel );
224
226 bool equals ( const This &other ) const
227 {
228 return entity_.impl().equals( other.entity_.impl() );
229 }
230
232 Entity &dereference () const
233 {
234 return entity_;
235 }
236
238 int level () const
239 {
240 return entity_.impl().level();
241 }
242
244 void increment();
245
246 protected:
248 const GridImp &grid () const
249 {
250 return entity_.impl().grid();
251 }
252
253 private:
254 void nextElement ( ElementInfo &elementInfo );
255 void nextElementStop (ElementInfo &elementInfo );
256 bool stopAtElement ( const ElementInfo &elementInfo ) const;
257
258 void goNext ( ElementInfo &elementInfo );
259 void goNext ( const std::integral_constant< int, 0 > cdVariable,
260 ElementInfo &elementInfo );
261 void goNext ( const std::integral_constant< int, 1 > cdVariable,
262 ElementInfo &elementInfo );
263 template< int cd >
264 void goNext ( const std::integral_constant< int, cd > cdVariable,
265 ElementInfo &elementInfo );
266
267 mutable Entity entity_;
268
270 int level_;
271
273 int subEntity_;
274
275 MacroIterator macroIterator_;
276
277 // knows on which element a point,edge,face is viewed
278 const MarkerVector *marker_;
279 };
280
281
282
283 // Implementation of AlbertaMarkerVector
284 // -------------------------------------
285
286 template< int dim, int dimworld >
287 template< int codim >
289 ::subEntityOnElement ( const ElementInfo &elementInfo, int subEntity ) const
290 {
291 assert( marker_[ codim ] != 0 );
292
293 const int subIndex = dofNumbering_( elementInfo, codim, subEntity );
294 const int markIndex = marker_[ codim ][ subIndex ];
295 assert( (markIndex >= 0) );
296
297 const int index = dofNumbering_( elementInfo, 0, 0 );
298 return (markIndex == index);
299 }
300
301
302 template< int dim, int dimworld >
303 template< int firstCodim, class Iterator >
305 ::markSubEntities ( const Iterator &begin, const Iterator &end )
306 {
307 clear();
308 std::conditional< (firstCodim <= dimension), MarkSubEntities<true>, NoMarkSubEntities<false> >::type
309 ::template mark< firstCodim, Iterator >( dofNumbering_, marker_, begin, end );
310 }
311
312
313 template< int dim, int dimworld >
314 inline void AlbertaMarkerVector< dim, dimworld >::print ( std::ostream &out ) const
315 {
316 for( int codim = 1; codim <= dimension; ++codim )
317 {
318 int *marker = marker_[ codim ];
319 if( marker != 0 )
320 {
321 const int size = dofNumbering_.size( codim );
322 out << std::endl;
323 out << "Codimension " << codim << " (" << size << " entries)" << std::endl;
324 for( int i = 0; i < size; ++i )
325 out << "subentity " << i << " visited on Element " << marker[ i ] << std::endl;
326 }
327 }
328 }
329
330
331
332 // Implementation of AlbertaGridTreeIterator
333 // -----------------------------------------
334
335 template< int codim, class GridImp, bool leafIterator >
338 : entity_(),
339 level_( -1 ),
340 subEntity_( -1 ),
341 macroIterator_(),
342 marker_( NULL )
343 {}
344
345 template< int codim, class GridImp, bool leafIterator >
346 inline AlbertaGridTreeIterator< codim, GridImp, leafIterator >
347 ::AlbertaGridTreeIterator ( const GridImp &grid,
348 const MarkerVector *marker,
349 int travLevel )
350 : entity_( EntityImp( grid ) ),
351 level_( travLevel ),
352 subEntity_( (codim == 0 ? 0 : -1) ),
353 macroIterator_( grid.meshPointer().begin() ),
354 marker_( marker )
355 {
356 ElementInfo elementInfo = *macroIterator_;
357 nextElementStop( elementInfo );
358 if( codim > 0 )
359 goNext( elementInfo );
360 // it is ok to set the invalid ElementInfo
361 entity_.impl().setElement( elementInfo, subEntity_ );
362 }
363
364
365 // Make LevelIterator with point to element from previous iterations
366 template< int codim, class GridImp, bool leafIterator >
368 ::AlbertaGridTreeIterator ( const GridImp &grid,
369 int travLevel )
370 : entity_( EntityImp( grid ) ),
371 level_( travLevel ),
372 subEntity_( -1 ),
373 macroIterator_( grid.meshPointer().end() ),
374 marker_( 0 )
375 {}
376
377
378 // Make LevelIterator with point to element from previous iterations
379 template< int codim, class GridImp, bool leafIterator >
382 : entity_( other.entity_ ),
383 level_( other.level_ ),
384 subEntity_( other.subEntity_ ),
385 macroIterator_( other.macroIterator_ ),
386 marker_( other.marker_ )
387 {}
388
389
390 // Make LevelIterator with point to element from previous iterations
391 template< int codim, class GridImp, bool leafIterator >
394 {
395 entity_ = other.entity_;
396 level_ = other.level_;
397 subEntity_ = other.subEntity_;
398 macroIterator_ = other.macroIterator_;
399 marker_ = other.marker_;
400
401 return *this;
402 }
403
404
405 template< int codim, class GridImp, bool leafIterator >
407 {
408 ElementInfo elementInfo = entity_.impl().elementInfo_;
409 goNext ( elementInfo );
410 // it is ok to set the invalid ElementInfo
411 entity_.impl().setElement( elementInfo, subEntity_ );
412 }
413
414
415 template< int codim, class GridImp, bool leafIterator >
417 ::nextElement ( ElementInfo &elementInfo )
418 {
419 if( elementInfo.isLeaf() || (elementInfo.level() >= level_) )
420 {
421 while( (elementInfo.level() > 0) && (elementInfo.indexInFather() == 1) )
422 elementInfo = elementInfo.father();
423 if( elementInfo.level() == 0 )
424 {
425 ++macroIterator_;
426 elementInfo = *macroIterator_;
427 }
428 else
429 elementInfo = elementInfo.father().child( 1 );
430 }
431 else
432 elementInfo = elementInfo.child( 0 );
433 }
434
435
436 template< int codim, class GridImp, bool leafIterator >
437 inline void AlbertaGridTreeIterator< codim, GridImp, leafIterator >
438 ::nextElementStop ( ElementInfo &elementInfo )
439 {
440 while( !(!elementInfo || stopAtElement( elementInfo )) )
441 nextElement( elementInfo );
442 }
443
444
445 template< int codim, class GridImp, bool leafIterator >
446 inline bool AlbertaGridTreeIterator< codim, GridImp, leafIterator >
447 ::stopAtElement ( const ElementInfo &elementInfo ) const
448 {
449 if( !elementInfo )
450 return true;
451 return (leafIterator ? elementInfo.isLeaf() : (level_ == elementInfo.level()));
452 }
453
454
455 template< int codim, class GridImp, bool leafIterator >
456 inline void AlbertaGridTreeIterator< codim, GridImp, leafIterator >
457 ::goNext ( ElementInfo &elementInfo )
458 {
459 std::integral_constant< int, codim > codimVariable;
460 goNext( codimVariable, elementInfo );
461 }
462
463 template< int codim, class GridImp, bool leafIterator >
464 inline void AlbertaGridTreeIterator< codim, GridImp, leafIterator >
465 ::goNext ( const std::integral_constant< int, 0 > cdVariable,
466 ElementInfo &elementInfo )
467 {
468 assert( stopAtElement( elementInfo ) );
469
470 nextElement( elementInfo );
471 nextElementStop( elementInfo );
472 }
473
474 template< int codim, class GridImp, bool leafIterator >
475 inline void AlbertaGridTreeIterator< codim, GridImp, leafIterator >
476 ::goNext ( const std::integral_constant< int, 1 > cdVariable,
477 ElementInfo &elementInfo )
478 {
479 assert( stopAtElement( elementInfo ) );
480
481 ++subEntity_;
482 if( subEntity_ >= numSubEntities )
483 {
484 subEntity_ = 0;
485 nextElement( elementInfo );
486 nextElementStop( elementInfo );
487 if( !elementInfo )
488 return;
489 }
490
491 if( leafIterator )
492 {
493 const int face = (dimension == 1 ? (numSubEntities-1)-subEntity_ : subEntity_);
494
495 const ALBERTA EL *neighbor = elementInfo.elInfo().neigh[ face ];
496 if( (neighbor != NULL) && !elementInfo.isBoundary( face ) )
497 {
498 // face is reached from element with largest number
499 const int elIndex = grid().dofNumbering() ( elementInfo, 0, 0 );
500 const int nbIndex = grid().dofNumbering() ( neighbor, 0, 0 );
501 if( elIndex < nbIndex )
502 goNext( cdVariable, elementInfo );
503 }
504 // uncomment this assertion only if codimension 1 entities are marked
505 // assert( marker_->template subEntityOnElement< 1 >( elementInfo, subEntity_ ) );
506 }
507 else
508 {
509 assert( marker_ != 0 );
510 if( !marker_->template subEntityOnElement< 1 >( elementInfo, subEntity_ ) )
511 goNext( cdVariable, elementInfo );
512 }
513 }
514
515 template< int codim, class GridImp, bool leafIterator >
516 template< int cd >
517 inline void AlbertaGridTreeIterator< codim, GridImp, leafIterator >
518 ::goNext ( const std::integral_constant< int, cd > cdVariable,
519 ElementInfo &elementInfo )
520 {
521 assert( stopAtElement( elementInfo ) );
522
523 ++subEntity_;
524 if( subEntity_ >= numSubEntities )
525 {
526 subEntity_ = 0;
527 nextElement( elementInfo );
528 nextElementStop( elementInfo );
529 if( !elementInfo )
530 return;
531 }
532
533 assert( marker_ != 0 );
534 if( !marker_->template subEntityOnElement< cd >( elementInfo, subEntity_ ) )
535 goNext( cdVariable, elementInfo );
536 }
537
538}
539
540#endif // #if HAVE_ALBERTA
541
542#endif // #ifndef DUNE_ALBERTA_TREEITERATOR_HH
Definition: treeiterator.hh:184
bool equals(const This &other) const
equality
Definition: treeiterator.hh:226
AlbertaGridTreeIterator(const GridImp &grid, const MarkerVector *marker, int travLevel)
Constructor making begin iterator.
Definition: treeiterator.hh:347
AlbertaGridTreeIterator(const GridImp &grid, int travLevel)
Constructor making end iterator.
Definition: treeiterator.hh:368
This & operator=(const This &other)
Constructor making end iterator.
Definition: treeiterator.hh:393
Entity & dereference() const
dereferencing
Definition: treeiterator.hh:232
int level() const
ask for level of entities
Definition: treeiterator.hh:238
void increment()
increment
Definition: treeiterator.hh:406
const GridImp & grid() const
obtain a reference to the grid
Definition: treeiterator.hh:248
AlbertaGridTreeIterator(const This &other)
Constructor making end iterator.
Definition: treeiterator.hh:381
marker assigning subentities to one element containing them
Definition: treeiterator.hh:32
AlbertaMarkerVector(const DofNumbering &dofNumbering)
create AlbertaMarkerVector with empty vectors
Definition: treeiterator.hh:51
bool up2Date() const
return true if marking is up to date
Definition: treeiterator.hh:92
bool subEntityOnElement(const ElementInfo &elementInfo, int subEntity) const
visit subentity on this element?
Definition: treeiterator.hh:289
void print(std::ostream &out=std::cout) const
print for debugin' only
Definition: treeiterator.hh:314
Wrapper class for entities.
Definition: entity.hh:64
Traits for type conversions and type information.
provides a wrapper for ALBERTA's el_info structure
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:46
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:267
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
provides a wrapper for ALBERTA's mesh structure
Dune namespace.
Definition: alignedallocator.hh:14
Static tag representing a codimension.
Definition: dimension.hh:22
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)