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