Dune Core Modules (2.9.0)

brezzidouglasmarini1simplex2dlocalinterpolation.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_LOCALINTERPOLATION_HH
6#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
7
8#include <vector>
9
11#include <dune/localfunctions/common/localinterpolation.hh>
12
13namespace Dune
14{
23 template<class LB>
25 {
26
27 public:
30 {
31 sign0 = sign1 = sign2 = 1.0;
32 }
33
40 {
41 using std::sqrt;
42 sign0 = sign1 = sign2 = 1.0;
43 if (s & 1)
44 {
45 sign0 = -1.0;
46 }
47 if (s & 2)
48 {
49 sign1 = -1.0;
50 }
51 if (s & 4)
52 {
53 sign2 = -1.0;
54 }
55
56 n0[0] = 0.0;
57 n0[1] = -1.0;
58 n1[0] = -1.0;
59 n1[1] = 0.0;
60 n2[0] = 1.0/sqrt(2.0);
61 n2[1] = 1.0/sqrt(2.0);
62 c0 = 0.5*n0[0] - 1.0*n0[1];
63 c1 = -1.0*n1[0] + 0.5*n1[1];
64 c2 = 0.5*n2[0] + 0.5*n2[1];
65 }
66
75 template<typename F, typename C>
76 void interpolate (const F& ff, std::vector<C>& out) const
77 {
78 // f gives v*outer normal at a point on the edge!
79 typedef typename LB::Traits::RangeFieldType Scalar;
80
81 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
82
83 out.resize(6);
84 fill(out.begin(), out.end(), 0.0);
85
86 const int qOrder = 4;
88
89 for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
90 {
91 Scalar qPos = it->position();
92 typename LB::Traits::DomainType localPos;
93
94 localPos[0] = qPos;
95 localPos[1] = 0.0;
96 auto y = f(localPos);
97 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
98 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
99
100 localPos[0] = 0.0;
101 localPos[1] = qPos;
102 y = f(localPos);
103 out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
104 out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
105
106 localPos[0] = 1.0 - qPos;
107 localPos[1] = qPos;
108 y = f(localPos);
109 out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
110 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
111 }
112 }
113
114 private:
115 typename LB::Traits::RangeFieldType sign0,sign1,sign2;
116 typename LB::Traits::DomainType n0,n1,n2;
117 typename LB::Traits::RangeFieldType c0,c1,c2;
118 };
119}
120
121#endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:25
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:76
BDM1Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:39
BDM1Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:29
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:154
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:266
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:463
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 (Jul 15, 22:36, 2024)