Dune Core Modules (2.7.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 #include <dune/localfunctions/common/localinterpolation.hh>
10 
11 namespace Dune
12 {
13 
22  template<class LB>
24  {
25 
26  public:
29  {
30  sign0 = sign1 = sign2 = sign3 = 1.0;
31  }
32 
38  RT4Cube2DLocalInterpolation (unsigned int s)
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<typename F, typename C>
77  void interpolate (const F& ff, 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 
83  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
84 
85  out.resize(60);
86  fill(out.begin(), out.end(), 0.0);
87 
88  const int qOrder = 12;
90 
91  for (typename QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.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  auto y = f(localPos);
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  out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0;
102  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();
103  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;
104 
105  localPos[0] = 1.0;
106  localPos[1] = qPos;
107  y = f(localPos);
108  out[5] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
109  out[6] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
110  out[7] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
111  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();
112  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;
113 
114  localPos[0] = qPos;
115  localPos[1] = 0.0;
116  y = f(localPos);
117  out[10] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
118  out[11] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
119  out[12] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
120  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();
121  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;
122 
123  localPos[0] = qPos;
124  localPos[1] = 1.0;
125  y = f(localPos);
126  out[15] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
127  out[16] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
128  out[17] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
129  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();
130  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;
131  }
132 
134 
135  for (typename QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
136  it != rule2.end(); ++it)
137  {
138  FieldVector<double,2> qPos = it->position();
139 
140  auto y = f(qPos);
141  std::vector<std::vector<double> > l(2,std::vector<double> (5));
142  l[0][0]=1.0;
143  l[1][0]=1.0;
144  l[0][1]=2.0*qPos[0]-1.0;
145  l[1][1]=2.0*qPos[1]-1.0;
146  l[0][2]=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
147  l[1][2]=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
148  l[0][3]=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
149  l[1][3]=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
150  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);
151  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);
152 
153  for (int i=0;i<4;i++)
154  for (int j=0;j<5;j++)
155  out[20+i*5+j] +=y[0]*l[0][i]*l[1][j]*it->weight();
156 
157  for (int i=0;i<5;i++)
158  for (int j=0;j<4;j++)
159  out[40+i*4+j] +=y[1]*l[0][i]*l[1][j]*it->weight();
160  }
161  }
162 
163  private:
164  typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
165  typename LB::Traits::DomainType n0, n1, n2, n3;
166  };
167 }
168 
169 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS3_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
Second order Raviart-Thomas shape functions on the reference triangle.
Definition: raviartthomas4cube2dlocalinterpolation.hh:24
RT4Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas4cube2dlocalinterpolation.hh:28
RT4Cube2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas4cube2dlocalinterpolation.hh:38
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas4cube2dlocalinterpolation.hh:77
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.80.0 (May 16, 22:29, 2024)