Dune Core Modules (2.7.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 #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:
27 
33  RT4Cube2DLocalInterpolation (unsigned int s = 0)
34  {
35  sign0 = sign1 = sign2 = sign3 = 1.0;
36  if (s & 1)
37  {
38  sign0 *= -1.0;
39  }
40  if (s & 2)
41  {
42  sign1 *= -1.0;
43  }
44  if (s & 4)
45  {
46  sign2 *= -1.0;
47  }
48  if (s & 8)
49  {
50  sign3 *= -1.0;
51  }
52 
53  n0[0] = -1.0;
54  n0[1] = 0.0;
55  n1[0] = 1.0;
56  n1[1] = 0.0;
57  n2[0] = 0.0;
58  n2[1] = -1.0;
59  n3[0] = 0.0;
60  n3[1] = 1.0;
61  }
62 
71  template<typename F, typename C>
72  void interpolate (const F& ff, std::vector<C>& out) const
73  {
74  // f gives v*outer normal at a point on the edge!
75  typedef typename LB::Traits::RangeFieldType Scalar;
76  typedef typename LB::Traits::DomainFieldType Vector;
77 
78  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
79 
80  out.resize(60);
81  fill(out.begin(), out.end(), 0.0);
82 
83  const int qOrder = 12;
85 
86  for (typename QuadratureRule<Scalar,1>::const_iterator it=rule.begin(); it!=rule.end(); ++it)
87  {
88  Scalar qPos = it->position();
89  typename LB::Traits::DomainType localPos;
90 
91  localPos[0] = 0.0;
92  localPos[1] = qPos;
93  auto y = f(localPos);
94  out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
95  out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
96  out[2] += (y[0]*n0[0] + y[1]*n0[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign0;
97  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();
98  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;
99 
100  localPos[0] = 1.0;
101  localPos[1] = qPos;
102  y = f(localPos);
103  out[5] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
104  out[6] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
105  out[7] += (y[0]*n1[0] + y[1]*n1[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign1;
106  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();
107  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;
108 
109  localPos[0] = qPos;
110  localPos[1] = 0.0;
111  y = f(localPos);
112  out[10] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
113  out[11] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
114  out[12] += (y[0]*n2[0] + y[1]*n2[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign2;
115  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();
116  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;
117 
118  localPos[0] = qPos;
119  localPos[1] = 1.0;
120  y = f(localPos);
121  out[15] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
122  out[16] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
123  out[17] += (y[0]*n3[0] + y[1]*n3[1])*(6.0*qPos*qPos - 6.0*qPos + 1.0)*it->weight()*sign3;
124  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();
125  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;
126  }
127 
129 
130  for (typename QuadratureRule<Vector,2>::const_iterator it = rule2.begin();
131  it != rule2.end(); ++it)
132  {
133  FieldVector<double,2> qPos = it->position();
134 
135  auto y = f(qPos);
136  std::vector<std::vector<double> > l(2,std::vector<double> (5));
137  l[0][0]=1.0;
138  l[1][0]=1.0;
139  l[0][1]=2.0*qPos[0]-1.0;
140  l[1][1]=2.0*qPos[1]-1.0;
141  l[0][2]=6.0*qPos[0]*qPos[0]-6.0*qPos[0]+1.0;
142  l[1][2]=6.0*qPos[1]*qPos[1]-6.0*qPos[1]+1.0;
143  l[0][3]=20.0*qPos[0]*qPos[0]*qPos[0] - 30.0*qPos[0]*qPos[0] + 12.0*qPos[0] - 1.0;
144  l[1][3]=20.0*qPos[1]*qPos[1]*qPos[1] - 30.0*qPos[1]*qPos[1] + 12.0*qPos[1] - 1.0;
145  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);
146  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);
147 
148  for (int i=0;i<4;i++)
149  for (int j=0;j<5;j++)
150  out[20+i*5+j] +=y[0]*l[0][i]*l[1][j]*it->weight();
151 
152  for (int i=0;i<5;i++)
153  for (int j=0;j<4;j++)
154  out[40+i*4+j] +=y[1]*l[0][i]*l[1][j]*it->weight();
155  }
156  }
157 
158  private:
159  typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
160  typename LB::Traits::DomainType n0, n1, n2, n3;
161  };
162 }
163 
164 #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(unsigned int s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas4cube2dlocalinterpolation.hh:33
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas4cube2dlocalinterpolation.hh:72
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)