3 #include <dune/common/typeutilities.hh> 5 #include <amdis/LinearAlgebra.hpp> 6 #include <amdis/OperatorList.hpp> 7 #include <amdis/common/Concepts.hpp> 8 #include <amdis/common/ConceptsBase.hpp> 9 #include <amdis/common/FlatVector.hpp> 10 #include <amdis/common/TypeTraits.hpp> 11 #include <amdis/typetree/TreePath.hpp> 24 template <
class GB,
class T =
double,
class Traits = BackendTraits<GB>>
47 auto const localSize = basis->localView().maxSize();
53 REQUIRES(Concepts::Similar<GB_,GB>)>
58 GlobalBasis const& basis()
const {
return *basis_; }
78 template <
class ContextTag,
class Operator,
class TreePath>
79 void addOperator(ContextTag contextTag, Operator
const& op, TreePath path);
82 template <
class Operator,
class TreePath = RootTreePath>
83 void addOperator(Operator
const& op, TreePath path = {})
90 template <
class Operator,
class TreePath = RootTreePath>
91 void addIntersectionOperator(Operator
const& op, TreePath path = {})
93 using I =
typename GlobalBasis::LocalView::GridView::Intersection;
101 REQUIRES(Concepts::LocalView<LocalView>)>
102 void assemble(LocalView
const& localView);
112 VectorOperators<GlobalBasis,ElementVector> operators_;
114 std::shared_ptr<GlobalBasis const> basis_;
124 template <
class T =
double,
class GB>
125 auto makeLinearForm(GB&& basis)
128 return LF{FWD(basis)};
133 #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:82
The basic container that stores a base vector and a corresponding basis.
Definition: LinearSolverInterface.hpp:23
Definition: AdaptiveGrid.hpp:373
Contains all classes needed for solving linear and non linear equation systems.
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
Definition: OperatorList.hpp:14
Definition: OperatorList.hpp:15