DUNE PDELab (2.7)

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#include <dune/localfunctions/common/localinterpolation.hh>
10
11namespace Dune
12{
13
22 template<class LB>
24 {
25
26 public:
27
33 RT12DLocalInterpolation (unsigned int s = 0)
34 {
35 sign0 = sign1 = sign2 = 1.0;
36 if (s & 1)
37 {
38 sign0 = -1.0;
39 }
40 if (s & 2)
41 {
42 sign1 = -1.0;
43 }
44 if (s & 4)
45 {
46 sign2 = -1.0;
47 }
48 n0[0] = 0.0;
49 n0[1] = -1.0;
50 n1[0] = -1.0;
51 n1[1] = 0.0;
52 n2[0] = 1.0/sqrt(2.0);
53 n2[1] = 1.0/sqrt(2.0);
54 c0 = 0.5*n0[0] - 1.0*n0[1];
55 c1 = -1.0*n1[0] + 0.5*n1[1];
56 c2 = 0.5*n2[0] + 0.5*n2[1];
57 }
58
67 template<typename F, typename C>
68 void interpolate (const F& ff, std::vector<C>& out) const
69 {
70 // f gives v*outer normal at a point on the edge!
71 typedef typename LB::Traits::RangeFieldType Scalar;
72 typedef typename LB::Traits::DomainFieldType Vector;
73
74 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
75
76 out.resize(8);
77 fill(out.begin(), out.end(), 0.0);
78
79 const int qOrder1 = 4;
81
82 for (typename Dune::QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
83 it != rule1.end(); ++it)
84 {
85 Scalar qPos = it->position();
86 typename LB::Traits::DomainType localPos;
87
88 localPos[0] = qPos;
89 localPos[1] = 0.0;
90 auto y = f(localPos);
91 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0/c0;
92 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight()/c0;
93
94 localPos[0] = 0.0;
95 localPos[1] = qPos;
96 y = f(localPos);
97 out[1] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1/c1;
98 out[4] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight()/c1;
99
100 localPos[0] = 1.0 - qPos;
101 localPos[1] = qPos;
102 y = f(localPos);
103 out[2] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2/c2;
104 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(2.0*qPos - 1.0)*it->weight()/c2;
105 }
106
107 const int qOrder2 = 8;
109
110 for (typename Dune::QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
111 it != rule2.end(); ++it)
112 {
113 Dune::FieldVector<double,2> qPos = it->position();
114
115 auto y = f(qPos);
116 out[6] += y[0]*it->weight();
117 out[7] += y[1]*it->weight();
118 }
119 }
120
121 private:
122 typename LB::Traits::RangeFieldType sign0,sign1,sign2;
123 typename LB::Traits::DomainType n0,n1,n2;
124 typename LB::Traits::RangeFieldType c0,c1,c2;
125 };
126}
127#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALINTERPOLATION_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:126
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:254
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas12dlocalinterpolation.hh:24
RT12DLocalInterpolation(unsigned int s=0)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalinterpolation.hh:33
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas12dlocalinterpolation.hh:68
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:767
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:233
Dune namespace.
Definition: alignedallocator.hh:14
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)