AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
LinearForm.inc.hpp
1#pragma once
2
3#include <utility>
4
5#include <amdis/GridFunctionOperator.hpp>
6#include <amdis/LocalOperator.hpp>
7#include <amdis/OperatorAdapter.hpp>
8#include <amdis/functions/EntitySet.hpp>
9#include <amdis/operations/Assigner.hpp>
10#include <amdis/typetree/Traversal.hpp>
11#include <amdis/typetree/TreePath.hpp>
12
13namespace AMDiS {
14
15template <class GB, class T, class Traits>
16 template <class ContextTag, class Expr, class TreePath>
18addOperator(ContextTag contextTag, Expr const& expr, TreePath path)
19{
20 static_assert( Concepts::PreTreePath<TreePath>,
21 "path must be a valid treepath, or an integer/index-constant");
22
23 using LocalContext = typename Impl::ContextTagType<ContextTag>::type;
24 auto op = makeOperator<LocalContext>(expr, this->basis().gridView());
25 if constexpr(std::is_same_v<TreePath,RootTreePath>)
26 assembler_.push(contextTag, std::move(op));
27 else
28 assembler_.push(contextTag, OperatorAdapter{std::move(op), makeTreePath(path)});
29}
30
31
32template <class GB, class T, class Traits>
33 template <class LocalView,
34 REQUIRES_(Concepts::LocalView<LocalView>)>
36assemble(LocalView const& localView)
37{
38 assert(localView.bound());
39
40 elementVector_.resize(localView.size());
41 elementVector_ = 0;
42
43 auto const& gv = this->basis().gridView();
44 auto const& element = localView.element();
45 GlobalContext<TYPEOF(gv)> context{gv, element, element.geometry()};
46
47 auto rhsOp = localAssembler(assembler_);
48 rhsOp.bind(element);
49 rhsOp.assemble(context, localView.treeCache(), elementVector_);
50 rhsOp.unbind();
51
52 this->scatter(localView, elementVector_, Assigner::plus_assign{});
53}
54
55
56template <class GB, class T, class Traits>
59{
60 auto localView = this->basis().localView();
61
62 this->init(this->basis(), true);
63 for (auto const& element : entitySet(this->basis())) {
64 localView.bind(element);
65 this->assemble(localView);
66 localView.unbind();
67 }
68 this->finish();
69}
70
71} // end namespace AMDiS
Definition: ContextGeometry.hpp:162
void assemble()
Assemble all vector operators added by addOperator().
Definition: LinearForm.inc.hpp:58
void addOperator(ContextTag contextTag, Operator const &op, TreePath path)
Associate a local operator with this LinearForm.
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:16
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:36
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:181
Definition: Assigner.hpp:17