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