#ifndef DUNE_CURVEDGRID_GRIDVIEW_HH #define DUNE_CURVEDGRID_GRIDVIEW_HH #include <cassert> #include <utility> #include <dune/common/referencehelper.hh> #include <dune/curvedgrid/datahandle.hh> #include <dune/curvedgrid/declaration.hh> #include <dune/curvedgrid/indexsets.hh> #include <dune/curvedgrid/intersection.hh> #include <dune/curvedgrid/intersectioniterator.hh> #include <dune/curvedgrid/iterator.hh> #include <dune/grid/common/gridview.hh> namespace Dune { namespace Curved { // Internal Forward Declarations // ----------------------------- template <class HGV, class GF, bool useInterpolation> class GridView; // GridViewTraits // -------------- template <class HGV, class GF, bool useInterpolation> class GridViewTraits { friend class GridView<HGV, GF, useInterpolation>; using HostGridView = HGV; using GridFunction = ResolveReference_t<GF>; using HostGrid = typename HostGridView::Grid; using HostIntersection = typename HostGridView::Intersection; using HostIntersectionIterator = typename HostGridView::IntersectionIterator; public: using GridViewImp = GridView<HostGridView, GF, useInterpolation>; using Grid = Dune::CurvedGrid<HostGrid, GF, useInterpolation>; using IndexSet = Curved::IndexSet<const Grid, typename HostGridView::IndexSet>; using Intersection = Dune::Intersection<const Grid, Curved::Intersection<const Grid, HostIntersection>>; using IntersectionIterator = Dune::IntersectionIterator<const Grid, Curved::IntersectionIterator<const Grid, HostIntersectionIterator>, Curved::Intersection<const Grid, HostIntersection>>; using CollectiveCommunication = typename HostGridView::CollectiveCommunication; template <int codim> struct Codim { using IteratorImp = Curved::Iterator<HostGridView, codim, All_Partition, const Grid>; using Iterator = Dune::EntityIterator<codim, const Grid, IteratorImp>; using Entity = typename Grid::Traits::template Codim<codim>::Entity; using Geometry = typename Grid::template Codim<codim>::Geometry; using LocalGeometry = typename Grid::template Codim<codim>::LocalGeometry; template <PartitionIteratorType pit> struct Partition { using IteratorImp = Curved::Iterator<HostGridView, codim, pit, const Grid>; using Iterator = Dune::EntityIterator<codim, const Grid, IteratorImp>; }; }; static const bool conforming = HostGridView::conforming; }; // GridView // -------- template <class HGV, class GF, bool useInterpolation> class GridView { using Self = GridView; public: using Traits = GridViewTraits<HGV, GF, useInterpolation>; using HostGridView = typename Traits::HostGridView; using Grid = typename Traits::Grid; using GridFunction = ResolveReference_t<GF>; using IndexSet = typename Traits::IndexSet; using Intersection = typename Traits::Intersection; using IntersectionIterator = typename Traits::IntersectionIterator; using CollectiveCommunication = typename Traits::CollectiveCommunication; template <int codim> struct Codim : public Traits::template Codim<codim> {}; static const bool conforming = Traits::conforming; public: /// \brief construct the GridView from a hostGridView GridView (const Grid* grid, const HostGridView& hostGridView) : grid_(grid) , hostGridView_(hostGridView) {} /// \brief return a const reference to the CurvedGrid const Grid& grid () const { assert( !!grid_ ); return *grid_; } /// \brief return the geometry mapping of the grid const GridFunction& gridFunction () const { return grid().gridFunction(); } /// \brief obtain the index set /** * The lifetime of the returned index set is bound to the lifetime of the * grid view. Keep a copy of the grid view to prevent the index set from * becoming a dangling reference. */ const IndexSet& indexSet () const { if (!indexSet_) indexSet_.emplace(&hostGridView().indexSet()); return *indexSet_; } /// \brief obtain number of entities in a given codimension in the hostGridView. int size (int codim) const { return hostGridView().size(codim); } /// \brief obtain number of entities with a given geometry type in the hostGridView. int size (const GeometryType& type) const { return hostGridView().size(type); } /// \brief return whether the hostEntity associated to `e` is contained in the hostGridView. template<class EntityType> bool contains (const EntityType& e) const { return hostGridView().contains(e.impl().hostEntity()); } /// \brief obtain begin iterator for this view. template <int codim, PartitionIteratorType pit = All_Partition> typename Codim<codim>::template Partition<pit>::Iterator begin () const { using IteratorImp = typename Traits::template Codim<codim>::template Partition<pit>::IteratorImp; return IteratorImp::begin(gridFunction(), hostGridView(), grid().order()); } /// \brief obtain end iterator for this view. template <int codim, PartitionIteratorType pit = All_Partition> typename Codim<codim>::template Partition<pit>::Iterator end () const { using IteratorImp = typename Traits::template Codim<codim>::template Partition<pit>::IteratorImp; return IteratorImp::end(gridFunction(), hostGridView(), grid().order()); } /// \brief return begin intersection iterator with respect to this view IntersectionIterator ibegin (const typename Codim<0>::Entity& entity) const { using IteratorImpl = Curved::IntersectionIterator<const Grid, typename HostGridView::IntersectionIterator>; return IteratorImpl(entity, hostGridView().ibegin(entity.impl().hostEntity()), grid().order()); } /// \brief return end intersection iterator with respect to this view IntersectionIterator iend (const typename Codim<0>::Entity& entity) const { using IteratorImpl = Curved::IntersectionIterator<const Grid, typename HostGridView::IntersectionIterator>; return IteratorImpl(entity, hostGridView().iend(entity.impl().hostEntity()), grid().order()); } /// \brief obtain collective communication object of the hostGridView. const CollectiveCommunication& comm () const { return hostGridView().comm(); } ///\brief return the size of the overlap region of the hostGridView. int overlapSize (int codim) const { return hostGridView().overlapSize(codim); } ///\brief return the size of the ghost region of the hostGridView. int ghostSize (int codim) const { return hostGridView().ghostSize(codim); } /// \brief communicate via the hostGridView communication method using a wrapped dataHandle. template <class DataHandle, class Data> void communicate (CommDataHandleIF<DataHandle, Data>& dataHandle, InterfaceType interface, CommunicationDirection direction) const { using DataHandleIF = CommDataHandleIF<DataHandle, Data>; using WrappedDataHandle = Curved::CommDataHandle<Grid, DataHandleIF>; WrappedDataHandle wrappedDataHandle(grid(), dataHandle); hostGridView().communicate(wrappedDataHandle, interface, direction); } /// \brief return the GridView of the hostGrid. const HostGridView& hostGridView () const { return hostGridView_; } private: const Grid* grid_; HostGridView hostGridView_; mutable std::optional<IndexSet> indexSet_; }; } // end namespace Curved } // end namespace Dune #endif // DUNE_CURVEDGRID_GRIDVIEW_HH