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>;
50 template <
class NIS,
class Iter>
51 using hasIndices = decltype(std::declval<NIS>().indices(std::declval<Iter>()));
56 : globalBasis_(&globalBasis)
57 , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
59 , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
61 static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>,
Tree>());
62 initializeTree(tree_);
69 : globalBasis_(other.globalBasis_)
72 , nodeIndexSet_(other.nodeIndexSet_)
73 , element_(other.element_)
74 , indices_(other.indices_)
81 : globalBasis_(std::move(other.globalBasis_))
82 , tree_(std::move(other.tree_))
84 , nodeIndexSet_(std::move(other.nodeIndexSet_))
85 , element_(std::move(other.element_))
86 , indices_(std::move(other.indices_))
98 bindTree(tree_, *element_);
99 nodeIndexSet_.bind(tree_);
100 indices_.resize(
size());
102 if constexpr (Dune::Std::is_detected_v<hasIndices, NodeIndexSet, TYPEOF(indices_.begin())>)
103 nodeIndexSet_.indices(indices_.begin());
106 indices_[i] = nodeIndexSet_.index(i);
112 return bool(element_);
127 nodeIndexSet_.unbind();
155 return globalBasis_->preBasis().maxNodeSize();
168 return *globalBasis_;
181 NodeIndexSet nodeIndexSet_;
183 std::optional<Element> element_;
184 std::vector<MultiIndex> indices_;
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:144
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:153
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:32
auto makeNodeCache(Node const &node)
Construct a new local-basis cache from a basis-node.
Definition: NodeCache.hpp:35
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:172
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:110
Contains all classes needed for solving linear and non linear equation systems.
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:160
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
typename Impl::NodeCacheFactory< Node >::type NodeCache_t
Defines the type of a node cache associated to a given Node.
Definition: NodeCache.hpp:31
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:138
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:132
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:166
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:116
LocalView(LocalView const &other)
Copy constructor.
Definition: LocalView.hpp:68
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:95
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:66
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:125
typename NodeIndexSet::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:47
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:55
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