3#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
4#define DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
10#include <dune/typetree/treecontainer.hh>
12#include <dune/functions/functionspacebases/hierarchicnodetorangemap.hh>
13#include <dune/functions/functionspacebases/flatvectorview.hh>
14#include <dune/functions/gridfunctions/gridviewentityset.hh>
15#include <dune/functions/gridfunctions/gridfunction.hh>
16#include <dune/functions/backends/concepts.hh>
17#include <dune/functions/backends/istlvectorbackend.hh>
66template<
typename B,
typename V,
67 typename NTRE = HierarchicNodeToRangeMap,
68 typename R =
typename V::value_type>
79 using GridView =
typename Basis::GridView;
81 using Tree =
typename Basis::LocalView::Tree;
82 using NodeToRangeEntry = NTRE;
84 using Domain =
typename EntitySet::GlobalCoordinate;
94 using LocalView =
typename Basis::LocalView;
95 using size_type =
typename Tree::size_type;
98 using LocalBasisRange =
typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
101 using NodeData =
typename std::vector<LocalBasisRange<Node>>;
108 using Domain = LocalDomain;
109 using Range = GlobalFunction::Range;
110 using Element = GlobalFunction::Element;
113 : globalFunction_(&globalFunction)
114 , localView_(globalFunction.basis().localView())
115 , evaluationBuffer_(localView_.tree())
118 localDoFs_.reserve(localView_.maxSize());
121 LocalFunction(
const LocalFunction& other)
122 : globalFunction_(other.globalFunction_)
123 , localView_(other.localView_)
124 , evaluationBuffer_(localView_.tree())
125 , bound_(other.bound_)
127 localDoFs_.reserve(localView_.maxSize());
129 localDoFs_ = other.localDoFs_;
132 LocalFunction& operator=(
const LocalFunction& other)
134 globalFunction_ = other.globalFunction_;
135 localView_ = other.localView_;
136 bound_ = other.bound_;
138 localDoFs_ = other.localDoFs_;
148 void bind(
const Element& element)
150 localView_.bind(element);
163 localDoFs_.resize(localView_.size());
164 for (size_type i = 0; i < localView_.tree().size(); ++i)
169 size_t localIndex = localView_.tree().localIndex(i);
170 localDoFs_[localIndex] = globalFunction_->dofs()[localView_.index(localIndex)];
196 Range operator()(
const Domain& x)
const
201 const auto& nodeToRangeEntry = globalFunction_->nodeToRangeEntry();
202 const auto& fe = node.finiteElement();
203 const auto& localBasis = fe.localBasis();
204 auto& shapeFunctionValues = evaluationBuffer_[treePath];
206 localBasis.evaluateFunction(x, shapeFunctionValues);
209 auto re = flatVectorView(nodeToRangeEntry(node, treePath, y));
211 for (size_type i = 0; i < localBasis.size(); ++i)
214 auto c = flatVectorView(localDoFs_[node.localIndex(i)]);
217 auto v = flatVectorView(shapeFunctionValues[i]);
223 auto&& dimC = c.size();
224 auto dimV = v.size();
225 assert(dimC*dimV == re.size());
226 for(size_type j=0; j<dimC; ++j)
229 for(size_type k=0; k<dimV; ++k)
230 re[j*dimV + k] += c_j*v[k];
238 const Element& localContext()
const
240 return localView_.element();
243 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(
const LocalFunction& t)
251 LocalView localView_;
252 mutable PerNodeEvaluationBuffer evaluationBuffer_;
254 std::vector<Coefficient> localDoFs_;
257 template<
class B_T,
class V_T,
class NTRE_T>
259 entitySet_(basis.gridView()),
262 nodeToRangeEntry_(
wrap_or_move(
std::forward<NTRE_T>(nodeToRangeEntry)))
265 DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry) :
266 entitySet_(basis->gridView()),
268 coefficients_(coefficients),
269 nodeToRangeEntry_(nodeToRangeEntry)
272 const Basis& basis()
const
277 const V& dofs()
const
279 return *coefficients_;
282 const NodeToRangeEntry& nodeToRangeEntry()
const
284 return *nodeToRangeEntry_;
288 Range operator() (
const Domain& x)
const
293 friend typename Traits::DerivativeInterface
derivative(
const DiscreteGlobalBasisFunction& t)
315 return LocalFunction(t);
328 EntitySet entitySet_;
329 std::shared_ptr<const Basis> basis_;
330 std::shared_ptr<const V> coefficients_;
331 std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry_;
346template<
typename... TT>
351template<
typename R,
typename B,
typename V>
352auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
354 using Basis = std::decay_t<B>;
359 auto toConstVectorBackend = [&](
auto&& v) ->
decltype(
auto) {
360 if constexpr (models<Concept::ConstVectorBackend<Basis>,
decltype(v)>()) {
361 return std::forward<decltype(v)>(v);
367 using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
369 std::forward<B>(basis),
370 toConstVectorBackend(std::forward<V>(vector)),
A grid function induced by a global basis and a coefficient vector.
Definition: discreteglobalbasisfunction.hh:70
const EntitySet & entitySet() const
Get associated EntitySet.
Definition: discreteglobalbasisfunction.hh:321
GridView::template Codim< codim >::Entity Element
Type of Elements contained in this EntitySet.
Definition: gridviewentityset.hh:32
Element::Geometry::LocalCoordinate LocalCoordinate
Type of local coordinates with respect to the Element.
Definition: gridviewentityset.hh:35
Default exception for dummy implementations.
Definition: exceptions.hh:261
Traits for type conversions and type information.
typename AutonomousValueType< T >::type AutonomousValue
Type free of internal references that T can be converted to.
Definition: typetraits.hh:558
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:313
void localFunction(DiscreteGlobalBasisFunction< TT... > &&t)=delete
Construction of local functions from a temporary DiscreteGlobalBasisFunction (forbidden)
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:39
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:345
constexpr HybridTreePath< T... > treePath(const T &... t)
Constructs a new HybridTreePath from the given indices.
Definition: treepath.hh:188
void forEachLeafNode(Tree &&tree, LeafFunc &&leafFunc)
Traverse tree and visit each leaf node.
Definition: traversal.hh:304
std::decay_t< decltype(makeTreeContainer(std::declval< const Tree & >(), std::declval< Detail::LeafToDefaultConstructibleValue< LeafToValue > >()))> TreeContainer
Alias to container type generated by makeTreeContainer for give tree type and when using LeafToValue ...
Definition: treecontainer.hh:304
Dune namespace.
Definition: alignedallocator.hh:11
auto wrap_or_move(T &&t)
Capture R-value reference to shared_ptr.
Definition: shared_ptr.hh:94
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:37
A simple node to range map using the nested tree indices.
Definition: hierarchicnodetorangemap.hh:30