5 #include <amdis/GridFunctionOperator.hpp> 6 #include <amdis/LocalOperator.hpp> 7 #include <amdis/operations/Assigner.hpp> 8 #include <amdis/typetree/Traversal.hpp> 9 #include <amdis/typetree/TreePath.hpp> 13 template <
class GB,
class T,
class Traits>
14 template <
class ContextTag,
class Expr,
class TreePath>
16 addOperator(ContextTag contextTag, Expr
const& expr, TreePath path)
18 static_assert( Concepts::PreTreePath<TreePath>,
19 "path must be a valid treepath, or an integer/index-constant");
21 auto i = makeTreePath(path);
23 using LocalContext =
typename ContextTag::type;
24 auto op = makeOperator<LocalContext>(expr, this->basis().gridView());
25 operators_[i].push(contextTag, std::move(op));
29 template <
class GB,
class T,
class Traits>
30 template <
class LocalView,
class LocalOperators,
31 REQUIRES_(Concepts::LocalView<LocalView>)>
37 elementVector_.resize(localView.
size());
40 auto const& gv = this->basis().gridView();
41 auto const& element = localView.
element();
44 Traversal::forEachNode(localView.
treeCache(), [&](
auto const& node,
auto tp) {
45 auto& rhsOp = localOperators[tp];
47 rhsOp.assemble(context, node, elementVector_);
55 template <
class GB,
class T,
class Traits>
59 auto localView = this->basis().localView();
62 this->init(sizeInfo(this->basis()),
true);
63 for (
auto const& element : elements(this->basis().gridView(),
typename Traits::PartitionSet{})) {
64 localView.bind(element);
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:141
auto localOperators(Container const &c)
Definition: Operator.hpp:125
Definition: Assigner.hpp:16
Definition: ContextGeometry.hpp:168
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:106
Definition: AdaptBase.hpp:6
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:182
Geometry const & geometry() const
Return the geometry of the Element.
Definition: ContextGeometry.hpp:202
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:20
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:135
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:112