BiLinearForm.inc.hpp 2.66 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
#pragma once

#include <utility>

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

namespace AMDiS {

template <class RB, class CB, class T>
  template <class ContextTag, class Expr, class RowTreePath, class ColTreePath>
void BiLinearForm<RB,CB,T>::
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");

  auto i = child(this->rowBasis()->localView().tree(), makeTreePath(row));
  auto j = child(this->colBasis()->localView().tree(), makeTreePath(col));

  using LocalContext = typename ContextTag::type;
  using Traits = DefaultAssemblerTraits<LocalContext, ElementMatrix>;
  auto op = makeLocalOperator<LocalContext>(expr, this->rowBasis()->gridView());
  auto localAssembler = makeUniquePtr(makeAssembler<Traits>(std::move(op), i, j));

  operators_[i][j].push(contextTag, std::move(localAssembler));
}


template <class RB, class CB, class T>
void BiLinearForm<RB,CB,T>::
assemble(RowLocalView const& rowLocalView, ColLocalView const& colLocalView)
{
  elementMatrix_.resize(rowLocalView.size(), colLocalView.size());
  elementMatrix_ = 0;

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

  for_each_node(rowLocalView.tree(), [&](auto const& rowNode, auto) {
    for_each_node(colLocalView.tree(), [&](auto const& colNode, auto) {
      auto& matOp = operators_[rowNode][colNode];
      if (matOp) {
        matOp.bind(element, geometry);
        assembleOperators(gv, element, matOp, makeMatrixAssembler(rowNode, colNode, elementMatrix_));
        matOp.unbind();
      }
    });
  });

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


template <class RB, class CB, class T>
void BiLinearForm<RB,CB,T>::
assemble(SymmetryStructure symmetry)
{
  auto rowLocalView = this->rowBasis()->localView();
  auto colLocalView = this->colBasis()->localView();

  this->init(symmetry);
  for (auto const& element : elements(this->rowBasis()->gridView(), typename Traits::PartitionSet{})) {
    rowLocalView.bind(element);
    if (this->rowBasis() == this->colBasis())
      this->assemble(rowLocalView, rowLocalView);
    else {
      colLocalView.bind(element);
      this->assemble(rowLocalView, colLocalView);
      colLocalView.unbind();
    }
    rowLocalView.unbind(element);
  }
  this->finish();
}


} // end namespace AMDiS