Dune Core Modules (2.7.1)

raviartthomas1cube2dlocalinterpolation.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_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_HH
5
6#include <vector>
7
9#include <dune/localfunctions/common/localinterpolation.hh>
10
11
12namespace Dune
13{
14
23 template<class LB>
25 {
26
27 public:
33 RT1Cube2DLocalInterpolation (unsigned int s = 0)
34 {
35 sign0 = sign1 = sign2 = sign3 = 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 if (s & 8)
49 {
50 sign3 = -1.0;
51 }
52
53 n0[0] = -1.0;
54 n0[1] = 0.0;
55 n1[0] = 1.0;
56 n1[1] = 0.0;
57 n2[0] = 0.0;
58 n2[1] = -1.0;
59 n3[0] = 0.0;
60 n3[1] = 1.0;
61 }
62
71 template<class F, class C>
72 void interpolate (const F& ff, 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
78 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
79
80 out.resize(12);
81 fill(out.begin(), out.end(), 0.0);
82
83 const int qOrder = 3;
85
86 for (typename QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
87 it != rule1.end(); ++it)
88 {
89 Scalar qPos = it->position();
90 typename LB::Traits::DomainType localPos;
91
92 localPos[0] = 0.0;
93 localPos[1] = qPos;
94 auto y = f(localPos);
95 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
96 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
97
98 localPos[0] = 1.0;
99 localPos[1] = qPos;
100 y = f(localPos);
101 out[2] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
102 out[3] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
103
104 localPos[0] = qPos;
105 localPos[1] = 0.0;
106 y = f(localPos);
107 out[4] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
108 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
109
110 localPos[0] = qPos;
111 localPos[1] = 1.0;
112 y = f(localPos);
113 out[6] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
114 out[7] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
115 }
116
118
119 for (typename QuadratureRule<Vector,2>::const_iterator it=rule2.begin(); it!=rule2.end(); ++it)
120 {
121 FieldVector<double,2> qPos = it->position();
122
123 auto y = f(qPos);
124 out[8] += y[0]*it->weight();
125 out[9] += y[1]*it->weight();
126 out[10] += y[0]*qPos[1]*it->weight();
127 out[11] += y[1]*qPos[0]*it->weight();
128 }
129 }
130
131 private:
132 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
133 typename LB::Traits::DomainType n0, n1, n2, n3;
134 };
135}
136#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_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: raviartthomas1cube2dlocalinterpolation.hh:25
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube2dlocalinterpolation.hh:72
RT1Cube2DLocalInterpolation(unsigned int s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalinterpolation.hh:33
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:775
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 (Nov 13, 23:29, 2024)