7 #include <dune/common/concept.hh> 8 #include <dune/functions/functionspacebases/concepts.hh> 10 #include <amdis/common/OptionalNoCopy.hpp> 11 #include <amdis/common/TypeTraits.hpp> 12 #include <amdis/functions/NodeCache.hpp> 13 #include <amdis/functions/Nodes.hpp> 23 using PrefixPath = Dune::TypeTree::HybridTreePath<>;
26 using NodeIndexSet = NodeIndexSet_t<typename GB::PreBasis, PrefixPath>;
36 using Element =
typename GridView::template Codim<0>::Entity;
42 using Tree = Node_t<typename GlobalBasis::PreBasis, PrefixPath>;
53 : globalBasis_(&globalBasis)
54 , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
55 , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
57 static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>,
Tree>());
58 initializeTree(tree_);
70 bindTree(tree_, element_);
71 indices_.resize(tree_.size());
72 bindNodeIndexSet(globalBasis_->preBasis(), nodeIndexSet_, tree_, indices_.begin(), Dune::PriorityTag<5>{});
79 static_assert(not std::is_rvalue_reference_v<Element&&>);
101 unbindNodeIndexSet(nodeIndexSet_, Dune::PriorityTag<5>{});
115 treeCache_.emplace(makeNodeCache(tree_));
133 return globalBasis_->preBasis().maxNodeSize();
147 return *globalBasis_;
157 template <
class PB,
class NIS,
class Tree,
class Iter,
158 class = decltype(std::declval<PB>().indices(std::declval<Tree>(), std::declval<Iter>()))>
159 static void bindNodeIndexSet(PB
const& preBasis, NIS& , Tree
const&
tree, Iter it, Dune::PriorityTag<3>)
161 preBasis.indices(tree, it);
164 template <
class PB,
class NIS,
class Tree,
class Iter,
165 class = decltype(std::declval<NIS>().
bind(std::declval<Tree>())),
166 class = decltype(std::declval<NIS>().indices(std::declval<Iter>()))>
167 static void bindNodeIndexSet(PB
const& , NIS& nodeIndexSet, Tree
const& tree, Iter it, Dune::PriorityTag<2>)
169 nodeIndexSet.bind(tree);
170 nodeIndexSet.indices(it);
173 template <
class PB,
class NIS,
class Tree,
class Iter,
174 class = decltype(std::declval<NIS>().
bind(std::declval<Tree>())),
175 class = decltype(std::declval<NIS>().
index(0u))>
176 static void bindNodeIndexSet(PB
const& , NIS& nodeIndexSet, Tree
const& tree, Iter it, Dune::PriorityTag<1>)
178 nodeIndexSet.bind(tree);
179 for (
size_type i = 0; i < tree.size(); ++i)
180 *it++ = nodeIndexSet.index(i);
183 template <
class PB,
class NIS,
class Tree,
class Iter>
184 static void bindNodeIndexSet(PB
const& , NIS& , Tree
const& , Iter , Dune::PriorityTag<0>)
190 class = decltype(std::declval<NIS>().
unbind())>
191 static void unbindNodeIndexSet(NIS& nodeIndexSet, Dune::PriorityTag<1>)
193 nodeIndexSet.unbind();
197 static void unbindNodeIndexSet(NIS& , Dune::PriorityTag<0>)
206 NodeIndexSet nodeIndexSet_;
209 std::vector<MultiIndex> indices_ = {};
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:121
constexpr bool MultiIndex
A multi-index type.
Definition: Concepts.hpp:150
size_type maxSize() const
Maximum local size for any element on the GridView.
Definition: LocalView.hpp:131
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:33
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:151
typename GlobalBasis::PreBasis::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:48
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:83
Definition: AdaptBase.hpp:6
MultiIndex index(size_type i) const
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis...
Definition: LocalView.hpp:138
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:36
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:21
Node_t< typename GlobalBasis::PreBasis, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: LocalView.hpp:42
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:112
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:106
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:145
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:89
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:67
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:62
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:99
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:52
GB GlobalBasis
The global FE basis that this is a view on.
Definition: LocalView.hpp:30
NodeCache_t< Tree > TreeCache
Cached basis-tree.
Definition: LocalView.hpp:45
std::size_t size_type
The type used for sizes.
Definition: LocalView.hpp:39