3 #include <dune/common/typeutilities.hh> 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> 25 template <
class GB,
class T =
double,
class Traits = BackendTraits<GB>>
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>
93 template <
class Operator,
class TreePath = RootTreePath>
94 void addIntersectionOperator(
Operator const& op, TreePath path = {})
96 using I =
typename GlobalBasis::LocalView::GridView::Intersection;
103 template <
class LocalView,
class LocalOperators,
104 REQUIRES(Concepts::LocalView<LocalView>)>
110 auto& operators() {
return operators_; }
116 Recursive::forEach(operators_, [&](
auto& op) { op.
update(basis_->gridView()); });
125 VectorOperators<GlobalBasis,ElementVector> operators_;
127 std::shared_ptr<GlobalBasis const> basis_;
137 template <
class T =
double,
class GB>
138 auto makeLinearForm(GB&& basis)
141 return LF{FWD(basis)};
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
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
The base class for an operator to be used in an Assembler.
Definition: Operator.hpp:79
Definition: AdaptiveGrid.hpp:373
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
Definition: Observer.hpp:25
Definition: OperatorList.hpp:15
void update(GridView const &gv)
Update the operator data on a GridView.
Definition: Operator.hpp:103
Definition: OperatorList.hpp:16