Dune Core Modules (2.6.0)

p1localbasis.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_P1_LOCALBASIS_HH
4#define DUNE_P1_LOCALBASIS_HH
5
6#include <array>
7#include <numeric>
8
10
11#include <dune/localfunctions/common/localbasis.hh>
12
13namespace Dune
14{
26 template<class D, class R, int dim>
28 {
29 public:
33
35 unsigned int size () const
36 {
37 return dim+1;
38 }
39
41 inline void evaluateFunction (const typename Traits::DomainType& in,
42 std::vector<typename Traits::RangeType>& out) const
43 {
44 out.resize(size());
45 out[0] = 1.0;
46 for (size_t i=0; i<dim; i++) {
47 out[0] -= in[i];
48 out[i+1] = in[i];
49 }
50 }
51
53 inline void
54 evaluateJacobian (const typename Traits::DomainType& in, // position
55 std::vector<typename Traits::JacobianType>& out) const // return value
56 {
57 out.resize(size());
58
59 for (int i=0; i<dim; i++)
60 out[0][0][i] = -1;
61
62 for (int i=0; i<dim; i++)
63 for (int j=0; j<dim; j++)
64 out[i+1][0][j] = (i==j);
65
66 }
67
73 inline void partial(const std::array<unsigned int,dim>& order,
74 const typename Traits::DomainType& in,
75 std::vector<typename Traits::RangeType>& out) const
76 {
77 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
78
79 if (totalOrder==0)
80 evaluateFunction(in, out);
81 else if (totalOrder==1)
82 {
83 auto direction = std::find(order.begin(), order.end(), 1);
84 out.resize(size());
85
86 out[0] = -1;
87 for (int i=0; i<dim; i++)
88 out[i+1] = (i==(direction-order.begin()));
89 }
90 else // all higher order derivatives are zero
91 {
92 out.resize(size());
93
94 for (int i=0; i<dim+1; i++)
95 out[i] = 0;
96 }
97 }
98
100 unsigned int order () const
101 {
102 return 1;
103 }
104 };
105}
106#endif
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Linear Lagrange shape functions on the simplex.
Definition: p1localbasis.hh:28
unsigned int order() const
Polynomial order of the shape functions.
Definition: p1localbasis.hh:100
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.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: p1localbasis.hh:32
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: p1localbasis.hh:54
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 any order of all shape functions.
Definition: p1localbasis.hh:73
unsigned int size() const
number of shape functions
Definition: p1localbasis.hh:35
Implements a matrix constructed from a given type representing a field and compile-time given number ...
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 (Jul 15, 22:36, 2024)