AMDiS 2.10
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/OperatorAdapter.hpp>
9#include <amdis/functions/LocalViewPair.hpp>
10#include <amdis/typetree/Traversal.hpp>
11
12namespace AMDiS {
13
14template <class RB, class CB, class T, class Traits>
15 template <class ContextTag, class Expr, class RowTreePath, class ColTreePath>
17addOperator(ContextTag contextTag, Expr const& expr,
18 RowTreePath row, ColTreePath col)
19{
20 static_assert( Concepts::PreTreePath<RowTreePath>,
21 "row must be a valid treepath, or an integer/index-constant");
22 static_assert( Concepts::PreTreePath<ColTreePath>,
23 "col must be a valid treepath, or an integer/index-constant");
24
25 using LocalContext = typename Impl::ContextTagType<ContextTag>::type;
26 auto op = makeOperator<LocalContext>(expr, this->rowBasis().gridView());
27 if constexpr(std::is_same_v<RowTreePath,RootTreePath> && std::is_same_v<ColTreePath,RootTreePath>)
28 assembler_.push(contextTag, std::move(op));
29 else
30 assembler_.push(contextTag, OperatorAdapter{std::move(op),makeTreePath(row),makeTreePath(col)});
31 updatePattern_ = true;
32}
33
34
35template <class RB, class CB, class T, class Traits>
36 template <class RowLocalView, class ColLocalView,
37 REQUIRES_(Concepts::LocalView<RowLocalView>),
38 REQUIRES_(Concepts::LocalView<ColLocalView>)>
40assemble(RowLocalView const& rowLocalView, ColLocalView const& colLocalView)
41{
42 assert(rowLocalView.bound());
43 assert(colLocalView.bound());
44
45 elementMatrix_.resize(rowLocalView.size(), colLocalView.size());
46 elementMatrix_ = 0;
47
48 auto const& gv = this->rowBasis().gridView();
49 auto const& element = rowLocalView.element();
50 GlobalContext<TYPEOF(gv)> context{gv, element, element.geometry()};
51
52 auto matOp = localAssembler(assembler_);
53 matOp.bind(element);
54 matOp.assemble(context, rowLocalView.treeCache(), colLocalView.treeCache(), elementMatrix_);
55 matOp.unbind();
56
57 this->scatter(rowLocalView, colLocalView, elementMatrix_);
58}
59
60
61template <class RB, class CB, class T, class Traits>
64{
65 LocalViewPair<RB,CB> localViews{this->rowBasis(), this->colBasis()};
66
67 this->init();
68 for (auto const& element : entitySet(this->rowBasis())) {
69 localViews.bind(element);
70 this->assemble(localViews.row(), localViews.col());
71 localViews.unbind();
72 }
73 this->finish();
74}
75
76
77} // end namespace AMDiS
void addOperator(ContextTag contextTag, Operator const &op, Row row, Col col)
Associate a local operator with this BiLinearForm.
void assemble()
Assemble all matrix operators, TODO: incorporate boundary conditions.
Definition: BiLinearForm.inc.hpp:63
Definition: ContextGeometry.hpp:162
The restriction of a finite element basis to a single element.
Definition: LocalViewPair.hpp:10