Dune Core Modules (2.6.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#ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
5
6#include <vector>
7
9
10namespace Dune
11{
20 template<class LB>
22 {
23
24 public:
27 {
28 sign0 = sign1 = sign2 = 1.0;
29 }
30
37 {
38 sign0 = sign1 = sign2 = 1.0;
39 if (s & 1)
40 {
41 sign0 = -1.0;
42 }
43 if (s & 2)
44 {
45 sign1 = -1.0;
46 }
47 if (s & 4)
48 {
49 sign2 = -1.0;
50 }
51
52 n0[0] = 0.0;
53 n0[1] = -1.0;
54 n1[0] = -1.0;
55 n1[1] = 0.0;
56 n2[0] = 1.0/sqrt(2.0);
57 n2[1] = 1.0/sqrt(2.0);
58 c0 = 0.5*n0[0] - 1.0*n0[1];
59 c1 = -1.0*n1[0] + 0.5*n1[1];
60 c2 = 0.5*n2[0] + 0.5*n2[1];
61 }
62
71 template<typename F, typename C>
72 void interpolate (const F& f, std::vector<C>& out) const
73 {
74 // f gives v*outer normal at a point on the edge!
75 typedef typename LB::Traits::RangeFieldType Scalar;
76 typename F::Traits::RangeType y;
77
78 out.resize(6);
79 fill(out.begin(), out.end(), 0.0);
80
81 const int qOrder = 4;
83
84 for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
85 {
86 Scalar qPos = it->position();
87 typename LB::Traits::DomainType localPos;
88
89 localPos[0] = qPos;
90 localPos[1] = 0.0;
91 f.evaluate(localPos, y);
92 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
93 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
94
95 localPos[0] = 0.0;
96 localPos[1] = qPos;
97 f.evaluate(localPos, y);
98 out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
99 out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
100
101 localPos[0] = 1.0 - qPos;
102 localPos[1] = qPos;
103 f.evaluate(localPos, y);
104 out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
105 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
106 }
107 }
108
109 private:
110 typename LB::Traits::RangeFieldType sign0,sign1,sign2;
111 typename LB::Traits::DomainType n0,n1,n2;
112 typename LB::Traits::RangeFieldType c0,c1,c2;
113 };
114}
115
116#endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_SIMPLEX2D_LOCALINTERPOLATION_HH
First order Brezzi-Douglas-Marini shape functions on the reference triangle.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:22
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:72
BDM1Simplex2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:36
BDM1Simplex2DLocalInterpolation()
Standard constructor.
Definition: brezzidouglasmarini1simplex2dlocalinterpolation.hh:26
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:97
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:225
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:696
Dune namespace.
Definition: alignedallocator.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)