Dune Core Modules (2.7.0)

raviartthomas1cube2dlocalinterpolation.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_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALINTERPOLATION_HH
5 
6 #include <vector>
7 
9 #include <dune/localfunctions/common/localinterpolation.hh>
10 
11 
12 namespace Dune
13 {
14 
23  template<class LB>
25  {
26 
27  public:
30  {
31  sign0 = sign1 = sign2 = sign3 = 1.0;
32  }
33 
39  RT1Cube2DLocalInterpolation (unsigned int s)
40  {
41  sign0 = sign1 = sign2 = sign3 = 1.0;
42  if (s & 1)
43  {
44  sign0 = -1.0;
45  }
46  if (s & 2)
47  {
48  sign1 = -1.0;
49  }
50  if (s & 4)
51  {
52  sign2 = -1.0;
53  }
54  if (s & 8)
55  {
56  sign3 = -1.0;
57  }
58 
59  n0[0] = -1.0;
60  n0[1] = 0.0;
61  n1[0] = 1.0;
62  n1[1] = 0.0;
63  n2[0] = 0.0;
64  n2[1] = -1.0;
65  n3[0] = 0.0;
66  n3[1] = 1.0;
67  }
68 
77  template<class F, class C>
78  void interpolate (const F& ff, std::vector<C>& out) const
79  {
80  // f gives v*outer normal at a point on the edge!
81  typedef typename LB::Traits::RangeFieldType Scalar;
82  typedef typename LB::Traits::DomainFieldType Vector;
83 
84  auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
85 
86  out.resize(12);
87  fill(out.begin(), out.end(), 0.0);
88 
89  const int qOrder = 3;
91 
92  for (typename QuadratureRule<Scalar,1>::const_iterator it = rule1.begin();
93  it != rule1.end(); ++it)
94  {
95  Scalar qPos = it->position();
96  typename LB::Traits::DomainType localPos;
97 
98  localPos[0] = 0.0;
99  localPos[1] = qPos;
100  auto y = f(localPos);
101  out[0] += (y[0]*n0[0] + y[1]*n0[1])*it->weight()*sign0;
102  out[1] += (y[0]*n0[0] + y[1]*n0[1])*(2.0*qPos - 1.0)*it->weight();
103 
104  localPos[0] = 1.0;
105  localPos[1] = qPos;
106  y = f(localPos);
107  out[2] += (y[0]*n1[0] + y[1]*n1[1])*it->weight()*sign1;
108  out[3] += (y[0]*n1[0] + y[1]*n1[1])*(1.0 - 2.0*qPos)*it->weight();
109 
110  localPos[0] = qPos;
111  localPos[1] = 0.0;
112  y = f(localPos);
113  out[4] += (y[0]*n2[0] + y[1]*n2[1])*it->weight()*sign2;
114  out[5] += (y[0]*n2[0] + y[1]*n2[1])*(1.0 - 2.0*qPos)*it->weight();
115 
116  localPos[0] = qPos;
117  localPos[1] = 1.0;
118  y = f(localPos);
119  out[6] += (y[0]*n3[0] + y[1]*n3[1])*it->weight()*sign3;
120  out[7] += (y[0]*n3[0] + y[1]*n3[1])*(2.0*qPos - 1.0)*it->weight();
121  }
122 
124 
125  for (typename QuadratureRule<Vector,2>::const_iterator it=rule2.begin(); it!=rule2.end(); ++it)
126  {
127  FieldVector<double,2> qPos = it->position();
128 
129  auto y = f(qPos);
130  out[8] += y[0]*it->weight();
131  out[9] += y[1]*it->weight();
132  out[10] += y[0]*qPos[1]*it->weight();
133  out[11] += y[1]*qPos[0]*it->weight();
134  }
135  }
136 
137  private:
138  typename LB::Traits::RangeFieldType sign0, sign1, sign2, sign3;
139  typename LB::Traits::DomainType n0, n1, n2, n3;
140  };
141 }
142 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_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
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas1cube2dlocalinterpolation.hh:25
void interpolate(const F &ff, std::vector< C > &out) const
Interpolate a given function with shape functions.
Definition: raviartthomas1cube2dlocalinterpolation.hh:78
RT1Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas1cube2dlocalinterpolation.hh:29
RT1Cube2DLocalInterpolation(unsigned int s)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalinterpolation.hh:39
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)