DUNE PDELab (2.8)

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