Dune Core Modules (2.10.0)

dualp1localbasis.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_DUAL_P1_LOCALBASIS_HH
6#define DUNE_DUAL_P1_LOCALBASIS_HH
7
8#include <numeric>
9
12#include <dune/localfunctions/common/localbasis.hh>
13
14namespace Dune
15{
33 template<class D, class R, int dim, bool faceDualT=false>
35 {
36 public:
38 static const bool faceDual = faceDualT;
42
44 unsigned int size () const
45 {
46 return dim+1;
47 }
48
50 inline void evaluateFunction (const typename Traits::DomainType& in,
51 std::vector<typename Traits::RangeType>& out) const
52 {
53 // evaluate P1 basis functions
54 std::vector<typename Traits::RangeType> p1Values(size());
55
56 p1Values[0] = 1.0;
57
58 for (int i=0; i<dim; i++) {
59 p1Values[0] -= in[i];
60 p1Values[i+1] = in[i];
61 }
62
63 // compute dual basis function values as a linear combination of the Lagrange values
64 out.resize(size());
65
66 for (int i=0; i<=dim; i++) {
67 out[i] = (dim+!faceDual)*p1Values[i];
68 for (int j=0; j<i; j++)
69 out[i] -= p1Values[j];
70
71 for (int j=i+1; j<=dim; j++)
72 out[i] -= p1Values[j];
73 }
74 }
75
77 inline void
79 std::vector<typename Traits::JacobianType>& out) const
80 {
81 // evaluate P1 jacobians
82 std::vector<typename Traits::JacobianType> p1Jacs(size());
83
84 for (int i=0; i<dim; i++)
85 p1Jacs[0][0][i] = -1;
86
87 for (int i=0; i<dim; i++)
88 for (int j=0; j<dim; j++)
89 p1Jacs[i+1][0][j] = (i==j);
90
91 // compute dual basis jacobians as linear combination of the Lagrange jacobians
92 out.resize(size());
93
94 for (size_t i=0; i<=dim; i++) {
95 out[i][0] = 0;
96 out[i][0].axpy(dim+!faceDual,p1Jacs[i][0]);
97
98 for (size_t j=0; j<i; j++)
99 out[i][0] -= p1Jacs[j][0];
100
101 for (int j=i+1; j<=dim; j++)
102 out[i][0] -= p1Jacs[j][0];
103 }
104 }
105
107 void partial (const std::array<unsigned int, dim>& order,
108 const typename Traits::DomainType& in, // position
109 std::vector<typename Traits::RangeType>& out) const // return value
110 {
111 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
112 if (totalOrder == 0) {
113 evaluateFunction(in, out);
114 } else {
115 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
116 }
117 }
118
120 unsigned int order () const
121 {
122 return 1;
123 }
124 };
125}
126#endif
Dual Lagrange shape functions on the simplex.
Definition: dualp1localbasis.hh:35
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualp1localbasis.hh:120
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualp1localbasis.hh:50
static const bool faceDual
Determines if the basis is only biorthogonal on adjacent faces.
Definition: dualp1localbasis.hh:38
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: dualp1localbasis.hh:107
unsigned int size() const
number of shape functions
Definition: dualp1localbasis.hh:44
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: dualp1localbasis.hh:41
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualp1localbasis.hh:78
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:91
Default exception for dummy implementations.
Definition: exceptions.hh:263
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.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:279
Dune namespace.
Definition: alignedallocator.hh:13
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.111.3 (Dec 27, 23:30, 2024)