Dune Core Modules (2.10.0)

raviartthomas0prismlocalinterpolation.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_LOCALFUNCTIONS_RAVIARTTHOMAS0_PRISM_LOCALINTERPOLATION_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PRISM_LOCALINTERPOLATION_HH
7
8#include <vector>
9
10namespace Dune
11{
20 template<class LB>
22 {
23
24 public:
25
31 RT0PrismLocalInterpolation (std::bitset<5> s = 0)
32 {
33 typedef typename LB::Traits::RangeFieldType Scalar;
34
35 for (size_t i=0; i<5; i++)
36 sign[i] = (s[i]) ? -1.0 : 1.0;
37
38 Scalar r = 1/std::sqrt(2);
39
40 n[0] = { 0.0, -1.0, 0.0};
41 n[1] = {-1.0, 0.0, 0.0};
42 n[2] = { r, r, 0.0};
43 n[3] = { 0.0, 0.0, -1.0};
44 n[4] = { 0.0, 0.0, 1.0};
45
46 c[0] = 1.0;
47 c[1] = 1.0;
48 c[2] = std::sqrt(2);
49 c[3] = 1/2.0;
50 c[4] = 1/2.0;
51
52 m[0] = { 0.5, 0.0, 0.5};
53 m[1] = { 0.0, 0.5, 0.5};
54 m[2] = { 0.5, 0.5, 0.5};
55 m[3] = { 1/3.0, 1/3.0, 0.0};
56 m[4] = { 1/3.0, 1/3.0, 1.0};
57 }
58
67 template<class F, class C>
68 void interpolate (const F& f, std::vector<C>& out) const
69 {
70 out.resize(5);
71 for(int i=0; i<5; i++)
72 out[i] = f(m[i]).dot(n[i]) * c[i] * sign[i];
73 }
74
75 private:
76 // Facet orientations
77 std::array<typename LB::Traits::RangeFieldType, 5> sign;
78 // Facet area
79 std::array<typename LB::Traits::RangeFieldType, 5> c;
80
81 // Facet normals
82 std::array<typename LB::Traits::DomainType, 5> n;
83 // Facet midpoints
84 std::array<typename LB::Traits::DomainType, 5> m;
85 };
86}
87#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PRISM_LOCALINTERPOLATION_HH
First order Raviart-Thomas shape functions on the reference prism.
Definition: raviartthomas0prismlocalinterpolation.hh:22
RT0PrismLocalInterpolation(std::bitset< 5 > s=0)
Make set number s, where 0 <= s < 32.
Definition: raviartthomas0prismlocalinterpolation.hh:31
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas0prismlocalinterpolation.hh:68
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 27, 23:30, 2024)