AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
LinearForm.inc.hpp
1 #pragma once
2 
3 #include <utility>
4 
5 #include <amdis/Assembler.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 #include <amdis/utility/AssembleOperators.hpp>
11 
12 namespace AMDiS {
13 
14 template <class GB, class T, class Traits>
15  template <class ContextTag, class Expr, class TreePath>
17 addOperator(ContextTag contextTag, Expr const& expr, TreePath path)
18 {
19  static_assert( Concepts::PreTreePath<TreePath>,
20  "path must be a valid treepath, or an integer/index-constant");
21 
22  auto i = makeTreePath(path);
23  auto node = child(this->basis().localView().treeCache(), i);
24 
25  using LocalContext = typename ContextTag::type;
26  using Tr = DefaultAssemblerTraits<LocalContext, ElementVector>;
27  auto op = makeLocalOperator<LocalContext>(expr, this->basis().gridView());
28  auto localAssembler = makeUniquePtr(makeAssembler<Tr>(std::move(op), node));
29 
30  operators_[i].push(contextTag, std::move(localAssembler));
31 }
32 
33 
34 template <class GB, class T, class Traits>
35  template <class LocalView, REQUIRES_(Concepts::LocalView<LocalView>)>
37 assemble(LocalView const& localView)
38 {
39  elementVector_.resize(localView.size());
40  elementVector_ = 0;
41 
42  auto const& gv = this->basis().gridView();
43  auto const& element = localView.element();
44  auto geometry = element.geometry();
45 
46  Traversal::forEachNode(localView.treeCache(), [&](auto const& node, auto tp) {
47  auto& rhsOp = operators_[tp];
48  if (rhsOp) {
49  rhsOp.bind(element, geometry);
50  assembleOperators(gv, element, rhsOp, makeVectorAssembler(node, elementVector_));
51  rhsOp.unbind();
52  }
53  });
54 
55  this->scatter(localView, elementVector_, Assigner::plus_assign{});
56 }
57 
58 
59 template <class GB, class T, class Traits>
62 {
63  auto localView = this->basis().localView();
64 
65  this->init(sizeInfo(*this->basis()), true);
66  for (auto const& element : elements(this->basis().gridView(), typename Traits::PartitionSet{})) {
67  localView.bind(element);
68  this->assemble(localView);
69  localView.unbind();
70  }
71  this->finish();
72 }
73 
74 } // end namespace AMDiS
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:144
void addOperator(ContextTag contextTag, Operator const &op, TreePath path)
Associate a local operator with this LinearForm.
Definition: Assigner.hpp:16
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
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:138
void assemble()
Assemble all vector operators added by addOperator().
Definition: LinearForm.inc.hpp:61
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:116
auto makeUniquePtr(Obj &&obj)
Create a unique_ptr by copy/move construction.
Definition: TypeTraits.hpp:107
void init(int &argc, char **&argv, std::string const &initFileName="")
Initialized the Environment for MPI.
Definition: AMDiS.hpp:29