DUNE-FEM (unstable)

adaptiveleafindexset.hh
1#ifndef DUNE_FEM_ADAPTIVELEAFINDEXSET_HH
2#define DUNE_FEM_ADAPTIVELEAFINDEXSET_HH
3
4#include <cassert>
5#include <cstddef>
6
7#include <algorithm>
8#include <string>
9#include <vector>
10#include <memory>
11#include <type_traits>
12
13#include <dune/fem/common/forloop.hh>
14
15#include <dune/fem/gridpart/codimindexset.hh>
16#include <dune/fem/gridpart/common/gridpart.hh>
17#include <dune/fem/gridpart/common/persistentindexset.hh>
18#include <dune/fem/io/file/iointerface.hh>
19#include <dune/fem/version.hh>
20
21namespace Dune
22{
23
24 namespace Fem
25 {
26
27 // Internal forward declations
28 // ---------------------------
29
30 template < class GridPartImp >
31 class AdaptiveLeafIndexSet;
32 template < class GridPartImp >
33 class IntersectionAdaptiveLeafIndexSet;
34 template < class GridPartImp >
35 class DGAdaptiveLeafIndexSet;
36
37
38
39
41 //
42 // --AdaptiveIndexSetBaseTraits
43 //
45
46 template< class GridPart, class IndexSet >
47 class AdaptiveIndexSetBaseTraits
48 {
49 public:
50 // index set type derived from AdaptiveIndexSetBase
51 typedef IndexSet IndexSetType;
52
53 // grid part type
54 typedef GridPart GridPartType;
55 // grid type
56 typedef typename GridPartType :: GridType GridType;
57
58 // grid (part) dimension
59 static const int dimension = GridPartType :: dimension;
60
61 template< int codim >
62 struct Codim
63 {
64 // entity type
65 typedef typename GridPartType :: template Codim< codim > :: EntityType Entity;
66 };
67
68 // type of codim index set
69 typedef CodimIndexSet< GridType > CodimIndexSetType;
70 // index type used
71 typedef typename CodimIndexSetType :: IndexType IndexType;
72
73 // container of geometry types
74 typedef std::vector< GeometryType > Types;
75 };
76
77
78
80 //
81 // --AdaptiveIndexSetBase
82 //
84
96 template <class TraitsImp >
98 : public PersistentAdaptiveIndexSet< TraitsImp >
99 {
101 typedef PersistentAdaptiveIndexSet< TraitsImp > BaseType;
102
103 protected:
104 typedef typename TraitsImp :: GridPartType GridPartType;
105 typedef typename GridPartType :: GridType GridType;
106
107 typedef typename TraitsImp :: CodimIndexSetType CodimIndexSetType ;
108
109 typedef typename GridType::template Codim< 0 >::Entity GridElementType;
110
111 public:
113 static const int dimension = BaseType::dimension;
114
116 static const int numCodimensions = TraitsImp :: numCodimensions ;
117
119 static const int intersectionCodimension = TraitsImp :: intersectionCodimension ;
120
123
125 typedef typename BaseType :: IndexType IndexType;
126
128 typedef typename BaseType :: Types Types;
129
131 typedef typename BaseType :: template Codim< 0 > :: Entity ElementType;
132
134 typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
135
137 typedef typename GridPartType :: IntersectionType IntersectionType;
138
139
140 private:
141
142 template< int codim , bool gridHasCodim >
143 struct CountElementsBase
144 {
145 static void apply ( const ThisType &indexSet, const GeometryType &type, IndexType& count )
146 {
147 if( type.dim() == dimension - codim )
148 count = indexSet.template countElements< codim >( type, std::integral_constant<bool,true>() );
149 }
150 };
151
152 template< int codim >
153 struct CountElementsBase< codim, false >
154 {
155 static void apply ( const ThisType &indexSet, const GeometryType &type, IndexType& count )
156 {
157 if( type.dim() == dimension - codim )
158 count = indexSet.template countElements< codim >( type, std::integral_constant<bool,false>() );
159 }
160 };
161
162 template< int codim >
163 struct CountElements // TODO: hasEntity should be replaced by hasEntityIterator here
164 : public CountElementsBase< codim, Dune::Fem::GridPartCapabilities::hasEntity< GridPartType, codim > :: v >
165 {
166 };
167
168
169 template< int codim >
170 struct InsertSubEntities
171 {
172 static void apply ( const ThisType &indexSet, const GridElementType &element )
173 {
174 // if codimension is not available return
175 if( ! indexSet.codimAvailable( codim ) ) return ;
176
177 // if codimension is not used return
178 if( !indexSet.codimUsed( codim ) ) return;
179
180 CodimIndexSetType &codimSet = indexSet.codimLeafSet( codim );
181
182 const int count = element.subEntities( codim );
183 for( int i = 0; i < count; ++i )
184 {
185 codimSet.insertSubEntity( element, i );
186 }
187 }
188 };
189
190 template< int codim , bool gridHasCodim >
191 struct InsertGhostSubEntitiesBase
192 {
193 static void apply ( ThisType &indexSet, const GridElementType &element,
194 const bool skipGhosts )
195 {
196 // if codimension is not available return
197 if( ! indexSet.codimAvailable( codim ) ) return ;
198
199 // if codimension is not used return
200 if( !indexSet.codimUsed( codim ) ) return;
201
202 CodimIndexSetType &codimSet = indexSet.codimLeafSet( codim );
203
204 for( unsigned int i = 0; i < element.subEntities( codim ); ++i )
205 {
206 if( !skipGhosts || (element.partitionType() != GhostEntity) )
207 codimSet.insertGhost( element.template subEntity< codim >( i ) );
208 }
209 }
210 };
211
212 template< int codim >
213 struct InsertGhostSubEntitiesBase< codim, false >
214 {
215 static void apply ( ThisType &indexSet, const GridElementType &element,
216 const bool skipGhosts )
217 {}
218 };
219
220 template< int codim >
221 struct InsertGhostSubEntities
222 : public InsertGhostSubEntitiesBase< codim, Dune::Capabilities::hasEntity < GridType, codim > :: v >
223 {
224 };
225
226 template< int codim , bool gridHasCodim >
227 struct CallSetUpCodimSetBase
228 {
229 static void apply ( const int cd, const ThisType &indexSet )
230 {
231 // if codimension is not available return
232 if( ! indexSet.codimAvailable( codim ) ) return ;
233
234 if( cd == codim && ! indexSet.codimUsed( cd ) )
235 indexSet.template setupCodimSet< codim >(std::integral_constant<bool,true>());
236 }
237 };
238
239 template< int codim >
240 struct CallSetUpCodimSetBase< codim, false >
241 {
242 static void apply ( const int cd, const ThisType &indexSet )
243 {
244 if( cd == codim && ! indexSet.codimUsed( cd ) )
245 indexSet.template setupCodimSet< codim >(std::integral_constant<bool,false>());
246 }
247 };
248
249 template< int codim >
250 struct CallSetUpCodimSet
251 : public CallSetUpCodimSetBase< codim, Dune::Capabilities::hasEntity < GridType, codim > :: v >
252 {
253 };
254
255
257 // subentity extractor
259
260 template < int codim, bool gridHasCodim >
261 struct GetSubEntityBase
262 {
263 typedef typename GridPartType :: template Codim< codim > :: EntityType Entity;
264 static Entity subEntity( const ElementType& element, const int subEn )
265 {
266 return element.template subEntity< codim > ( subEn );
267 }
268 };
269
270 template < int codim >
271 struct GetSubEntityBase< codim, false >
272 {
273 typedef typename GridPartType :: template Codim< 0 > :: EntityType Entity;
274 static Entity subEntity( const ElementType& element, const int subEn )
275 {
276 DUNE_THROW(NotImplemented,"stupid grid without entities of codim 1 used");
277 }
278 };
279
280 struct GetFaceEntity
281 : public GetSubEntityBase< 1, Dune::Capabilities::hasEntity < GridType, 1 > :: v >
282 {
283 };
284
286 typedef typename GetFaceEntity :: Entity FaceType;
287
289 enum { CartesianNonAdaptiveGrid = Dune::Capabilities::isCartesian<GridType>::v &&
290 ! Capabilities::isLocallyAdaptive<GridType>::v };
291
292 // my type, to be revised
293 enum { myType = ( numCodimensions == 1 ) ? ( (CartesianNonAdaptiveGrid) ? -1 : 665 ) : 6 };
294
295 // max num of codimension (to avoid compiler warnings)
296 enum { maxNumCodimension = ((dimension + 1) > numCodimensions) ? dimension + 2 : numCodimensions+1 };
297
299 static const PartitionIteratorType pitype = GridPartType :: indexSetPartitionType ;
300
301 // grid part (for iterator access, no index set)
302 const GridPartType& gridPart_;
303
304 // pointer storage in case index set is created by passing a grid
305 std::unique_ptr< const GridPartType > gridPartPtr_;
306
307 // Codimension leaf index sets
308 mutable std::unique_ptr< CodimIndexSetType > codimLeafSet_[ numCodimensions ];
309
310 // vector holding geometry types
311 std::vector< std::vector< GeometryType > > geomTypes_;
312
313 // actual sequence number
314 int sequence_;
315
317 mutable bool compressed_;
318
319 protected:
320 using BaseType::grid_;
321 using BaseType::dofManager_;
322
323 // return true if codim is supported
324 bool codimAvailable( const int codim ) const
325 {
326 return (codim < numCodimensions && codim >= 0);
327 }
328
329 // return true if indices for this codim exist
330 bool codimUsed( const int codim ) const
331 {
332 return codimAvailable( codim ) && codimLeafSet_[ codim ] ;
333 }
334
335 CodimIndexSetType& codimLeafSet( const int codim ) const
336 {
337 assert( codimUsed( codim ) );
338 return *codimLeafSet_[ codim ];
339 }
340
341 public:
342 void requestCodimensions ( const std::vector< int >& codimensions ) const
343 {
344 // enable requested codimensions and rebuild index set
345 for( const auto& codim : codimensions )
346 {
347 // start loop from 1 since codim 0 is always created
348 Fem::ForLoop< CallSetUpCodimSet, 1, dimension >::apply( codim, *this );
349 }
350 }
351
353 AdaptiveIndexSetBase ( const GridType* grid )
354 : AdaptiveIndexSetBase( *(new GridPartType( const_cast< GridType& > (*grid), typename GridPartType::NoIndexSetType() )) )
355 {
356 // store pointer to avoid memory leaks
357 gridPartPtr_.reset( &gridPart_ );
358 }
359
361 AdaptiveIndexSetBase ( const GridPartType& gridPart )
362 : BaseType( gridPart.grid() )
363 , gridPart_( gridPart )
364 , sequence_( dofManager_.sequence() )
365 , compressed_(true) // at start the set is compressed
366 {
367 // codim 0 is used by default
368 codimLeafSet_[ 0 ].reset( new CodimIndexSetType( grid_, 0 ) );
369
370 // set the codim of each Codim Set.
371 for(int codim = 0; codim < numCodimensions; ++codim )
372 {
373 if( codim == intersectionCodimension )
374 codimLeafSet_[ codim ].reset( new CodimIndexSetType( grid_, 1 ) );
375 }
376
377 {
378 // get level-0 view, this is already used in GridPtr (DFG parser)
379 auto macroView = grid_.levelGridView( 0 );
380 setupGeomTypes( macroView.indexSet() );
381 }
382
383 // build index set
385 }
386
387 protected:
388 template <class MacroIndexSet>
389 void setupGeomTypes(const MacroIndexSet& indexSet)
390 {
392 // resize vector of geometry types
393 geomTypes_.resize( dimension+1 );
394 for(int codim=0; codim <= dimension; ++codim )
395 {
396 const int size = indexSet.types( codim ).size();
397 // copy geometry types
398 geomTypes_[ codim ].resize( size );
399 std::copy_n( indexSet.types( codim ).begin(), size, geomTypes_[ codim ].begin() );
400 }
401 }
402
403 public:
405 int type () const
406 {
407 return myType;
408 }
409
411 virtual std::string name () const
412 {
413 return "AdaptiveIndexSetBase";
414 }
415
416 //****************************************************************
417 //
418 // INTERFACE METHODS for DUNE INDEX SETS
419 // --size
420 //
421 //****************************************************************
424 {
425 const int codim = dimension - type.dim();
426
427 // true if only one geometry type is present
428 const bool onlySingleGeometryType = hasSingleGeometryType || ( geomTypes( codim ).size() == 1 ) ;
429 // use size of codim index set if possible
430 if( codimAvailable( codim ) && onlySingleGeometryType )
431 {
432 if( codimUsed( codim ) )
433 return type == geomTypes( codim )[ 0 ] ? codimLeafSet( codim ).size() : 0;
434 }
435
436 // count entities for given geometry type
437 IndexType count = 0 ;
438 Fem::ForLoop< CountElements, 0, dimension > :: apply( *this, type, count );
439 return count;
440 }
441
443 IndexType size ( int codim ) const
444 {
445 assert( codim < numCodimensions );
447 {
448 return codimLeafSet( codim ).size();
449 }
450
451 // count size for all geometry types
452 IndexType count = 0 ;
453 const size_t types = geomTypes( codim ).size();
454 for( size_t i=0; i<types; ++i )
455 {
456 count += size( geomTypes( codim )[ i ] );
457 }
458 return count ;
459 }
460
462 const std::vector <GeometryType> & geomTypes (const int codim) const
463 {
464 assert( codim >= 0 && codim < int(geomTypes_.size()) );
465 return geomTypes_[ codim ];
466 }
467
469 Types types( const int codim ) const
470 {
471 return geomTypes( codim );
472 }
473
475 template <class EntityType>
476 bool contains (const EntityType & en) const
477 {
478 enum { codim = EntityType::codimension };
479 if( codimAvailable( codim ) )
480 {
481 assert( codimUsed( codim ) );
482 return codimLeafSet( codim ).exists( gridEntity( en ) );
483 }
484 else
485 return false;
486 }
487
488 //****************************************************************
489 //
490 // METHODS for Adaptation with DofManger
491 //
492 //****************************************************************
493
495 void insertEntity( const GridElementType &entity )
496 {
497 // here we have to add the support of higher codims
499 insertIndex( entity );
500 }
501
503 void removeEntity( const GridElementType &entity )
504 {
505 removeIndex( entity );
506 }
507
510
512 void resize ()
513 {
515
516 #if HAVE_MPI
517 if( CartesianNonAdaptiveGrid &&
518 grid_.comm().size() > 1 )
519 {
520 // only done for structured grids
521 clear();
522
523 // this should only be the case of YaspGrid
524 markAllBelowOld<Interior_Partition>();
525 if( pitype > Interior_Partition )
526 {
527 markAllBelowOld< pitype >();
528 }
529 compressed_ = true;
530 }
531 else
532 #endif
533 {
534 // use a hierarchic walk to mark new elements
535 markAllBelowOld< pitype > ();
536
537 #if HAVE_MPI
538 // only if ghost are really supported
539 if( pitype == All_Partition )
540 {
541 if( grid_.comm().size() > 1 )
542 {
543 // make sure that also ghosts have indices
544 markAllUsed<Ghost_Partition>();
545 }
546 }
547 #endif
548 }
549 }
550
552 bool compress ();
553
554 public:
556 // index methods
557 // --index
560 template< class Entity >
561 IndexType index ( const Entity &entity ) const
562 {
563 return index< Entity :: codimension >( entity );
564 }
565
567 template< int codim >
569 index ( const typename GridPartType::template Codim< codim >::EntityType &entity ) const
570 {
571 if( codimAvailable( codim ) )
572 {
573 if( (codim != 0) && ! codimUsed( codim ) )
574 setupCodimSet< codim >(std::integral_constant<bool,Dune::Capabilities::hasEntity < GridType, codim > :: v>());
575
576 return codimLeafSet( codim ).index( gridEntity( entity ) );
577 }
578 else
579 {
580 DUNE_THROW( NotImplemented, (name() + " does not support indices for codim = ") << codim );
581 return -1;
582 }
583 }
584
585 /* \brief return index for intersection */
586 IndexType index ( const IntersectionType &intersection ) const
587 {
588 enum { codim = intersectionCodimension };
589 if( codimAvailable( codim ) )
590 {
591 // this in only done on first call
592 setupIntersections();
593
594 // get corresponding face entity pointer
595 FaceType face = getIntersectionFace( intersection );
596
597 return codimLeafSet( codim ).index( gridEntity( face ) );
598 }
599 else
600 {
601 DUNE_THROW( NotImplemented, (name() + " does not support indices for intersections, intersectionCodim = ") << codim );
602 return -1;
603 }
604 }
605
606 /* \brief return index for sub entity of given intersection and subEntityNumber */
608 subIndex ( const IntersectionType &intersection,
609 int subNumber, unsigned int codim ) const
610 {
611 DUNE_THROW( NotImplemented, (name() + " does not support subIndices for intersections, intersectionCodim = ") << codim );
612 return -1;
613 }
614
616 template< class Entity >
617 IndexType subIndex ( const Entity &entity, int subNumber, unsigned int codim ) const
618 {
619 return subIndex< Entity::codimension >( entity, subNumber, codim );
620 }
621
623 template< int cd >
624 IndexType subIndex ( const typename GridPartType::template Codim< cd >::EntityType &entity,
625 int subNumber, unsigned int codim ) const
626 {
627 assert( (int( codim ) >= cd) && (int( codim ) <= dimension) );
628 if( !codimAvailable( codim ) )
629 DUNE_THROW( NotImplemented, (name() + " does not support indices for codim = ") << codim );
630
631 if( (codim != 0) && ! codimUsed( codim ) )
632 {
633 Fem::ForLoop< CallSetUpCodimSet, 1, dimension >::apply( codim, *this );
634 }
635
636 const CodimIndexSetType &codimSet = codimLeafSet( codim );
637 const IndexType idx = codimSet.subIndex( gridEntity( entity ), subNumber );
638 assert( (idx >= 0) && (idx < IndexType( codimSet.size() )) );
639 return idx;
640 }
641
643 //
644 // DoF adjustment methods, needed by AdaptiveDofMapper interface
645 //
647
650 {
651 const int codim = dimension - type.dim();
652 // this index set works only with single geometry types adaptively
653 assert( hasSingleGeometryType || geomTypes( codim ).size() == 1 );
654
655 // if mapper is asking for types that are
656 // not in the index set then 0 should be returned
657 if( geomTypes( codim )[ 0 ] != type )
658 {
659 return 0;
660 }
661
662 return numberOfHoles( codim );
663 }
664
666 IndexType numberOfHoles ( const int codim ) const
667 {
668 if( codimAvailable( codim ) && codimUsed( codim ) )
669 {
670 return codimLeafSet( codim ).numberOfHoles();
671 }
672 else
673 return 0;
674 }
675
678 {
679 const int codim = dimension - type.dim();
680 assert( hasSingleGeometryType || geomTypes( codim ).size() == 1 );
681 return oldIndex( hole, codim );
682 }
683
685 IndexType oldIndex (const IndexType hole, const int codim ) const
686 {
687 if( codimAvailable( codim ) )
688 {
689 assert( codimUsed( codim ) );
690 return codimLeafSet( codim ).oldIndex( hole );
691 }
692 else
693 {
694 DUNE_THROW( NotImplemented, (name() + " does not support indices for codim = ") << codim );
695 return IndexType(-1);
696 }
697 }
698
701 {
702 const int codim = dimension - type.dim();
703 assert( hasSingleGeometryType || geomTypes( codim ).size() == 1 );
704 return newIndex( hole, codim );
705 }
706
708 IndexType newIndex (const IndexType hole , const int codim ) const
709 {
710 if( codimAvailable( codim ) )
711 {
712 assert( codimUsed( codim ) );
713 return codimLeafSet( codim ).newIndex( hole );
714 }
715 else
716 {
717 DUNE_THROW( NotImplemented, (name() + " does not support indices for codim = ") << codim );
718 return IndexType(-1);
719 }
720 }
721
722 protected:
723 // Note: The following methods forward to Dune::Fem::CodimIndexSet, which
724 // expects a Dune::Grid as template argument; all arguments passed to
725 // members of Dune::Fem::CodimIndexSet must be compatible with the
726 // template grid type. If entities returned by the grid and the grid part
727 // respectively differ in type, Dune::Fem::AdaptiveLeafIndexSetBase will
728 // call the necessary operations from grid part entities to grid entites.
729
730 // memorise index
731 void insertIndex ( const GridElementType &entity );
732
733 // memorise indices for all intersections
734 void insertIntersections ( const GridElementType &entity ) const;
735
736 // insert index temporarily
737 void insertTemporary ( const GridElementType &entity );
738
739 // set indices to unsed so that they are cleaned on compress
740 void removeIndex ( const GridElementType &entity );
741
742 // check whether entity can be inserted or not
743 void checkHierarchy ( const GridElementType &entity, bool wasNew );
744
745 // mark indices that are still used (and give new indices to new elements)
746 template <PartitionIteratorType pt>
747 void markAllUsed ();
748
750 void clear();
751
754
755 // give all entities that lie below the old entities new numbers
756 // here we need the hierarchic iterator because for example for some
757 // grid more the one level of new elements can be created during adaption
758 // there for we start to give new number for all elements below the old
759 // element
760 template <PartitionIteratorType pt>
761 void markAllBelowOld ();
762
763 // mark indices that are still used (and give new indices to new elements)
764 template< int codim >
765 void setupCodimSet (const std::integral_constant<bool,true> &hasEntities) const;
766 template< int codim >
767 void setupCodimSet (const std::integral_constant<bool,false> &hasEntities) const;
768
769 // mark indices that are still used (and give new indices to new intersections)
770 void setupIntersections () const;
771
772 // count elements by iterating over grid and compare
773 // entities of given codim with given type
774 template< int codim >
775 inline IndexType countElements ( GeometryType type, const std::integral_constant<bool,true> &hasEntities ) const;
776 template< int codim >
777 inline IndexType countElements ( GeometryType type, const std::integral_constant<bool,false> &hasEntities ) const;
778
779 public:
781 template< class StreamTraits >
783
785 template< class StreamTraits >
787
788 protected:
789 FaceType getIntersectionFace( const IntersectionType& intersection ) const
790 {
791 ElementType inside = intersection.inside();
792 return getIntersectionFace( intersection, inside );
793 }
794
795 FaceType getIntersectionFace( const IntersectionType& intersection,
796 const ElementType& inside ) const
797 {
798 if( ! intersection.conforming() && intersection.neighbor() )
799 {
800 ElementType outside = intersection.outside();
801 // only if outside is more refined then inside
802 if( inside.level() < outside.level() )
803 return GetFaceEntity :: subEntity( outside, intersection.indexInOutside() );
804 }
805
806 // default: get subentity of inside
807 return GetFaceEntity :: subEntity( inside, intersection.indexInInside() );
808 }
809 };
810
811 template< class TraitsImp >
812 inline void
814 {
815 codimLeafSet( 0 ).resize();
816
817 // if more than one codimension is supported
818 if( numCodimensions > 1 )
819 {
820 for( int codim = 1; codim < numCodimensions; ++codim )
821 {
822 if( codimUsed( codim ) )
823 codimLeafSet( codim ).resize();
824 }
825 }
826 }
827
828
829 // --compress
830 template< class TraitsImp >
831 inline bool
833 {
834 // reset list of holes in any case
835 for( int codim = 0; codim < numCodimensions; ++codim )
836 {
837 if( codimUsed( codim ) )
838 codimLeafSet( codim ).clearHoles();
839 }
840
841 if( compressed_ )
842 {
843 // if set already compress, do noting for serial runs
844 // in parallel runs check sequence number of dof manager
845 if( (grid_.comm().size() == 1) || (sequence_ == dofManager_.sequence()) )
846 return false;
847 }
848
849 // prepare index sets for setup
850 for( int codim = 0; codim < numCodimensions; ++codim )
851 {
852 if( codimUsed( codim ) )
853 codimLeafSet( codim ).prepareCompress();
854 }
855
856 // mark all indices still needed
857 setupIndexSet();
858
859 // true if a least one index is moved
860 bool haveToCopy = codimLeafSet( 0 ).compress();
861 for( int codim = 1; codim < numCodimensions; ++codim )
862 {
863 if( codimUsed( codim ) )
864 haveToCopy |= codimLeafSet( codim ).compress();
865 }
866
867 // now status is compressed
868 compressed_ = true;
869 // update sequence number
870 sequence_ = dofManager_.sequence();
871
872 return haveToCopy;
873 }
874
875
876 template< class TraitsImp >
877 inline void
878 AdaptiveIndexSetBase< TraitsImp >::insertIndex ( const GridElementType &entity )
879 {
880#if HAVE_MPI
881 // we need special treatment for ghosts
882 // ghosts should not be inlcuded in holes list
883 if( entity.partitionType() == GhostEntity )
884 {
885 codimLeafSet( 0 ).insertGhost( entity );
886 const bool skipGhosts = (pitype != All_Partition);
887 // only for index sets upporting more than one codim
888 if( numCodimensions > 1 )
889 Fem::ForLoop< InsertGhostSubEntities, 1, dimension >::apply( *this, entity, skipGhosts );
890 }
891 else
892#endif // HAVE_MPI
893 {
894 codimLeafSet( 0 ).insert( entity );
895 // only for index sets supporting more than one codim
896 if( numCodimensions > 1 )
897 Fem::ForLoop< InsertSubEntities, 1, dimension >::apply( *this, entity );
898
899 }
900
901 assert( codimLeafSet( 0 ).exists( entity ) );
902
903 // insert intersections if this is enabled
904 if( intersectionCodimension > 0 )
905 {
906 insertIntersections( entity );
907 }
908
909 // now consecutivity is no longer guaranteed
910 compressed_ = false;
911 }
912
913 template< class TraitsImp >
914 inline void
915 AdaptiveIndexSetBase< TraitsImp >::insertIntersections ( const GridElementType &gridElement ) const
916 {
917 codimLeafSet( intersectionCodimension ).resize();
918
919 const ElementType &element = gridPart_.convert( gridElement );
920 const IntersectionIteratorType endiit = gridPart_.iend( element );
921 for( IntersectionIteratorType iit = gridPart_.ibegin( element );
922 iit != endiit ; ++ iit )
923 {
924 // get intersection
925 const IntersectionType& intersection = *iit ;
926
927 // get correct face pointer
928 FaceType face = getIntersectionFace( intersection, element );
929
930 // insert face into index set
931 codimLeafSet( intersectionCodimension ).insert( gridEntity( face ) );
932 }
933 }
934
935 template< class TraitsImp >
936 inline void
937 AdaptiveIndexSetBase< TraitsImp >::insertTemporary( const GridElementType &entity )
938 {
939 insertIndex( entity );
940 codimLeafSet( 0 ).markForRemoval( entity );
941 }
942
943 template< class TraitsImp >
944 inline void
945 AdaptiveIndexSetBase< TraitsImp >::removeIndex( const GridElementType &entity )
946 {
947 // remove entities (only mark them as unused)
948 codimLeafSet( 0 ).markForRemoval( entity );
949
950 // don't remove higher codim indices (will be done on compression
951
952 // now consecutivity is no longer guaranteed
953 compressed_ = false;
954 }
955
956
957 template< class TraitsImp >
958 inline void
959 AdaptiveIndexSetBase< TraitsImp >
960 ::checkHierarchy ( const GridElementType &entity, bool wasNew )
961 {
962 // for leaf entities, just insert the index
963 // this certainly the case for grids without hierarchy
964 if( entity.isLeaf() )
965 {
966 insertIndex( entity );
967 return;
968 }
969
970 bool isNew = wasNew ;
971 typedef typename GridType::HierarchicIterator HierarchicIterator;
972
973 if( isNew )
974 {
975 // this is a new entity, so insert it,
976 // but only temporarily because it's not a leaf entity
977 insertTemporary( entity );
978 }
979 else
980 {
981 // if we were a leaf entity, all children are new
982 isNew = codimLeafSet( 0 ).validIndex( entity );
983 }
984
985 // entity has children and we need to go deeper
986 const int childLevel = entity.level() + 1;
987 const HierarchicIterator end = entity.hend( childLevel );
988 for( HierarchicIterator it = entity.hbegin( childLevel ); it != end; ++it )
989 checkHierarchy( *it, isNew );
990 }
991
992
993 template< class TraitsImp >
994 template< PartitionIteratorType pt >
995 inline void
996 AdaptiveIndexSetBase< TraitsImp >::markAllUsed ()
997 {
998 // make correct size of vectors
999 resizeVectors();
1000
1001 // mark all indices as unused
1002 for( int codim = 0; codim < numCodimensions; ++codim )
1003 {
1004 if( codimUsed( codim ) )
1005 codimLeafSet( codim ).resetUsed();
1006 }
1007
1008 typedef typename GridPartType
1009 ::template Codim< 0 > :: template Partition< pt > :: IteratorType Iterator;
1010
1011 const Iterator end = gridPart_.template end< 0, pt >();
1012 for( Iterator it = gridPart_.template begin< 0, pt >(); it != end; ++it )
1013 insertIndex( gridEntity( *it ) );
1014 }
1015
1016 template< class TraitsImp >
1017 inline void
1019 {
1020 // for structured grids clear all information
1021 // this in only done when setting up grids or after
1022 // read of parallel data on serial grids
1023 if( CartesianNonAdaptiveGrid )
1024 {
1025 // mark all indices as unused
1026 for( int codim = 0; codim < numCodimensions; ++codim )
1027 {
1028 if( codimUsed( codim ) )
1029 {
1030 // clear all information
1031 codimLeafSet( codim ).clear();
1032 }
1033 }
1034 }
1035 }
1036
1037 template< class TraitsImp >
1038 inline void
1040 {
1041 // only done for structured grids
1042 clear();
1043
1044#if HAVE_MPI
1045 // for YaspGrid we need all interior indices first
1046 if( CartesianNonAdaptiveGrid &&
1047 grid_.comm().size() > 1 )
1048 {
1049 // we should only get here for YaspGrid
1050 markAllUsed<Interior_Partition> ();
1051 if( pitype > Interior_Partition )
1052 markAllUsed< pitype >();
1053 }
1054 else
1055#endif
1056 {
1057 // give all entities that lie on the leaf level new numbers
1058 markAllUsed< pitype > ();
1059 }
1060 }
1061
1062 template< class TraitsImp >
1063 template< PartitionIteratorType pt >
1064 inline void
1066 {
1067 // mark all indices as unused
1068 for( int codim = 0; codim < numCodimensions; ++codim )
1069 {
1070 if( codimUsed( codim ) )
1071 codimLeafSet( codim ).resetUsed();
1072 }
1073
1074 // get macro iterator
1075 auto macroView = grid_.levelGridView( 0 );
1076
1077 const auto macroend = macroView.template end< 0, pt >();
1078 for( auto macroit = macroView.template begin< 0, pt >();
1079 macroit != macroend; ++macroit )
1080 {
1081 checkHierarchy( *macroit, false );
1082 }
1083 }
1084
1085
1086 template< class TraitsImp >
1087 template< int codim >
1088 inline void
1089 AdaptiveIndexSetBase< TraitsImp >::setupCodimSet (const std::integral_constant<bool,true>&) const
1090 {
1091 // if codim is not available do nothing
1092 if( ! codimAvailable( codim ) ) return ;
1093
1094 // create codimLeafSet if not existing
1095 if( ! codimLeafSet_[ codim ] )
1096 {
1097 codimLeafSet_[ codim ].reset( new CodimIndexSetType( grid_, codim ) );
1098 }
1099
1100 // resize if necessary
1101 codimLeafSet( codim ).resize();
1102
1103 // walk over grid parts entity set and insert entities
1104 typedef typename GridPartType
1105 ::template Codim< codim >::template Partition< pitype > :: IteratorType Iterator;
1106
1107 const Iterator end = gridPart_.template end< codim, pitype >();
1108 for( Iterator it = gridPart_.template begin< codim, pitype >(); it != end; ++it )
1109 codimLeafSet( codim ).insert( gridEntity( *it ) );
1110 }
1111
1112 template< class TraitsImp >
1113 template< int codim >
1114 inline void
1115 AdaptiveIndexSetBase< TraitsImp >::setupCodimSet (const std::integral_constant<bool,false>&) const
1116 {
1117 // if codim is not available do nothing
1118 if( ! codimAvailable( codim ) ) return ;
1119
1120 // create codimLeafSet if not existing
1121 if( ! codimLeafSet_[ codim ] )
1122 {
1123 codimLeafSet_[ codim ].reset( new CodimIndexSetType( grid_, codim ) );
1124 }
1125
1126 // resize if necessary
1127 codimLeafSet( codim ).resize();
1128
1129 typedef typename GridPartType
1130 ::template Codim< 0 >::template Partition< pitype > :: IteratorType Iterator;
1131
1132 const Iterator end = gridPart_.template end< 0, pitype >();
1133 for( Iterator it = gridPart_.template begin< 0, pitype >(); it != end; ++it )
1134 {
1135 const ElementType& element = *it ;
1136 const GridElementType &gridElement = gridEntity( element );
1137 InsertSubEntities< codim >::apply( *this, gridElement );
1138 }
1139 }
1140
1141
1142 template< class TraitsImp >
1143 inline void
1144 AdaptiveIndexSetBase< TraitsImp >::setupIntersections() const
1145 {
1146 // if intersectionCodimension < 0 then this feature is disabled
1147 if( intersectionCodimension < 0 ) return ;
1148
1149 // do nothing if insections are already available
1150 if( codimUsed( intersectionCodimension ) ) return ;
1151
1152 // resize if necessary
1153 codimLeafSet( intersectionCodimension ).resize();
1154
1155 // walk over grid parts entity set and insert entities
1156 typedef typename GridPartType
1157 ::template Codim< 0 >::template Partition< pitype > :: IteratorType Iterator;
1158
1159 const Iterator end = gridPart_.template end< 0, pitype >();
1160 for( Iterator it = gridPart_.template begin< 0, pitype >(); it != end; ++it )
1161 {
1162 // insert all intersections of this entity
1163 insertIntersections( gridEntity( *it ) );
1164 }
1165 }
1166
1167 template< class TraitsImp >
1168 template< int codim >
1170 AdaptiveIndexSetBase< TraitsImp >::countElements ( GeometryType type, const std::integral_constant<bool,true>& ) const
1171 {
1172 typedef typename GridPartType
1173 ::template Codim< codim > :: template Partition< pitype > :: IteratorType Iterator;
1174
1175 const Iterator begin = gridPart_.template begin< codim, pitype >();
1176 const Iterator end = gridPart_.template end< codim, pitype >();
1177 IndexType count = 0;
1178 for( Iterator it = begin; it != end; ++it )
1179 {
1180 if( it->type() == type )
1181 {
1182 ++count;
1183 }
1184 }
1185 return count;
1186 }
1187
1188 template< class TraitsImp >
1189 template< int codim >
1191 AdaptiveIndexSetBase< TraitsImp >::countElements ( GeometryType type, const std::integral_constant<bool,false>& ) const
1192 {
1193 if( ! codimLeafSet_[ codim ] )
1194 return 0;
1195
1196 // make sure codimension is enabled
1197 assert( codimAvailable( codim ) );
1198
1199 // resize if necessary
1200 codimLeafSet( codim ).resize();
1201
1202 typedef typename GridPartType
1203 ::template Codim< 0 >::template Partition< pitype > :: IteratorType Iterator;
1204
1205 typedef typename GridPartType::ctype ctype;
1206
1207 const Iterator end = gridPart_.template end< 0, pitype >();
1208 IndexType count = 0;
1209 for( Iterator it = gridPart_.template begin< 0, pitype >(); it != end; ++it )
1210 {
1211 const ElementType& element = *it ;
1212 const GridElementType &gridElement = gridEntity( element );
1213 const int subEntities = gridElement.subEntities( codim );
1214 for (int i=0; i < subEntities; ++i)
1215 {
1216 if (! codimLeafSet( codim ).exists( gridElement, i) )
1217 {
1218 codimLeafSet( codim ).insertSubEntity( gridElement,i );
1220 general( gridElement.type() ).type( i, codim ) == type )
1221 {
1222 ++count;
1223 }
1224 }
1225 }
1226 }
1227
1228 return count;
1229 }
1230
1231
1232 template< class TraitsImp >
1233 template< class StreamTraits >
1236 {
1237 // write name for identification
1238 const std::string myname( name() );
1239 out << myname;
1240
1241 // write number of codimensions
1242 out << numCodimensions ;
1243
1244 // write whether codim is used
1245 for( int i = 0; i < numCodimensions; ++i )
1246 out << codimUsed( i );
1247
1248 // write all sets
1249 for( int i = 0; i < numCodimensions; ++i )
1250 {
1251 if( codimUsed( i ) )
1252 codimLeafSet( i ).write( out );
1253 }
1254
1255 // if we got until here writing was sucessful
1256 return true;
1257 }
1258
1259
1260 template< class TraitsImp >
1261 template< class StreamTraits >
1264 {
1265 {
1266 // read name and check compatibility
1267 std::string storedName;
1268 in >> storedName;
1269
1270 std::string myname( name() );
1271
1272 if( myname != storedName )
1273 {
1274 size_t length = std::min( myname.size(), storedName.size() );
1275 // only print the first character of whatever storedName is
1276 std::string found = storedName.substr(0, length-1 );
1278 "AdaptiveIndexSetBase::read: got " << found
1279 << " (expected " << myname << ")." );
1280 }
1281 }
1282
1283 // read number of codimensions
1284 int numCodim;
1285 in >> numCodim;
1286
1287 // make sure everything is correct
1288 // assert( numCodim == numCodimensions );
1289 if( numCodim != numCodimensions )
1290 DUNE_THROW(InvalidStateException,"AdaptiveIndexSetBase::read: got wrong number of codimensions" << numCodim << " instead of " << numCodimensions);
1291
1292 // read codim used
1293 for( int i = 0; i < numCodimensions; ++i )
1294 {
1295 bool codimInUse = false ;
1296 in >> codimInUse;
1297 if( codimInUse && ! codimLeafSet_[ i ] )
1298 {
1299 codimLeafSet_[ i ].reset( new CodimIndexSetType( grid_, (i == intersectionCodimension ) ? 1 : i ) );
1300 }
1301 }
1302
1303 for( int i = 0; i < numCodimensions; ++i )
1304 {
1305 if( codimUsed( i ) )
1306 codimLeafSet( i ).read( in );
1307 }
1308
1309 // in parallel runs we have to compress here
1310 if( grid_.comm().size() > 1 )
1311 compressed_ = false;
1312
1313 // if we got until here reading was sucessful
1314 return true;
1315 }
1316
1317
1318
1320 //
1321 // --AdaptiveLeafIndexSet
1322 //
1324
1325 template< class GridPartImp >
1326 struct AdaptiveLeafIndexSetTraits
1327 : public AdaptiveIndexSetBaseTraits< GridPartImp, AdaptiveLeafIndexSet< GridPartImp > >
1328 {
1329 // number of codimensions
1330 enum { numCodimensions = GridPartImp :: dimension + 1 };
1331 // first comdimension that is supported (not yet supported)
1332 enum { startingCodimension = 0 };
1333 // intersection codimensions (this is usually dimension + 1 )
1334 enum { intersectionCodimension = -1 };
1335 };
1336
1348 template < class GridPartImp >
1350 : public AdaptiveIndexSetBase< AdaptiveLeafIndexSetTraits< GridPartImp > >
1351 {
1353 typedef AdaptiveLeafIndexSetTraits< GridPartImp > Traits;
1354
1355 public:
1356 typedef typename BaseType :: GridType GridType;
1357 typedef typename BaseType :: GridPartType GridPartType;
1358
1360 AdaptiveLeafIndexSet (const GridType* grid)
1361 : BaseType(grid)
1362 {}
1363
1365 AdaptiveLeafIndexSet (const GridPartType& gridPart)
1366 : BaseType(gridPart)
1367 {}
1368
1370 virtual std::string name () const
1371 {
1372 return "AdaptiveLeafIndexSet";
1373 }
1374
1375 bool compress ()
1376 {
1377 const bool compressed = BaseType::compress();
1378
1379#ifndef NDEBUG
1380 if( this->grid_.comm().size() == 1 )
1381 {
1382 for( int codim = Traits::startingCodimension; codim < Traits::numCodimensions; ++codim )
1383 assert( this->codimUsed( codim ) ?
1384 size_t(this->size( codim )) == size_t(this->grid_.size( codim )) : true );
1385 }
1386#endif // #ifndef NDEBUG
1387
1388 return compressed;
1389 }
1390 };
1391
1392
1394 //
1395 // --IntersectionAdaptiveLeafIndexSet
1396 //
1398
1399 template< class GridPartImp >
1400 struct IntersectionAdaptiveLeafIndexSetTraits
1401 : public AdaptiveIndexSetBaseTraits< GridPartImp, IntersectionAdaptiveLeafIndexSet< GridPartImp > >
1402 {
1403 // number of codimensions
1404 enum { numCodimensions = GridPartImp :: dimension + 2 };
1405 // intersection codimensions (this is usually dimension + 1 )
1406 enum { intersectionCodimension = numCodimensions - 1 };
1407 // first comdimension that is supported (not yet supported)
1408 enum { startingCodimension = 0 };
1409 };
1410
1422 template < class GridPartImp >
1423 class IntersectionAdaptiveLeafIndexSet
1424 : public AdaptiveIndexSetBase< IntersectionAdaptiveLeafIndexSetTraits< GridPartImp > >
1425 {
1426 typedef AdaptiveIndexSetBase< IntersectionAdaptiveLeafIndexSetTraits< GridPartImp > > BaseType;
1427 typedef IntersectionAdaptiveLeafIndexSetTraits< GridPartImp > Traits;
1428
1429 public:
1430 typedef typename BaseType :: GridType GridType;
1431 typedef typename BaseType :: GridPartType GridPartType;
1433 IntersectionAdaptiveLeafIndexSet (const GridType* grid)
1434 : BaseType(grid)
1435 {}
1436
1438 IntersectionAdaptiveLeafIndexSet (const GridPartType& gridPart)
1439 : BaseType(gridPart)
1440 {}
1441
1443 virtual std::string name () const
1444 {
1445 return "IntersectionAdaptiveLeafIndexSet";
1446 }
1447
1448 bool compress ()
1449 {
1450 const bool compressed = BaseType::compress();
1451
1452#ifndef NDEBUG
1453 if( this->grid_.comm().size() == 1 )
1454 {
1455 for( int codim = Traits::startingCodimension; codim < Traits::numCodimensions; ++codim )
1456 if( codim != Traits::intersectionCodimension )
1457 assert( size_t(this->size( codim )) == size_t(this->grid_.size( codim )) );
1458 }
1459#endif // #ifndef NDEBUG
1460
1461 return compressed;
1462 }
1463 };
1464
1466 //
1467 // --DGAdaptiveLeafIndexSet
1468 //
1470
1471 template< class GridPartImp >
1472 struct DGAdaptiveLeafIndexSetTraits
1473 : public AdaptiveIndexSetBaseTraits< GridPartImp, DGAdaptiveLeafIndexSet< GridPartImp > >
1474 {
1475 // this index set only supports one codimension, codim zero
1476 enum { numCodimensions = 1 };
1477 // first comdimension that is supported (not yet supported)
1478 enum { startingCodimension = 0 };
1479 // intersection codimensions (this is usually dimension + 1 )
1480 enum { intersectionCodimension = -1 };
1481 };
1482
1493 template < class GridPartImp >
1495 : public AdaptiveIndexSetBase< DGAdaptiveLeafIndexSetTraits< GridPartImp > >
1496 {
1498 typedef DGAdaptiveLeafIndexSetTraits< GridPartImp > Traits;
1499
1500 public:
1501 typedef typename BaseType :: GridType GridType;
1502 typedef typename BaseType :: GridPartType GridPartType;
1504 DGAdaptiveLeafIndexSet (const GridType* grid)
1505 : BaseType(grid)
1506 {}
1507
1508 DGAdaptiveLeafIndexSet (const GridPartType& gridPart)
1509 : BaseType(gridPart)
1510 {}
1511
1513 virtual std::string name () const
1514 {
1515 return "DGAdaptiveLeafIndexSet";
1516 }
1517
1518 bool compress ()
1519 {
1520 const bool compressed = BaseType::compress();
1521
1522#ifndef NDEBUG
1523 if( this->grid_.comm().size() == 1 )
1524 assert( size_t(this->size( 0 )) == size_t(this->grid_.size( 0 )) );
1525#endif // #ifndef NDEBUG
1526
1527 return compressed;
1528 }
1529 };
1530
1531 } // namespace Fem
1532
1533} // namespace Dune
1534
1535#endif // #ifndef DUNE_FEM_ADAPTIVELEAFINDEXSET_HH
Wrapper class for entities.
Definition: entity.hh:66
consecutive, persistent index set for the leaf level based on the grid's hierarchy index set
Definition: adaptiveleafindexset.hh:99
IndexType subIndex(const typename GridPartType::template Codim< cd >::EntityType &entity, int subNumber, unsigned int codim) const
return index for given subentity *‍/
Definition: adaptiveleafindexset.hh:624
IndexType numberOfHoles(GeometryType type) const
return number of holes for given type *‍/
Definition: adaptiveleafindexset.hh:649
AdaptiveIndexSetBase(const GridPartType &gridPart)
Constructor.
Definition: adaptiveleafindexset.hh:361
virtual std::string name() const
return name of index set
Definition: adaptiveleafindexset.hh:411
void insertEntity(const GridElementType &entity)
please doc me *‍/
Definition: adaptiveleafindexset.hh:495
IndexType index(const Entity &entity) const
return number of entities of given type *‍/
Definition: adaptiveleafindexset.hh:561
bool read(InStreamInterface< StreamTraits > &in)
please doc me *‍/
Definition: adaptiveleafindexset.hh:1263
BaseType::template Codim< 0 >::Entity ElementType
type of codimension 0 Entity
Definition: adaptiveleafindexset.hh:131
IndexType subIndex(const Entity &entity, int subNumber, unsigned int codim) const
return index for given subentity *‍/
Definition: adaptiveleafindexset.hh:617
IndexType index(const typename GridPartType::template Codim< codim >::EntityType &entity) const
return number of entities of given type *‍/
Definition: adaptiveleafindexset.hh:569
void setupIndexSet()
mark all indices of interest
Definition: adaptiveleafindexset.hh:1039
void resize()
please doc me *‍/
Definition: adaptiveleafindexset.hh:512
IndexType numberOfHoles(const int codim) const
return number of holes of the sets indices
Definition: adaptiveleafindexset.hh:666
Types types(const int codim) const
return range of geometry types *‍/
Definition: adaptiveleafindexset.hh:469
AdaptiveIndexSetBase(const GridType *grid)
Constructor.
Definition: adaptiveleafindexset.hh:353
void removeEntity(const GridElementType &entity)
please doc me *‍/
Definition: adaptiveleafindexset.hh:503
static const int intersectionCodimension
intersection codimension (numCodim-1 if enabled, otherwise -1)
Definition: adaptiveleafindexset.hh:119
IndexType oldIndex(const IndexType hole, const int codim) const
return old index, for dof manager only
Definition: adaptiveleafindexset.hh:685
IndexType newIndex(IndexType hole, GeometryType type) const
return new index for given hole and type *‍/
Definition: adaptiveleafindexset.hh:700
IndexType newIndex(const IndexType hole, const int codim) const
return new index, for dof manager only returns index
Definition: adaptiveleafindexset.hh:708
void resizeVectors()
reallocate the vector for new size
Definition: adaptiveleafindexset.hh:813
GridPartType::IntersectionType IntersectionType
type of intersections
Definition: adaptiveleafindexset.hh:137
IndexType oldIndex(IndexType hole, GeometryType type) const
return old index for given hole and type *‍/
Definition: adaptiveleafindexset.hh:677
static const int dimension
grid dimension *‍/
Definition: adaptiveleafindexset.hh:113
bool contains(const EntityType &en) const
return true if entity has index *‍/
Definition: adaptiveleafindexset.hh:476
bool write(OutStreamInterface< StreamTraits > &out) const
please doc me *‍/
Definition: adaptiveleafindexset.hh:1235
void clear()
clear index set (only for structured grids)
Definition: adaptiveleafindexset.hh:1018
BaseType::IndexType IndexType
index type *‍/
Definition: adaptiveleafindexset.hh:125
BaseType::Types Types
geometry type range type *‍/
Definition: adaptiveleafindexset.hh:128
static const int numCodimensions
number of supported codimensions
Definition: adaptiveleafindexset.hh:116
IndexType size(int codim) const
return number of entities of given type *‍/
Definition: adaptiveleafindexset.hh:443
const std::vector< GeometryType > & geomTypes(const int codim) const
*‍/
Definition: adaptiveleafindexset.hh:462
void setupGeomTypes(const MacroIndexSet &indexSet)
Definition: adaptiveleafindexset.hh:389
IndexType size(GeometryType type) const
return number of entities of given type *‍/
Definition: adaptiveleafindexset.hh:423
bool compress()
please doc me *‍/
Definition: adaptiveleafindexset.hh:832
int type() const
return type of index set, for GrapeDataIO
Definition: adaptiveleafindexset.hh:405
static const bool hasSingleGeometryType
true if only one geometry type is available
Definition: adaptiveleafindexset.hh:122
GridPartType::IntersectionIteratorType IntersectionIteratorType
type of intersection iterator
Definition: adaptiveleafindexset.hh:134
consecutive, persistent index set for the leaf level based on the grid's hierarchy index set
Definition: adaptiveleafindexset.hh:1351
virtual std::string name() const
return name of index set
Definition: adaptiveleafindexset.hh:1370
AdaptiveLeafIndexSet(const GridType *grid)
Constructor.
Definition: adaptiveleafindexset.hh:1360
AdaptiveLeafIndexSet(const GridPartType &gridPart)
Constructor.
Definition: adaptiveleafindexset.hh:1365
consecutive, persistent index set for the leaf level based on the grid's hierarchy index set
Definition: adaptiveleafindexset.hh:1496
DGAdaptiveLeafIndexSet(const GridType *grid)
Constructor.
Definition: adaptiveleafindexset.hh:1504
virtual std::string name() const
return name of index set
Definition: adaptiveleafindexset.hh:1513
abstract interface for an input stream
Definition: streams.hh:190
Traits::IndexType IndexType
index type
Definition: indexset.hh:136
abstract interface for an output stream
Definition: streams.hh:48
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
Default exception for dummy implementations.
Definition: exceptions.hh:263
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
PartitionIteratorType
Parameter to be used for the parallel level- and leaf iterators.
Definition: gridenums.hh:136
@ All_Partition
all entities
Definition: gridenums.hh:141
@ Interior_Partition
only interior entities
Definition: gridenums.hh:137
@ GhostEntity
ghost entities
Definition: gridenums.hh:35
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:58
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:27
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: capabilities.hh:48
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:128
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)