3#include <dune/common/typeutilities.hh>
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>
25 template <
class GB,
class T =
double,
class Traits = BackendTraits>
27 :
public VectorFacade<T, Traits::template Vector<GB>::template Impl>
50 auto const localSize = basis->localView().maxSize();
56 REQUIRES(Concepts::Similar<GB_,GB>)>
61 GlobalBasis const& basis()
const {
return *basis_; }
81 template <
class ContextTag,
class Operator,
class TreePath>
85 template <
class Operator,
class TreePath = RootTreePath>
88 using E =
typename GlobalBasis::LocalView::Element;
93 template <
class Operator,
class TreePath = RootTreePath>
94 void addIntersectionOperator(Operator
const& op, TreePath path = {})
96 using I =
typename GlobalBasis::LocalView::GridView::Intersection;
97 addOperator(tag::intersection_operator<I>{}, op, path);
104 REQUIRES(Concepts::LocalView<LocalView>)>
110 auto& assembler() {
return assembler_; }
116 assembler_.update(basis_->gridView());
126 std::shared_ptr<GlobalBasis const> basis_;
132 LinearForm(GB&& basis)
133 -> LinearForm<Underlying_t<GB>>;
136 template <
class T =
double,
class GB>
137 auto makeLinearForm(GB&& basis)
139 using LF = LinearForm<Underlying_t<GB>, T>;
140 return LF{FWD(basis)};
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
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