DUNE PDELab (git)

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// 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_RT0TRIANGLELOCALBASIS_HH
6#define DUNE_RT0TRIANGLELOCALBASIS_HH
7
8#include <numeric>
9
11
12#include <dune/localfunctions/common/localbasis.hh>
13
14namespace Dune
15{
25 template<class D, class R>
27 {
28 public:
31
33 RT02DLocalBasis (std::bitset<3> s = 0)
34 {
35 for (int i=0; i<3; i++)
36 sign_[i] = s[i] ? -1.0 : 1.0;
37 }
38
40 unsigned int size () const
41 {
42 return 3;
43 }
44
46 inline void evaluateFunction (const typename Traits::DomainType& in,
47 std::vector<typename Traits::RangeType>& out) const
48 {
49 out.resize(3);
50 out[0] = {sign_[0]*in[0], sign_[0]*(in[1]-D(1))};
51 out[1] = {sign_[1]*(in[0]-D(1)), sign_[1]*in[1]};
52 out[2] = {sign_[2]*in[0], sign_[2]*in[1]};
53 }
54
56 inline void
57 evaluateJacobian (const typename Traits::DomainType& in, // position
58 std::vector<typename Traits::JacobianType>& out) const // return value
59 {
60 out.resize(3);
61 for (int i=0; i<3; i++)
62 {
63 out[i][0] = {sign_[i], 0};
64 out[i][1] = { 0, sign_[i]};
65 }
66 }
67
69 void partial (const std::array<unsigned int, 2>& order,
70 const typename Traits::DomainType& in, // position
71 std::vector<typename Traits::RangeType>& out) const // return value
72 {
73 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
74 if (totalOrder == 0) {
75 evaluateFunction(in, out);
76 } else if (totalOrder == 1) {
77 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
78 out.resize(size());
79
80 for (int i=0; i<3; i++)
81 {
82 out[i][direction] = sign_[i];
83 out[i][1-direction] = 0;
84 }
85 } else {
86 out.resize(size());
87 for (std::size_t i = 0; i < size(); ++i)
88 for (std::size_t j = 0; j < 2; ++j)
89 out[i][j] = 0;
90 }
91
92 }
93
95 unsigned int order () const
96 {
97 return 1;
98 }
99
100 private:
101
102 // Signs of the edge normals
103 std::array<R,3> sign_;
104 };
105}
106#endif
A dense n x m matrix.
Definition: fmatrix.hh:117
Definition: raviartthomas02dlocalbasis.hh:27
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas02dlocalbasis.hh:95
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas02dlocalbasis.hh:57
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas02dlocalbasis.hh:46
unsigned int size() const
number of shape functions
Definition: raviartthomas02dlocalbasis.hh:40
RT02DLocalBasis(std::bitset< 3 > s=0)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas02dlocalbasis.hh:33
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:69
Implements a matrix constructed from a given type representing a field and compile-time given number ...
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 (Nov 13, 23:29, 2024)