dune-localfunctions 2.11
Loading...
Searching...
No Matches
raviartthomas03dlocalbasis.hh
Go to the documentation of this file.
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_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
7
8#include <numeric>
9
10#include <dune/common/fmatrix.hh>
11
13
14namespace Dune
15{
25 template<class D, class R>
27 {
28 public:
29 typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,3,Dune::FieldVector<R,3>,
30 Dune::FieldMatrix<R,3,3> > Traits;
31
33 RT03DLocalBasis (std::bitset<4> s = 0)
34 {
35 for (int i=0; i<4; i++)
36 sign_[i] = s[i] ? -1.0 : 1.0;
37 }
38
40 unsigned int size () const
41 {
42 return 4;
43 }
44
46 inline void evaluateFunction (const typename Traits::DomainType& in,
47 std::vector<typename Traits::RangeType>& out) const
48 {
49 out.resize(4);
50 out[0] = {sign_[0]*2.0* in[0], sign_[0]*2.0* in[1], sign_[0]*2.0*(in[2]-D(1))};
51 out[1] = {sign_[1]*2.0* in[0], sign_[1]*2.0*(in[1]-D(1)), sign_[1]*2.0* in[2] };
52 out[2] = {sign_[2]*2.0*(in[0]-D(1)), sign_[2]*2.0* in[1], sign_[2]*2.0* in[2] };
53 out[3] = {sign_[3]*2.0* in[0], sign_[3]*2.0* in[1], sign_[3]*2.0* in[2] };
54 }
55
57 inline void
58 evaluateJacobian (const typename Traits::DomainType& in, // position
59 std::vector<typename Traits::JacobianType>& out) const // return value
60 {
61 out.resize(4);
62 for (int i=0; i<4; i++)
63 {
64 out[i][0] = {2.0*sign_[i], 0, 0};
65 out[i][1] = { 0,2.0*sign_[i], 0};
66 out[i][2] = { 0, 0,2.0*sign_[i]};
67 }
68 }
69
71 void partial (const std::array<unsigned int, 3>& order,
72 const typename Traits::DomainType& in, // position
73 std::vector<typename Traits::RangeType>& out) const // return value
74 {
75 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
76 if (totalOrder == 0) {
77 evaluateFunction(in, out);
78 } else if (totalOrder == 1) {
79 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
80 out.resize(size());
81
82 for (std::size_t i=0; i<size(); i++)
83 {
84 out[i][direction] = sign_[i]*2.0;
85 out[i][(direction+1)%3] = 0;
86 out[i][(direction+2)%3] = 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 < 3; ++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 face normals
106 std::array<R,4> sign_;
107 };
108}
109#endif
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
RT03DLocalBasis(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition raviartthomas03dlocalbasis.hh:33
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition raviartthomas03dlocalbasis.hh:71
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas03dlocalbasis.hh:98
unsigned int size() const
number of shape functions
Definition raviartthomas03dlocalbasis.hh:40
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition raviartthomas03dlocalbasis.hh:58
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition raviartthomas03dlocalbasis.hh:46
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > > Traits
Definition raviartthomas03dlocalbasis.hh:30