Dune Core Modules (2.9.1)

brezzidouglasmarini1simplex2dlocalbasis.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_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
7
8#include <array>
9#include <bitset>
10#include <numeric>
11#include <vector>
12
14
15#include "../../common/localbasis.hh"
16
17namespace Dune
18{
28 template<class D, class R>
30 {
31
32 public:
35
38 {
39 for (size_t i=0; i<3; i++)
40 sign_[i] = 1.0;
41 }
42
48 BDM1Simplex2DLocalBasis (std::bitset<3> s)
49 {
50 for (size_t i=0; i<3; i++)
51 sign_[i] = s[i] ? -1.0 : 1.0;
52 }
53
55 unsigned int size () const
56 {
57 return 6;
58 }
59
66 inline void evaluateFunction (const typename Traits::DomainType& in,
67 std::vector<typename Traits::RangeType>& out) const
68 {
69 out.resize(6);
70
71 out[0][0] = sign_[0]*in[0];
72 out[0][1] = sign_[0]*(in[1] - 1.0);
73 out[1][0] = sign_[1]*(in[0] - 1.0);
74 out[1][1] = sign_[1]*in[1];
75 out[2][0] = sign_[2]*in[0];
76 out[2][1] = sign_[2]*in[1];
77 out[3][0] = 3.0*in[0];
78 out[3][1] = 3.0 - 6.0*in[0] - 3.0*in[1];
79 out[4][0] = -3.0 + 3.0*in[0] + 6.0*in[1];
80 out[4][1] = -3.0*in[1];
81 out[5][0] = -3.0*in[0];
82 out[5][1] = 3.0*in[1];
83 }
84
91 inline void evaluateJacobian (const typename Traits::DomainType& in,
92 std::vector<typename Traits::JacobianType>& out) const
93 {
94 out.resize(6);
95
96 out[0][0][0] = sign_[0];
97 out[0][0][1] = 0.0;
98 out[0][1][0] = 0.0;
99 out[0][1][1] = sign_[0];
100
101 out[1][0][0] = sign_[1];
102 out[1][0][1] = 0.0;
103 out[1][1][0] = 0.0;
104 out[1][1][1] = sign_[1];
105
106 out[2][0][0] = sign_[2];
107 out[2][0][1] = 0.0;
108 out[2][1][0] = 0.0;
109 out[2][1][1] = sign_[2];
110
111 out[3][0][0] = 3.0;
112 out[3][0][1] = 0.0;
113 out[3][1][0] = -6.0;
114 out[3][1][1] = -3.0;
115
116 out[4][0][0] = 3.0;
117 out[4][0][1] = 6.0;
118 out[4][1][0] = 0.0;
119 out[4][1][1] = -3.0;
120
121 out[5][0][0] = -3.0;
122 out[5][0][1] = 0.0;
123 out[5][1][0] = 0.0;
124 out[5][1][1] = 3.0;
125 }
126
128 void partial (const std::array<unsigned int, 2>& order,
129 const typename Traits::DomainType& in, // position
130 std::vector<typename Traits::RangeType>& out) const // return value
131 {
132 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
133 if (totalOrder == 0) {
134 evaluateFunction(in, out);
135 } else if (totalOrder == 1) {
136 out.resize(size());
137 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
138
139 switch (direction) {
140 case 0:
141 out[0][0] = sign_[0];
142 out[0][1] = 0.0;
143
144 out[1][0] = sign_[1];
145 out[1][1] = 0.0;
146
147 out[2][0] = sign_[2];
148 out[2][1] = 0.0;
149
150 out[3][0] = 3.0;
151 out[3][1] = -6.0;
152
153 out[4][0] = 3.0;
154 out[4][1] = 0.0;
155
156 out[5][0] = -3.0;
157 out[5][1] = 0.0;
158 break;
159 case 1:
160 out[0][0] = 0.0;
161 out[0][1] = sign_[0];
162
163 out[1][0] = 0.0;
164 out[1][1] = sign_[1];
165
166 out[2][0] = 0.0;
167 out[2][1] = sign_[2];
168
169 out[3][0] = 0.0;
170 out[3][1] = -3.0;
171
172 out[4][0] = 6.0;
173 out[4][1] = -3.0;
174
175 out[5][0] = 0.0;
176 out[5][1] = 3.0;
177 break;
178 default:
179 DUNE_THROW(RangeError, "Component out of range.");
180 }
181 } else {
182 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
183 }
184 }
185
187 unsigned int order () const
188 {
189 return 1;
190 }
191
192 private:
193 std::array<R,3> sign_;
194 };
195}
196#endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALBASIS_HH
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:30
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:91
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:55
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:66
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: brezzidouglasmarini1simplex2dlocalbasis.hh:128
BDM1Simplex2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:37
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:187
BDM1Simplex2DLocalBasis(std::bitset< 3 > s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalbasis.hh:48
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Default exception for dummy implementations.
Definition: exceptions.hh:263
Default exception class for range errors.
Definition: exceptions.hh:254
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
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 (Nov 21, 23:30, 2024)