7 #include <dune/common/concept.hh> 8 #include <dune/functions/functionspacebases/concepts.hh> 10 #include <amdis/common/TypeTraits.hpp> 11 #include <amdis/functions/NodeCache.hpp> 12 #include <amdis/functions/Nodes.hpp> 22 using PrefixPath = Dune::TypeTree::HybridTreePath<>;
25 using NodeIndexSet = NodeIndexSet_t<typename GB::PreBasis, PrefixPath>;
35 using Element =
typename GridView::template Codim<0>::Entity;
41 using Tree = Node_t<typename GlobalBasis::PreBasis, PrefixPath>;
52 : globalBasis_(&globalBasis)
53 , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
54 , treeCache_(makeNodeCache(tree_))
55 , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
57 static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>,
Tree>());
58 initializeTree(tree_);
65 : globalBasis_(other.globalBasis_)
66 , element_(other.element_)
68 , treeCache_(makeNodeCache(tree_))
69 , nodeIndexSet_(other.nodeIndexSet_)
70 , indices_(other.indices_)
77 : globalBasis_(other.globalBasis_)
78 , element_(other.element_)
79 , tree_(std::move(other.tree_))
80 , treeCache_(makeNodeCache(tree_))
81 , nodeIndexSet_(std::move(other.nodeIndexSet_))
82 , indices_(std::move(other.indices_))
94 bindTree(tree_, element);
95 indices_.resize(tree_.size());
96 bindNodeIndexSet(globalBasis_->preBasis(), nodeIndexSet_, tree_, indices_.begin(), Dune::PriorityTag<5>{});
102 static_assert(not std::is_rvalue_reference_v<Element&&>);
108 return element_ !=
nullptr;
124 unbindNodeIndexSet(nodeIndexSet_, Dune::PriorityTag<5>{});
153 return globalBasis_->preBasis().maxNodeSize();
167 return *globalBasis_;
177 template <
class PB,
class NIS,
class Tree,
class Iter,
178 class = decltype(std::declval<PB>().indices(std::declval<Tree>(), std::declval<Iter>()))>
179 static void bindNodeIndexSet(PB
const& preBasis, NIS& , Tree
const&
tree, Iter it, Dune::PriorityTag<3>)
181 preBasis.indices(tree, it);
184 template <
class PB,
class NIS,
class Tree,
class Iter,
185 class = decltype(std::declval<NIS>().
bind(std::declval<Tree>())),
186 class = decltype(std::declval<NIS>().indices(std::declval<Iter>()))>
187 static void bindNodeIndexSet(PB
const& , NIS& nodeIndexSet, Tree
const& tree, Iter it, Dune::PriorityTag<2>)
189 nodeIndexSet.bind(tree);
190 nodeIndexSet.indices(it);
193 template <
class PB,
class NIS,
class Tree,
class Iter,
194 class = decltype(std::declval<NIS>().
bind(std::declval<Tree>())),
195 class = decltype(std::declval<NIS>().
index(0u))>
196 static void bindNodeIndexSet(PB
const& , NIS& nodeIndexSet, Tree
const& tree, Iter it, Dune::PriorityTag<1>)
198 nodeIndexSet.bind(tree);
199 for (
size_type i = 0; i < tree.size(); ++i)
200 *it++ = nodeIndexSet.index(i);
203 template <
class PB,
class NIS,
class Tree,
class Iter>
204 static void bindNodeIndexSet(PB
const& , NIS& , Tree
const& , Iter , Dune::PriorityTag<0>)
210 class = decltype(std::declval<NIS>().
unbind())>
211 static void unbindNodeIndexSet(NIS& nodeIndexSet, Dune::PriorityTag<1>)
213 nodeIndexSet.unbind();
217 static void unbindNodeIndexSet(NIS& , Dune::PriorityTag<0>)
224 Element const* element_ =
nullptr;
227 NodeIndexSet nodeIndexSet_;
228 std::vector<MultiIndex> indices_;
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:141
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:151
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:32
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:171
typename GlobalBasis::PreBasis::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:47
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:106
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:158
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:35
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:20
Node_t< typename GlobalBasis::PreBasis, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: LocalView.hpp:41
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:135
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:129
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:165
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:112
LocalView(LocalView const &other)
Copy constructor.
Definition: LocalView.hpp:64
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:91
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:63
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:122
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:51
GB GlobalBasis
The global FE basis that this is a view on.
Definition: LocalView.hpp:29
NodeCache_t< Tree > TreeCache
Cached basis-tree.
Definition: LocalView.hpp:44
std::size_t size_type
The type used for sizes.
Definition: LocalView.hpp:38