Dune Core Modules (2.6.0)

raviartthomas4cube2dlocalinterpolation.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_RAVIARTTHOMAS4_CUBE2D_LOCALINTERPOLATION_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS4_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(60);
84 fill(out.begin(), out.end(), 0.0);
85
86 const int qOrder = 12;
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 out[4] += (y[0]*n0[0] + y[1]*n0[1])*(1.0-20.0*qPos+90.0*pow(qPos,2)-140.0*pow(qPos,3)+70.0*pow(qPos,4))*it->weight()*sign0;
102
103 localPos[0] = 1.0;
104 localPos[1] = qPos;
105 f.evaluate(localPos, y);
106 out[5] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
107 out[6] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
108 out[7] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
109 out[8] += (y[0]*n1[0] + y[1]*n1[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
110 out[9] += (y[0]*n1[0] + y[1]*n1[1])*(1.0-20.0*qPos+90.0*pow(qPos,2)-140.0*pow(qPos,3)+70.0*pow(qPos,4))*it->weight()*sign1;
111
112 localPos[0] = qPos;
113 localPos[1] = 0.0;
114 f.evaluate(localPos, y);
115 out[10] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
116 out[11] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
117 out[12] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
118 out[13] += (y[0]*n2[0] + y[1]*n2[1])*(-20.0*qPos*qPos*qPos + 30.0*qPos*qPos - 12.0*qPos + 1.0)*it->weight();
119 out[14] += (y[0]*n2[0] + y[1]*n2[1])*(1.0-20.0*qPos+90.0*pow(qPos,2)-140.0*pow(qPos,3)+70.0*pow(qPos,4))*it->weight()*sign2;
120
121 localPos[0] = qPos;
122 localPos[1] = 1.0;
123 f.evaluate(localPos, y);
124 out[15] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
125 out[16] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
126 out[17] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
127 out[18] += (y[0]*n3[0] + y[1]*n3[1])*(20.0*qPos*qPos*qPos - 30.0*qPos*qPos + 12.0*qPos - 1.0)*it->weight();
128 out[19] += (y[0]*n3[0] + y[1]*n3[1])*(1.0-20.0*qPos+90.0*pow(qPos,2)-140.0*pow(qPos,3)+70.0*pow(qPos,4))*it->weight()*sign3;
129 }
130
132
133 for (typename QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
134 it != rule2.end(); ++it)
135 {
136 FieldVector<double,2> qPos = it->position();
137
138 f.evaluate(qPos, y);
139 std::vector<std::vector<double> > l(2,std::vector<double> (5));
140 l[0][0]=1.0;
141 l[1][0]=1.0;
142 l[0][1]=2.0*qPos[0]-1.0;
143 l[1][1]=2.0*qPos[1]-1.0;
144 l[0][2]=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
145 l[1][2]=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
146 l[0][3]=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
147 l[1][3]=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
148 l[0][4]=1.0-20.0*qPos[0]+90.0*pow(qPos[0],2)-140.0*pow(qPos[0],3)+70.0*pow(qPos[0],4);
149 l[1][4]=1.0-20.0*qPos[1]+90.0*pow(qPos[1],2)-140.0*pow(qPos[1],3)+70.0*pow(qPos[1],4);
150
151 for (int i=0;i<4;i++)
152 for (int j=0;j<5;j++)
153 out[20+i*5+j] +=y[0]*l[0][i]*l[1][j]*it->weight();
154
155 for (int i=0;i<5;i++)
156 for (int j=0;j<4;j++)
157 out[40+i*4+j] +=y[1]*l[0][i]*l[1][j]*it->weight();
158 }
159 }
160
161 private:
162 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
163 typename LB::Traits::DomainType n0, n1, n2, n3;
164 };
165}
166
167#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 triangle.
Definition: raviartthomas4cube2dlocalinterpolation.hh:23
RT4Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas4cube2dlocalinterpolation.hh:27
void interpolate(const F &f, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas4cube2dlocalinterpolation.hh:76
RT4Cube2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas4cube2dlocalinterpolation.hh:37
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 (Jul 15, 22:36, 2024)