Dune Core Modules (2.6.0)

raviartthomas12dlocalinterpolation.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_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
5
6#include <vector>
7
9
10namespace Dune
11{
12
21 template<class LB>
23 {
24
25 public:
28 {
29 sign0 = sign1 = sign2 = 1.0;
30 }
31
37 RT12DLocalInterpolation (unsigned int s)
38 {
39 sign0 = sign1 = sign2 = 1.0;
40 if (s & 1)
41 {
42 sign0 = -1.0;
43 }
44 if (s & 2)
45 {
46 sign1 = -1.0;
47 }
48 if (s & 4)
49 {
50 sign2 = -1.0;
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 typedef typename LB::Traits::DomainFieldType Vector;
77 typename F::Traits::RangeType y;
78
79 out.resize(8);
80 fill(out.begin(), out.end(), 0.0);
81
82 const int qOrder1 = 4;
84
85 for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
86 it != rule1.end(); ++it)
87 {
88 Scalar qPos = it->position();
89 typename LB::Traits::DomainType localPos;
90
91 localPos[0] = qPos;
92 localPos[1] = 0.0;
93 f.evaluate(localPos, y);
94 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
95 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
96
97 localPos[0] = 0.0;
98 localPos[1] = qPos;
99 f.evaluate(localPos, y);
100 out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
101 out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
102
103 localPos[0] = 1.0 - qPos;
104 localPos[1] = qPos;
105 f.evaluate(localPos, y);
106 out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
107 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
108 }
109
110 const int qOrder2 = 8;
112
113 for (typename Dune::QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
114 it != rule2.end(); ++it)
115 {
116 Dune::FieldVector<double,2> qPos = it->position();
117
118 f.evaluate(qPos, y);
119 out[6] += y[0]*it->weight();
120 out[7] += y[1]*it->weight();
121 }
122 }
123
124 private:
125 typename LB::Traits::RangeFieldType sign0,sign1,sign2;
126 typename LB::Traits::DomainType n0,n1,n2;
127 typename LB::Traits::RangeFieldType c0,c1,c2;
128 };
129}
130#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:93
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
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas12dlocalinterpolation.hh:23
RT12DLocalInterpolation()
Standard constructor.
Definition: raviartthomas12dlocalinterpolation.hh:27
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas12dlocalinterpolation.hh:72
RT12DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalinterpolation.hh:37
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 (Jul 15, 22:36, 2024)