AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
LinearForm.hpp
1 #pragma once
2 
3 #include <dune/common/typeutilities.hh>
4 
5 #include <amdis/LinearAlgebra.hpp>
6 #include <amdis/Observer.hpp>
7 #include <amdis/OperatorList.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 
14 namespace AMDiS
15 {
17 
25  template <class GB, class T = double, class Traits = BackendTraits<GB>>
26  class LinearForm
27  : public VectorFacade<T, Traits::template VectorImpl>
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;
98  }
100 
101 
103  template <class LocalView, class LocalOperators,
104  REQUIRES(Concepts::LocalView<LocalView>)>
105  void assemble(LocalView const& localView, LocalOperators& localOperators);
106 
108  void assemble();
109 
110  auto& operators() { return operators_; }
111 
113  void updateImpl(event::adapt e) override
114  {
115  if (e.value) {
116  Recursive::forEach(operators_, [&](auto& op) { op.update(basis_->gridView()); });
117  }
118  }
119 
120  private:
122  ElementVector elementVector_;
123 
125  VectorOperators<GlobalBasis,ElementVector> operators_;
126 
127  std::shared_ptr<GlobalBasis const> basis_;
128  };
129 
130 
131  // deduction guide
132  template <class GB>
133  LinearForm(GB&& basis)
135 
136 
137  template <class T = double, class GB>
138  auto makeLinearForm(GB&& basis)
139  {
140  using LF = LinearForm<Underlying_t<GB>, T>;
141  return LF{FWD(basis)};
142  }
143 
144 } // end namespace AMDiS
145 
146 #include <amdis/LinearForm.inc.hpp>
std::size_t localSize() const
Return the number of entries in the local part of the vector.
Definition: VectorFacade.hpp:85
void addOperator(ContextTag contextTag, Operator const &op, TreePath path)
Associate a local operator with this LinearForm.
The basic container that stores a base vector and a corresponding basis.
Definition: VectorFacade.hpp:39
auto localOperators(Container const &c)
Definition: Operator.hpp:125
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
The base class for an operator to be used in an Assembler.
Definition: Operator.hpp:79
Definition: AdaptiveGrid.hpp:373
The basic container that stores a base vector and a corresponding basis.
Definition: LinearForm.hpp:26
Definition: AdaptBase.hpp:6
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:182
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:35
Implementation of the ObserverInterface.
Definition: Observer.hpp:100
void assemble()
Assemble all vector operators added by addOperator().
Definition: LinearForm.inc.hpp:57
Definition: Observer.hpp:25
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
Definition: OperatorList.hpp:15
void update(GridView const &gv)
Update the operator data on a GridView.
Definition: Operator.hpp:103
T CoefficientType
The type of the elements of the DOFVector.
Definition: LinearForm.hpp:38
Definition: OperatorList.hpp:16