BiLinearForm.inc.hpp 2.93 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
#pragma once

#include <utility>

#include <amdis/Assembler.hpp>
#include <amdis/LocalOperator.hpp>
#include <amdis/typetree/Traversal.hpp>
#include <amdis/utility/AssembleOperators.hpp>

namespace AMDiS {

Praetorius, Simon's avatar
Praetorius, Simon committed
12
template <class RB, class CB, class T, class Traits>
13
  template <class ContextTag, class Expr, class RowTreePath, class ColTreePath>
Praetorius, Simon's avatar
Praetorius, Simon committed
14
void BiLinearForm<RB,CB,T,Traits>::
15
16
17
18
19
20
21
22
addOperator(ContextTag contextTag, Expr const& expr,
            RowTreePath row, ColTreePath col)
{
  static_assert( Concepts::PreTreePath<RowTreePath>,
      "row must be a valid treepath, or an integer/index-constant");
  static_assert( Concepts::PreTreePath<ColTreePath>,
      "col must be a valid treepath, or an integer/index-constant");

23
  auto i = makeTreePath(row);
24
  auto node_i = child(this->rowBasis().localView().treeCache(), i);
25
  auto j = makeTreePath(col);
26
  auto node_j = child(this->colBasis().localView().treeCache(), j);
27
28

  using LocalContext = typename ContextTag::type;
Praetorius, Simon's avatar
Praetorius, Simon committed
29
  using Tr = DefaultAssemblerTraits<LocalContext, ElementMatrix>;
30
  auto op = makeLocalOperator<LocalContext>(expr, this->rowBasis().gridView());
31
  auto localAssembler = makeUniquePtr(makeAssembler<Tr>(std::move(op), node_i, node_j));
32
33

  operators_[i][j].push(contextTag, std::move(localAssembler));
34
  updatePattern_ = true;
35
36
37
}


Praetorius, Simon's avatar
Praetorius, Simon committed
38
template <class RB, class CB, class T, class Traits>
39
40
41
  template <class RowLocalView, class ColLocalView,
    REQUIRES_(Concepts::LocalView<RowLocalView>),
    REQUIRES_(Concepts::LocalView<ColLocalView>)>
Praetorius, Simon's avatar
Praetorius, Simon committed
42
void BiLinearForm<RB,CB,T,Traits>::
43
44
45
46
47
assemble(RowLocalView const& rowLocalView, ColLocalView const& colLocalView)
{
  elementMatrix_.resize(rowLocalView.size(), colLocalView.size());
  elementMatrix_ = 0;

48
  auto const& gv = this->rowBasis().gridView();
49
50
51
  auto const& element = rowLocalView.element();
  auto geometry = element.geometry();

52
53
  Traversal::forEachNode(rowLocalView.treeCache(), [&](auto const& rowNode, auto rowTp) {
    Traversal::forEachNode(colLocalView.treeCache(), [&](auto const& colNode, auto colTp) {
54
      auto& matOp = operators_[rowTp][colTp];
55
56
57
58
59
60
61
62
63
64
65
66
      if (matOp) {
        matOp.bind(element, geometry);
        assembleOperators(gv, element, matOp, makeMatrixAssembler(rowNode, colNode, elementMatrix_));
        matOp.unbind();
      }
    });
  });

  this->scatter(rowLocalView, colLocalView, elementMatrix_);
}


Praetorius, Simon's avatar
Praetorius, Simon committed
67
68
template <class RB, class CB, class T, class Traits>
void BiLinearForm<RB,CB,T,Traits>::
69
assemble()
70
{
71
72
  auto rowLocalView = this->rowBasis().localView();
  auto colLocalView = this->colBasis().localView();
73

74
  this->init();
75
  for (auto const& element : elements(this->rowBasis().gridView(), typename Traits::PartitionSet{})) {
76
    rowLocalView.bind(element);
Praetorius, Simon's avatar
Praetorius, Simon committed
77
    if (this->rowBasis_ == this->colBasis_)
78
79
80
81
82
83
      this->assemble(rowLocalView, rowLocalView);
    else {
      colLocalView.bind(element);
      this->assemble(rowLocalView, colLocalView);
      colLocalView.unbind();
    }
84
    rowLocalView.unbind();
85
86
87
88
89
90
  }
  this->finish();
}


} // end namespace AMDiS