AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
BiLinearForm.inc.hpp
1 #pragma once
2 
3 #include <utility>
4 
5 #include <amdis/ContextGeometry.hpp>
6 #include <amdis/GridFunctionOperator.hpp>
7 #include <amdis/LocalOperator.hpp>
8 #include <amdis/typetree/Traversal.hpp>
9 
10 namespace AMDiS {
11 
12 template <class RB, class CB, class T, class Traits>
13  template <class ContextTag, class Expr, class RowTreePath, class ColTreePath>
15 addOperator(ContextTag contextTag, Expr const& expr,
16  RowTreePath row, ColTreePath col)
17 {
18  static_assert( Concepts::PreTreePath<RowTreePath>,
19  "row must be a valid treepath, or an integer/index-constant");
20  static_assert( Concepts::PreTreePath<ColTreePath>,
21  "col must be a valid treepath, or an integer/index-constant");
22 
23  auto i = makeTreePath(row);
24  auto j = makeTreePath(col);
25 
26  using LocalContext = typename ContextTag::type;
27  auto op = makeOperator<LocalContext>(expr, this->rowBasis().gridView());
28  operators_[i][j].push(contextTag, std::move(op));
29  updatePattern_ = true;
30 }
31 
32 
33 template <class RB, class CB, class T, class Traits>
34  template <class RowLocalView, class ColLocalView, class LocalOperators,
35  REQUIRES_(Concepts::LocalView<RowLocalView>),
36  REQUIRES_(Concepts::LocalView<ColLocalView>)>
38 assemble(RowLocalView const& rowLocalView, ColLocalView const& colLocalView,
39  LocalOperators& localOperators)
40 {
41  assert(rowLocalView.isBound());
42  assert(colLocalView.isBound());
43 
44  elementMatrix_.resize(rowLocalView.size(), colLocalView.size());
45  elementMatrix_ = 0;
46 
47  auto const& gv = this->rowBasis().gridView();
48  auto const& element = rowLocalView.element();
49  GlobalContext<TYPEOF(gv)> context{gv, element, element.geometry()};
50 
51  Traversal::forEachNode(rowLocalView.treeCache(), [&](auto const& rowCache, auto rowTp) {
52  Traversal::forEachNode(colLocalView.treeCache(), [&](auto const& colCache, auto colTp) {
53  auto& matOp = localOperators[rowTp][colTp];
54  matOp.bind(element);
55  matOp.assemble(context, rowCache, colCache, elementMatrix_);
56  matOp.unbind();
57  });
58  });
59 
60  this->scatter(rowLocalView, colLocalView, elementMatrix_);
61 }
62 
63 
64 template <class RB, class CB, class T, class Traits>
67 {
68  auto rowLocalView = this->rowBasis().localView();
69  auto colLocalView = this->colBasis().localView();
70 
71  this->init();
72  auto localOperators = AMDiS::localOperators(operators_);
73  for (auto const& element : elements(this->rowBasis().gridView(), typename Traits::PartitionSet{})) {
74  rowLocalView.bind(element);
75  if (this->rowBasis_ == this->colBasis_)
76  this->assemble(rowLocalView, rowLocalView, localOperators);
77  else {
78  colLocalView.bind(element);
79  this->assemble(rowLocalView, colLocalView, localOperators);
80  colLocalView.unbind();
81  }
82  rowLocalView.unbind();
83  }
84  this->finish();
85 }
86 
87 
88 } // end namespace AMDiS
auto localOperators(Container const &c)
Definition: Operator.hpp:125
Definition: ContextGeometry.hpp:168
Definition: AdaptBase.hpp:6
Geometry const & geometry() const
Return the geometry of the Element.
Definition: ContextGeometry.hpp:202
void assemble()
Assemble all matrix operators, TODO: incorporate boundary conditions.
Definition: BiLinearForm.inc.hpp:66
void addOperator(ContextTag contextTag, Operator const &op, Row row, Col col)
Associate a local operator with this BiLinearForm.