Dune Core Modules (2.9.1)

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 (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_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{
32 template<class D, class R, int dim, bool faceDualT=false>
34 {
35 public:
37 static const bool faceDual = faceDualT;
41
43 unsigned int size () const
44 {
45 return dim+1;
46 }
47
49 inline void evaluateFunction (const typename Traits::DomainType& in,
50 std::vector<typename Traits::RangeType>& out) const
51 {
52 // evaluate P1 basis functions
53 std::vector<typename Traits::RangeType> p1Values(size());
54
55 p1Values[0] = 1.0;
56
57 for (int i=0; i<dim; i++) {
58 p1Values[0] -= in[i];
59 p1Values[i+1] = in[i];
60 }
61
62 // compute dual basis function values as a linear combination of the Lagrange values
63 out.resize(size());
64
65 for (int i=0; i<=dim; i++) {
66 out[i] = (dim+!faceDual)*p1Values[i];
67 for (int j=0; j<i; j++)
68 out[i] -= p1Values[j];
69
70 for (int j=i+1; j<=dim; j++)
71 out[i] -= p1Values[j];
72 }
73 }
74
76 inline void
78 std::vector<typename Traits::JacobianType>& out) const
79 {
80 // evaluate P1 jacobians
81 std::vector<typename Traits::JacobianType> p1Jacs(size());
82
83 for (int i=0; i<dim; i++)
84 p1Jacs[0][0][i] = -1;
85
86 for (int i=0; i<dim; i++)
87 for (int j=0; j<dim; j++)
88 p1Jacs[i+1][0][j] = (i==j);
89
90 // compute dual basis jacobians as linear combination of the Lagrange jacobians
91 out.resize(size());
92
93 for (size_t i=0; i<=dim; i++) {
94 out[i][0] = 0;
95 out[i][0].axpy(dim+!faceDual,p1Jacs[i][0]);
96
97 for (size_t j=0; j<i; j++)
98 out[i][0] -= p1Jacs[j][0];
99
100 for (int j=i+1; j<=dim; j++)
101 out[i][0] -= p1Jacs[j][0];
102 }
103 }
104
106 void partial (const std::array<unsigned int, dim>& order,
107 const typename Traits::DomainType& in, // position
108 std::vector<typename Traits::RangeType>& out) const // return value
109 {
110 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
111 if (totalOrder == 0) {
112 evaluateFunction(in, out);
113 } else {
114 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
115 }
116 }
117
119 unsigned int order () const
120 {
121 return 1;
122 }
123 };
124}
125#endif
Dual Lagrange shape functions on the simplex.
Definition: dualp1localbasis.hh:34
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualp1localbasis.hh:119
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualp1localbasis.hh:49
static const bool faceDual
Determines if the basis is only biorthogonal on adjacent faces.
Definition: dualp1localbasis.hh:37
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:106
unsigned int size() const
number of shape functions
Definition: dualp1localbasis.hh:43
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:40
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualp1localbasis.hh:77
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
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:291
Dune namespace.
Definition: alignedallocator.hh:13
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:42
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)