Dune Core Modules (2.7.0)

raviartthomas02dlocalbasis.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_RT0TRIANGLELOCALBASIS_HH
4#define DUNE_RT0TRIANGLELOCALBASIS_HH
5
6#include <numeric>
7
9
10#include <dune/localfunctions/common/localbasis.hh>
11
12namespace Dune
13{
22 template<class D, class R>
24 {
25 public:
28
31 {
32 std::fill(sign_.begin(), sign_.end(), 1.0);
33 }
34
36 RT02DLocalBasis (std::bitset<3> s)
37 {
38 for (int i=0; i<3; i++)
39 sign_[i] = s[i] ? -1.0 : 1.0;
40 }
41
43 unsigned int size () const
44 {
45 return 3;
46 }
47
49 inline void evaluateFunction (const typename Traits::DomainType& in,
50 std::vector<typename Traits::RangeType>& out) const
51 {
52 out.resize(3);
53 out[0] = {sign_[0]*in[0], sign_[0]*(in[1]-D(1))};
54 out[1] = {sign_[1]*(in[0]-D(1)), sign_[1]*in[1]};
55 out[2] = {sign_[2]*in[0], sign_[2]*in[1]};
56 }
57
59 inline void
60 evaluateJacobian (const typename Traits::DomainType& in, // position
61 std::vector<typename Traits::JacobianType>& out) const // return value
62 {
63 out.resize(3);
64 for (int i=0; i<3; i++)
65 {
66 out[i][0] = {sign_[i], 0};
67 out[i][1] = { 0, sign_[i]};
68 }
69 }
70
72 void partial (const std::array<unsigned int, 2>& order,
73 const typename Traits::DomainType& in, // position
74 std::vector<typename Traits::RangeType>& out) const // return value
75 {
76 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
77 if (totalOrder == 0) {
78 evaluateFunction(in, out);
79 } else if (totalOrder == 1) {
80 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
81 out.resize(size());
82
83 for (int i=0; i<3; i++)
84 {
85 out[i][direction] = sign_[i];
86 out[i][1-direction] = 0;
87 }
88 } else {
89 out.resize(size());
90 for (std::size_t i = 0; i < size(); ++i)
91 for (std::size_t j = 0; j < 2; ++j)
92 out[i][j] = 0;
93 }
94
95 }
96
98 unsigned int order () const
99 {
100 return 1;
101 }
102
103 private:
104
105 // Signs of the edge normals
106 std::array<R,3> sign_;
107 };
108}
109#endif
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Lowest order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas02dlocalbasis.hh:24
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas02dlocalbasis.hh:98
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:60
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas02dlocalbasis.hh:49
unsigned int size() const
number of shape functions
Definition: raviartthomas02dlocalbasis.hh:43
RT02DLocalBasis()
Standard constructor.
Definition: raviartthomas02dlocalbasis.hh:30
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:72
RT02DLocalBasis(std::bitset< 3 > s)
Make set numer s, where 0<=s<4.
Definition: raviartthomas02dlocalbasis.hh:36
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:290
Dune namespace.
Definition: alignedallocator.hh:14
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)