AMDiS  0.3
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/operations/Assigner.hpp>
8 #include <amdis/typetree/Traversal.hpp>
9 #include <amdis/typetree/TreePath.hpp>
10 
11 namespace AMDiS {
12 
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)
17 {
18  static_assert( Concepts::PreTreePath<TreePath>,
19  "path must be a valid treepath, or an integer/index-constant");
20 
21  auto i = makeTreePath(path);
22 
23  using LocalContext = typename ContextTag::type;
24  auto op = makeOperator<LocalContext>(expr, this->basis().gridView());
25  operators_[i].push(contextTag, std::move(op));
26 }
27 
28 
29 template <class GB, class T, class Traits>
30  template <class LocalView, class LocalOperators,
31  REQUIRES_(Concepts::LocalView<LocalView>)>
33 assemble(LocalView const& localView, LocalOperators& localOperators)
34 {
35  assert(localView.isBound());
36 
37  elementVector_.resize(localView.size());
38  elementVector_ = 0;
39 
40  auto const& gv = this->basis().gridView();
41  auto const& element = localView.element();
42  GlobalContext<TYPEOF(gv)> context{gv, element, element.geometry()};
43 
44  Traversal::forEachNode(localView.treeCache(), [&](auto const& node, auto tp) {
45  auto& rhsOp = localOperators[tp];
46  rhsOp.bind(element);
47  rhsOp.assemble(context, node, elementVector_);
48  rhsOp.unbind();
49  });
50 
51  this->scatter(localView, elementVector_, Assigner::plus_assign{});
52 }
53 
54 
55 template <class GB, class T, class Traits>
58 {
59  auto localView = this->basis().localView();
60  auto localOperators = AMDiS::localOperators(operators_);
61 
62  this->init(sizeInfo(this->basis()), true);
63  for (auto const& element : elements(this->basis().gridView(), typename Traits::PartitionSet{})) {
64  localView.bind(element);
65  this->assemble(localView, localOperators);
66  localView.unbind();
67  }
68  this->finish();
69 }
70 
71 } // end namespace AMDiS
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:141
void addOperator(ContextTag contextTag, Operator const &op, TreePath path)
Associate a local operator with this LinearForm.
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
void assemble()
Assemble all vector operators added by addOperator().
Definition: LinearForm.inc.hpp:57
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:112