Dune Core Modules (2.9.1)

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// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS4_CUBE2D_LOCALINTERPOLATION_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS4_CUBE2D_LOCALINTERPOLATION_HH
7
8#include <vector>
9
11#include <dune/localfunctions/common/localinterpolation.hh>
12
13namespace Dune
14{
15
24 template<class LB>
26 {
27
28 public:
29
35 RT4Cube2DLocalInterpolation (unsigned int s = 0)
36 {
37 sign0 = sign1 = sign2 = sign3 = 1.0;
38 if (s & 1)
39 {
40 sign0 *= -1.0;
41 }
42 if (s & 2)
43 {
44 sign1 *= -1.0;
45 }
46 if (s & 4)
47 {
48 sign2 *= -1.0;
49 }
50 if (s & 8)
51 {
52 sign3 *= -1.0;
53 }
54
55 n0[0] = -1.0;
56 n0[1] = 0.0;
57 n1[0] = 1.0;
58 n1[1] = 0.0;
59 n2[0] = 0.0;
60 n2[1] = -1.0;
61 n3[0] = 0.0;
62 n3[1] = 1.0;
63 }
64
73 template<typename F, typename C>
74 void interpolate (const F& ff, std::vector<C>& out) const
75 {
76 // f gives v*outer normal at a point on the edge!
77 typedef typename LB::Traits::RangeFieldType Scalar;
78 typedef typename LB::Traits::DomainFieldType Vector;
79
80 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
81
82 out.resize(60);
83 fill(out.begin(), out.end(), 0.0);
84
85 const int qOrder = 12;
87
88 for (typename QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
89 {
90 Scalar qPos = it->position();
91 typename LB::Traits::DomainType localPos;
92
93 localPos[0] = 0.0;
94 localPos[1] = qPos;
95 auto y = f(localPos);
96 out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
97 out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
98 out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0;
99 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();
100 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;
101
102 localPos[0] = 1.0;
103 localPos[1] = qPos;
104 y = f(localPos);
105 out[5] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
106 out[6] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
107 out[7] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
108 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();
109 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;
110
111 localPos[0] = qPos;
112 localPos[1] = 0.0;
113 y = f(localPos);
114 out[10] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
115 out[11] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
116 out[12] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
117 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();
118 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;
119
120 localPos[0] = qPos;
121 localPos[1] = 1.0;
122 y = f(localPos);
123 out[15] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
124 out[16] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
125 out[17] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
126 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();
127 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;
128 }
129
131
132 for (typename QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
133 it != rule2.end(); ++it)
134 {
135 FieldVector<double,2> qPos = it->position();
136
137 auto y = f(qPos);
138 std::vector<std::vector<double> > l(2,std::vector<double> (5));
139 l[0][0]=1.0;
140 l[1][0]=1.0;
141 l[0][1]=2.0*qPos[0]-1.0;
142 l[1][1]=2.0*qPos[1]-1.0;
143 l[0][2]=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
144 l[1][2]=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
145 l[0][3]=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
146 l[1][3]=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
147 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);
148 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);
149
150 for (int i=0;i<4;i++)
151 for (int j=0;j<5;j++)
152 out[20+i*5+j] +=y[0]*l[0][i]*l[1][j]*it->weight();
153
154 for (int i=0;i<5;i++)
155 for (int j=0;j<4;j++)
156 out[40+i*4+j] +=y[1]*l[0][i]*l[1][j]*it->weight();
157 }
158 }
159
160 private:
161 typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
162 typename LB::Traits::DomainType n0, n1, n2, n3;
163 };
164}
165
166#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_CUBE2D_LOCALINTERPOLATION_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:154
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:266
Second order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas4cube2dlocalinterpolation.hh:26
RT4Cube2DLocalInterpolation(unsigned int s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas4cube2dlocalinterpolation.hh:35
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas4cube2dlocalinterpolation.hh:74
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:473
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)