Dune Core Modules (2.6.0)

raviartthomas3cube2dlocalinterpolation.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_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
5
6#include <vector>
7
9
10namespace Dune
11{
12
21 template<class LB>
23 {
24
25 public:
28 {
29 sign0 = sign1 = sign2 = sign3 = 1.0;
30 }
31
38 {
39 sign0 = sign1 = sign2 = sign3 = 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 if (s & 8)
53 {
54 sign3 *= -1.0;
55 }
56
57 n0[0] = -1.0;
58 n0[1] = 0.0;
59 n1[0] = 1.0;
60 n1[1] = 0.0;
61 n2[0] = 0.0;
62 n2[1] = -1.0;
63 n3[0] = 0.0;
64 n3[1] = 1.0;
65 }
66
75 template<typename F, typename C>
76 void interpolate (const F& f, std::vector<C>& out) const
77 {
78 // f gives v*outer normal at a point on the edge!
79 typedef typename LB::Traits::RangeFieldType Scalar;
80 typedef typename LB::Traits::DomainFieldType Vector;
81 typename F::Traits::RangeType y;
82
83 out.resize(40);
84 fill(out.begin(), out.end(), 0.0);
85
86 const int qOrder = 9;
88
89 for (typename QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
90 {
91 Scalar qPos = it->position();
92 typename LB::Traits::DomainType localPos;
93
94 localPos[0] = 0.0;
95 localPos[1] = qPos;
96 f.evaluate(localPos, y);
97 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
98 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
99 out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0;
100 out[3] += (y[0]*n0[0] + y[1]*n0[1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*it->weight();
101
102 localPos[0] = 1.0;
103 localPos[1] = qPos;
104 f.evaluate(localPos, y);
105 out[4] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
106 out[5] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
107 out[6] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
108 out[7] += (y[0]*n1[0] + y[1]*n1[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
109
110 localPos[0] = qPos;
111 localPos[1] = 0.0;
112 f.evaluate(localPos, y);
113 out[8] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
114 out[9] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
115 out[10] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
116 out[11] += (y[0]*n2[0] + y[1]*n2[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
117
118 localPos[0] = qPos;
119 localPos[1] = 1.0;
120 f.evaluate(localPos, y);
121 out[12] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
122 out[13] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
123 out[14] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
124 out[15] += (y[0]*n3[0] + y[1]*n3[1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*it->weight();
125 }
126
128
129 for (typename QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
130 it != rule2.end(); ++it)
131 {
132 FieldVector<double,2> qPos = it->position();
133
134 f.evaluate(qPos, y);
135 double l0_x=1.0;
136 double l1_x=2.0*qPos[0]-1.0;
137 double l2_x=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
138 double l3_x=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
139 double l0_y=1.0;
140 double l1_y=2.0*qPos[1]-1.0;
141 double l2_y=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
142 double l3_y=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
143
144 out[16] += y[0]*l0_x*l0_y*it->weight();
145 out[17] += y[0]*l0_x*l1_y*it->weight();
146 out[18] += y[0]*l0_x*l2_y*it->weight();
147 out[19] += y[0]*l0_x*l3_y*it->weight();
148 out[20] += y[0]*l1_x*l0_y*it->weight();
149 out[21] += y[0]*l1_x*l1_y*it->weight();
150 out[22] += y[0]*l1_x*l2_y*it->weight();
151 out[23] += y[0]*l1_x*l3_y*it->weight();
152 out[24] += y[0]*l2_x*l0_y*it->weight();
153 out[25] += y[0]*l2_x*l1_y*it->weight();
154 out[26] += y[0]*l2_x*l2_y*it->weight();
155 out[27] += y[0]*l2_x*l3_y*it->weight();
156
157 out[28] += y[1]*l0_x*l0_y*it->weight();
158 out[29] += y[1]*l0_x*l1_y*it->weight();
159 out[30] += y[1]*l0_x*l2_y*it->weight();
160 out[31] += y[1]*l1_x*l0_y*it->weight();
161 out[32] += y[1]*l1_x*l1_y*it->weight();
162 out[33] += y[1]*l1_x*l2_y*it->weight();
163 out[34] += y[1]*l2_x*l0_y*it->weight();
164 out[35] += y[1]*l2_x*l1_y*it->weight();
165 out[36] += y[1]*l2_x*l2_y*it->weight();
166 out[37] += y[1]*l3_x*l0_y*it->weight();
167 out[38] += y[1]*l3_x*l1_y*it->weight();
168 out[39] += y[1]*l3_x*l2_y*it->weight();
169 }
170 }
171
172 private:
173 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
174 typename LB::Traits::DomainType n0, n1, n2, n3;
175 };
176}
177
178#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_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
Second order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas3cube2dlocalinterpolation.hh:23
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas3cube2dlocalinterpolation.hh:76
RT3Cube2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas3cube2dlocalinterpolation.hh:37
RT3Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas3cube2dlocalinterpolation.hh:27
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)