Dune Core Modules (unstable)

refinedp0localbasis.hh
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 © 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 
13 #include <dune/localfunctions/common/localbasis.hh>
15 
16 namespace Dune
17 {
18 
38  template<class D, class R, int dim>
40  : public RefinedSimplexLocalBasis<D,dim>
41  {
42  // 2 to the k-th power
43  constexpr static int N = 1<<dim;
44  public:
47 
49  unsigned int size () const
50  {
51  return N;
52  }
53 
55  inline void evaluateFunction (const typename Traits::DomainType& in,
56  std::vector<typename Traits::RangeType>& out) const
57  {
58  int subElement = this->getSubElement(in);
59  out.resize(N);
60  for(int i=0; i<N; ++i)
61  out[i] = (i==subElement) ? 1 : 0;
62  }
63 
64  inline void
65  evaluateJacobian (const typename Traits::DomainType& in, // position
66  std::vector<typename Traits::JacobianType>& out) const // return value
67  {
68  out.resize(N);
69  for(int i=0; i<N; ++i)
70  out[i][0] = 0;
71  }
72 
74  void partial (const std::array<unsigned int, dim>& order,
75  const typename Traits::DomainType& in, // position
76  std::vector<typename Traits::RangeType>& out) const // return value
77  {
78  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
79  if (totalOrder == 0) {
80  evaluateFunction(in, out);
81  } else {
82  out.resize(size());
83  for (std::size_t i = 0; i < size(); ++i)
84  out[i] = 0;
85  }
86  }
87 
92  unsigned int order () const
93  {
94  return 0;
95  }
96 
97  };
98 
99 }
100 #endif
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Uniformly refined constant shape functions on a unit simplex in R^dim.
Definition: refinedp0localbasis.hh:41
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:74
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: refinedp0localbasis.hh:55
unsigned int order() const
Polynomial order of the shape functions.
Definition: refinedp0localbasis.hh:92
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:46
unsigned int size() const
number of shape functions
Definition: refinedp0localbasis.hh:49
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:279
Dune namespace.
Definition: alignedallocator.hh:13
Contains a base class for LocalBasis classes based on uniform refinement.
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:35
D DomainType
domain type
Definition: localbasis.hh:43
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)