dune-localfunctions 2.9.0
refinedp0localbasis.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_REFINED_P0_LOCALBASIS_HH
6#define DUNE_REFINED_P0_LOCALBASIS_HH
7
8#include <numeric>
9
10#include <dune/common/fvector.hh>
11#include <dune/common/fmatrix.hh>
12
15
16namespace Dune
17{
18
37 template<class D, class R, int dim>
39 : public RefinedSimplexLocalBasis<D,dim>
40 {
41 // 2 to the k-th power
42 constexpr static int N = 1<<dim;
43 public:
45 typedef LocalBasisTraits<D,dim,Dune::FieldVector<D,dim>,R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,dim> > Traits;
46
48 unsigned int size () const
49 {
50 return N;
51 }
52
54 inline void evaluateFunction (const typename Traits::DomainType& in,
55 std::vector<typename Traits::RangeType>& out) const
56 {
57 int subElement = this->getSubElement(in);
58 out.resize(N);
59 for(int i=0; i<N; ++i)
60 out[i] = (i==subElement) ? 1 : 0;
61 }
62
63 inline void
64 evaluateJacobian (const typename Traits::DomainType& in, // position
65 std::vector<typename Traits::JacobianType>& out) const // return value
66 {
67 out.resize(N);
68 for(int i=0; i<N; ++i)
69 out[i][0] = 0;
70 }
71
73 void partial (const std::array<unsigned int, dim>& order,
74 const typename Traits::DomainType& in, // position
75 std::vector<typename Traits::RangeType>& out) const // return value
76 {
77 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
78 if (totalOrder == 0) {
79 evaluateFunction(in, out);
80 } else {
81 out.resize(size());
82 for (std::size_t i = 0; i < size(); ++i)
83 out[i] = 0;
84 }
85 }
86
91 unsigned int order () const
92 {
93 return 0;
94 }
95
96 };
97
98}
99#endif
Contains a base class for LocalBasis classes based on uniform refinement.
Definition: bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:34
D DomainType
domain type
Definition: common/localbasis.hh:42
Definition: refinedsimplexlocalbasis.hh:20
Uniformly refined constant shape functions on a unit simplex in R^dim.
Definition: refinedp0localbasis.hh:40
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: refinedp0localbasis.hh:73
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: refinedp0localbasis.hh:54
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Definition: refinedp0localbasis.hh:64
unsigned int order() const
Polynomial order of the shape functions.
Definition: refinedp0localbasis.hh:91
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
export type traits for function signature
Definition: refinedp0localbasis.hh:45
unsigned int size() const
number of shape functions
Definition: refinedp0localbasis.hh:48