AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
LinearForm.hpp
1#pragma once
2
3#include <dune/common/typeutilities.hh>
4
5#include <amdis/Assembler.hpp>
6#include <amdis/LinearAlgebra.hpp>
7#include <amdis/Observer.hpp>
8#include <amdis/common/Concepts.hpp>
9#include <amdis/common/ConceptsBase.hpp>
10#include <amdis/common/FlatVector.hpp>
11#include <amdis/common/TypeTraits.hpp>
12#include <amdis/typetree/TreePath.hpp>
13
14namespace AMDiS
15{
17
25 template <class GB, class T = double, class Traits = BackendTraits>
27 : public VectorFacade<T, Traits::template Vector<GB>::template Impl>
28 , private Observer<event::adapt>
29 {
31 using Self = LinearForm;
32
33 public:
35 using GlobalBasis = GB;
36
38 using CoefficientType = T;
39
42
43 public:
45 explicit LinearForm(std::shared_ptr<GB> const& basis)
46 : Super(*basis)
47 , Observer<event::adapt>(*basis)
48 , basis_(basis)
49 {
50 auto const localSize = basis->localView().maxSize();
51 elementVector_.resize(localSize);
52 }
53
55 template <class GB_,
56 REQUIRES(Concepts::Similar<GB_,GB>)>
57 explicit LinearForm(GB_&& basis)
58 : LinearForm(Dune::wrap_or_move(FWD(basis)))
59 {}
60
61 GlobalBasis const& basis() const { return *basis_; }
62
64
81 template <class ContextTag, class Operator, class TreePath>
82 void addOperator(ContextTag contextTag, Operator const& op, TreePath path);
83
84 // Add an operator to be assembled on the elements of the grid
85 template <class Operator, class TreePath = RootTreePath>
86 void addOperator(Operator const& op, TreePath path = {})
87 {
88 using E = typename GlobalBasis::LocalView::Element;
90 }
91
92 // Add an operator to be assembled on the intersections of the grid
93 template <class Operator, class TreePath = RootTreePath>
94 void addIntersectionOperator(Operator const& op, TreePath path = {})
95 {
96 using I = typename GlobalBasis::LocalView::GridView::Intersection;
97 addOperator(tag::intersection_operator<I>{}, op, path);
98 }
100
101
103 template <class LocalView,
104 REQUIRES(Concepts::LocalView<LocalView>)>
105 void assemble(LocalView const& localView);
106
108 void assemble();
109
110 auto& assembler() { return assembler_; }
111
113 void updateImpl(event::adapt e) override
114 {
115 if (e.value)
116 assembler_.update(basis_->gridView());
117 }
118
119 private:
121 ElementVector elementVector_;
122
125
126 std::shared_ptr<GlobalBasis const> basis_;
127 };
128
129
130 // deduction guide
131 template <class GB>
132 LinearForm(GB&& basis)
133 -> LinearForm<Underlying_t<GB>>;
134
135
136 template <class T = double, class GB>
137 auto makeLinearForm(GB&& basis)
138 {
139 using LF = LinearForm<Underlying_t<GB>, T>;
140 return LF{FWD(basis)};
141 }
142
143} // end namespace AMDiS
144
145#include <amdis/LinearForm.inc.hpp>
An Assembler is the main driver for building the local element matrices and vectors by assembling ope...
Definition: Assembler.hpp:132
An adaptive container that stores a value per grid element.
Definition: ElementVector.hpp:25
Global basis defined on a pre-basis.
Definition: GlobalBasis.hpp:48
The basic container that stores a base vector and a corresponding basis.
Definition: LinearForm.hpp:29
T CoefficientType
The type of the elements of the DOFVector.
Definition: LinearForm.hpp:38
GB GlobalBasis
The type of the functionspace basis associated to this linearform.
Definition: LinearForm.hpp:35
LinearForm(GB_ &&basis)
Wraps the passed global basis into a (non-destroying) shared_ptr.
Definition: LinearForm.hpp:57
void assemble()
Assemble all vector operators added by addOperator().
Definition: LinearForm.inc.hpp:58
void updateImpl(event::adapt e) override
Update all operators on the updated gridView.
Definition: LinearForm.hpp:113
LinearForm(std::shared_ptr< GB > const &basis)
Constructor. Stores the shared_ptr of the basis and creates a new DataTransfer.
Definition: LinearForm.hpp:45
void addOperator(ContextTag contextTag, Operator const &op, TreePath path)
Associate a local operator with this LinearForm.
Implementation of the ObserverInterface.
Definition: Observer.hpp:104
The base class for an operator to be used in an Assembler.
Definition: Operator.hpp:80
The basic container that stores a base vector and a corresponding basis.
Definition: VectorFacade.hpp:39
std::size_t localSize() const
Return the number of entries in the local part of the vector.
Definition: VectorFacade.hpp:84
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:181
Definition: Observer.hpp:25
Definition: Assembler.hpp:15