dune-localfunctions  2.8.0
raviartthomas1cube3dlocalinterpolation.hh
Go to the documentation of this file.
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_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
5 
6 #include <vector>
7 
8 #include <dune/geometry/quadraturerules.hh>
10 
11 namespace Dune
12 {
21  template<class LB>
23  {
24 
25  public:
26 
32  RT1Cube3DLocalInterpolation (std::bitset<6> s = 0)
33  {
34  for (size_t i=0; i<6; i++)
35  sign_[i] = (s[i]) ? -1.0 : 1.0;
36 
37  n_[0] = {-1.0, 0.0, 0.0};
38  n_[1] = { 1.0, 0.0, 0.0};
39  n_[2] = { 0.0, -1.0, 0.0};
40  n_[3] = { 0.0, 1.0, 0.0};
41  n_[4] = { 0.0, 0.0, -1.0};
42  n_[5] = { 0.0, 0.0, 1.0};
43  }
44 
53  template<class F, class C>
54  void interpolate (const F& ff, std::vector<C>& out) const
55  {
56  // f gives v*outer normal at a point on the edge!
57  typedef typename LB::Traits::RangeFieldType Scalar;
58  typedef typename LB::Traits::DomainFieldType Vector;
59 
60  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
61 
62  out.resize(36);
63  fill(out.begin(), out.end(), 0.0);
64 
65  const int qOrder = 3;
66  const auto& rule1 = QuadratureRules<Scalar,2>::rule(GeometryTypes::cube(2), qOrder);
67 
68  for (auto&& qp : rule1)
69  {
70  Dune::FieldVector<Scalar,2> qPos = qp.position();
71  typename LB::Traits::DomainType localPos;
72 
73  localPos = {0.0, qPos[0], qPos[1]};
74  auto y = f(localPos);
75  out[0] += (y[0]*n_[0][0] + y[1]*n_[0][1] + y[2]*n_[0][2])*qp.weight()*sign_[0];
76  out[6] += (y[0]*n_[0][0] + y[1]*n_[0][1] + y[2]*n_[0][2])*(2.0*qPos[0] - 1.0)*qp.weight();
77  out[12] += (y[0]*n_[0][0] + y[1]*n_[0][1] + y[2]*n_[0][2])*(2.0*qPos[1] - 1.0)*qp.weight();
78  out[18] += (y[0]*n_[0][0] + y[1]*n_[0][1] + y[2]*n_[0][2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*qp.weight();
79 
80  localPos = {1.0, qPos[0], qPos[1]};
81  y = f(localPos);
82  out[1] += (y[0]*n_[1][0] + y[1]*n_[1][1] + y[2]*n_[1][2])*qp.weight()*sign_[1];
83  out[7] += (y[0]*n_[1][0] + y[1]*n_[1][1] + y[2]*n_[1][2])*(1.0 - 2.0*qPos[0])*qp.weight();
84  out[13] += (y[0]*n_[1][0] + y[1]*n_[1][1] + y[2]*n_[1][2])*(1.0 - 2.0*qPos[1])*qp.weight();
85  out[19] += (y[0]*n_[1][0] + y[1]*n_[1][1] + y[2]*n_[1][2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*qp.weight();
86 
87  localPos = {qPos[0], 0.0, qPos[1]};
88  y = f(localPos);
89  out[2] += (y[0]*n_[2][0] + y[1]*n_[2][1] + y[2]*n_[2][2])*qp.weight()*sign_[2];
90  out[8] += (y[0]*n_[2][0] + y[1]*n_[2][1] + y[2]*n_[2][2])*(1.0 - 2.0*qPos[0])*qp.weight();
91  out[14] += (y[0]*n_[2][0] + y[1]*n_[2][1] + y[2]*n_[2][2])*(2.0*qPos[1] - 1.0)*qp.weight();
92  out[20] += (y[0]*n_[2][0] + y[1]*n_[2][1] + y[2]*n_[2][2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*qp.weight();
93 
94  localPos = {qPos[0], 1.0, qPos[1]};
95  y = f(localPos);
96  out[3] += (y[0]*n_[3][0] + y[1]*n_[3][1] + y[2]*n_[3][2])*qp.weight()*sign_[3];
97  out[9] += (y[0]*n_[3][0] + y[1]*n_[3][1] + y[2]*n_[3][2])*(2.0*qPos[0] - 1.0)*qp.weight();
98  out[15] += (y[0]*n_[3][0] + y[1]*n_[3][1] + y[2]*n_[3][2])*(1.0 - 2.0*qPos[1])*qp.weight();
99  out[21] += (y[0]*n_[3][0] + y[1]*n_[3][1] + y[2]*n_[3][2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*qp.weight();
100 
101  localPos = {qPos[0], qPos[1], 0.0};
102  y = f(localPos);
103  out[4] += (y[0]*n_[4][0] + y[1]*n_[4][1] + y[2]*n_[4][2])*qp.weight()*sign_[4];
104  out[10] += (y[0]*n_[4][0] + y[1]*n_[4][1] + y[2]*n_[4][2])*(1.0 - 2.0*qPos[0])*qp.weight();
105  out[16] += (y[0]*n_[4][0] + y[1]*n_[4][1] + y[2]*n_[4][2])*(1.0 - 2.0*qPos[1])*qp.weight();
106  out[22] += (y[0]*n_[4][0] + y[1]*n_[4][1] + y[2]*n_[4][2])*(1.0 - 2.0*qPos[0])*(2.0*qPos[1] - 1.0)*qp.weight();
107 
108  localPos = {qPos[0], qPos[1], 1.0};
109  y = f(localPos);
110  out[5] += (y[0]*n_[5][0] + y[1]*n_[5][1] + y[2]*n_[5][2])*qp.weight()*sign_[5];
111  out[11] += (y[0]*n_[5][0] + y[1]*n_[5][1] + y[2]*n_[5][2])*(2.0*qPos[0] - 1.0)*qp.weight();
112  out[17] += (y[0]*n_[5][0] + y[1]*n_[5][1] + y[2]*n_[5][2])*(2.0*qPos[1] - 1.0)*qp.weight();
113  out[23] += (y[0]*n_[5][0] + y[1]*n_[5][1] + y[2]*n_[5][2])*(2.0*qPos[0] - 1.0)*(2.0*qPos[1] - 1.0)*qp.weight();
114  }
115 
116  const auto& rule2 = QuadratureRules<Vector,3>::rule(GeometryTypes::cube(3), qOrder);
117  for (auto&& qp : rule2)
118  {
119  FieldVector<double,3> qPos = qp.position();
120 
121  auto y = f(qPos);
122  out[24] += y[0]*qp.weight();
123  out[25] += y[1]*qp.weight();
124  out[26] += y[2]*qp.weight();
125  out[27] += y[0]*qPos[1]*qp.weight();
126  out[28] += y[0]*qPos[2]*qp.weight();
127  out[29] += y[1]*qPos[0]*qp.weight();
128  out[30] += y[1]*qPos[2]*qp.weight();
129  out[31] += y[2]*qPos[0]*qp.weight();
130  out[32] += y[2]*qPos[1]*qp.weight();
131  out[33] += y[0]*qPos[1]*qPos[2]*qp.weight();
132  out[34] += y[1]*qPos[0]*qPos[2]*qp.weight();
133  out[35] += y[2]*qPos[0]*qPos[1]*qp.weight();
134  }
135  }
136 
137  private:
138  // Facet orientations
139  std::array<typename LB::Traits::RangeFieldType, 6> sign_;
140 
141  // Facet normals
142  std::array<typename LB::Traits::DomainType, 6> n_;
143  };
144 }
145 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE3D_LOCALINTERPOLATION_HH
Definition: bdfmcube.hh:16
First order Raviart-Thomas shape functions on the reference hexahedron.
Definition: raviartthomas1cube3dlocalinterpolation.hh:23
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube3dlocalinterpolation.hh:54
RT1Cube3DLocalInterpolation(std::bitset< 6 > s=0)
Make set number s, where 0 <= s < 64.
Definition: raviartthomas1cube3dlocalinterpolation.hh:32