Dune Core Modules (2.6.0)

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
10
11namespace Dune
12{
13
22 template<class LB>
24 {
25
26 public:
29 {
30 sign0 = sign1 = sign2 = sign3 = 1.0;
31 }
32
39 {
40 sign0 = sign1 = sign2 = sign3 = 1.0;
41 if (s & 1)
42 {
43 sign0 = -1.0;
44 }
45 if (s & 2)
46 {
47 sign1 = -1.0;
48 }
49 if (s & 4)
50 {
51 sign2 = -1.0;
52 }
53 if (s & 8)
54 {
55 sign3 = -1.0;
56 }
57
58 n0[0] = -1.0;
59 n0[1] = 0.0;
60 n1[0] = 1.0;
61 n1[1] = 0.0;
62 n2[0] = 0.0;
63 n2[1] = -1.0;
64 n3[0] = 0.0;
65 n3[1] = 1.0;
66 }
67
76 template<class F, class C>
77 void interpolate (const F& f, std::vector<C>& out) const
78 {
79 // f gives v*outer normal at a point on the edge!
80 typedef typename LB::Traits::RangeFieldType Scalar;
81 typedef typename LB::Traits::DomainFieldType Vector;
82 typename F::Traits::RangeType y;
83
84 out.resize(12);
85 fill(out.begin(), out.end(), 0.0);
86
87 const int qOrder = 3;
89
90 for (typename QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
91 it != rule1.end(); ++it)
92 {
93 Scalar qPos = it->position();
94 typename LB::Traits::DomainType localPos;
95
96 localPos[0] = 0.0;
97 localPos[1] = qPos;
98 f.evaluate(localPos, y);
99 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
100 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
101
102 localPos[0] = 1.0;
103 localPos[1] = qPos;
104 f.evaluate(localPos, y);
105 out[2] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
106 out[3] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
107
108 localPos[0] = qPos;
109 localPos[1] = 0.0;
110 f.evaluate(localPos, y);
111 out[4] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
112 out[5] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
113
114 localPos[0] = qPos;
115 localPos[1] = 1.0;
116 f.evaluate(localPos, y);
117 out[6] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
118 out[7] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
119 }
120
122
123 for (typename QuadratureRule<Vector,2>::const_iterator it=rule2.begin(); it!=rule2.end(); ++it)
124 {
125 FieldVector<double,2> qPos = it->position();
126
127 f.evaluate(qPos, y);
128 out[8] += y[0]*it->weight();
129 out[9] += y[1]*it->weight();
130 out[10] += y[0]*qPos[1]*it->weight();
131 out[11] += y[1]*qPos[0]*it->weight();
132 }
133 }
134
135 private:
136 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
137 typename LB::Traits::DomainType n0, n1, n2, n3;
138 };
139}
140#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_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: raviartthomas1cube2dlocalinterpolation.hh:24
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube2dlocalinterpolation.hh:77
RT1Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas1cube2dlocalinterpolation.hh:28
RT1Cube2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalinterpolation.hh:38
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:705
Dune namespace.
Definition: alignedallocator.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)