dune-localfunctions  2.7.0
rannacherturek3dlocalbasis.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_RANNACHER_TUREK_3D_LOCALBASIS_HH
4 #define DUNE_RANNACHER_TUREK_3D_LOCALBASIS_HH
5 
6 #include <numeric>
7 #include <vector>
8 
9 #include <dune/common/fvector.hh>
10 #include <dune/common/fmatrix.hh>
11 
13 
14 namespace Dune
15 {
16 
17  template< class D, class R >
19  {
20  static const int coefficients[ 6 ][ 6 ];
21 
22  public:
24  R, 1, FieldVector< R, 1 >,
25  FieldMatrix< R, 1, 3 > > Traits;
26 
28  unsigned int size () const
29  {
30  return 6;
31  }
32 
34  inline void evaluateFunction ( const typename Traits::DomainType &in,
35  std::vector< typename Traits::RangeType > &out ) const
36  {
37  typedef typename Traits::RangeFieldType RangeFieldType;
38  RangeFieldType y[ 6 ] = { 1, in[ 0 ], in[ 1 ], in[ 2 ],
39  in[ 0 ]*in[ 0 ] - in[ 1 ]*in[ 1 ],
40  in[ 1 ]*in[ 1 ] - in[ 2 ]*in[ 2 ] };
41  out.resize( size() );
42  for( unsigned int i = 0; i < size(); ++i )
43  {
44  out[ i ] = RangeFieldType( 0 );
45  for( unsigned int j = 0; j < 6; ++j )
46  out[ i ] += coefficients[ i ][ j ]*y[ j ];
47  out[ i ] /= RangeFieldType( 3 );
48  }
49  }
50 
52  inline void evaluateJacobian ( const typename Traits::DomainType &in,
53  std::vector< typename Traits::JacobianType > &out ) const
54  {
55  typedef typename Traits::RangeFieldType RangeFieldType;
56  RangeFieldType y0[ 5 ] = { 1, 0, 0, 2*in[ 0 ], 0 };
57  RangeFieldType y1[ 5 ] = { 0, 1, 0, -2*in[ 1 ], 2*in[ 1 ] };
58  RangeFieldType y2[ 5 ] = { 0, 0, 1, 0, -2*in[ 2 ] };
59 
60  out.resize( size() );
61  for( unsigned int i = 0; i < size(); ++i )
62  {
63  out[ i ] = RangeFieldType( 0 );
64  for( unsigned int j = 0; j < 5; ++j )
65  {
66  out[ i ][ 0 ][ 0 ] += coefficients[ i ][ j+1 ]*y0[ j ];
67  out[ i ][ 0 ][ 1 ] += coefficients[ i ][ j+1 ]*y1[ j ];
68  out[ i ][ 0 ][ 2 ] += coefficients[ i ][ j+1 ]*y2[ j ];
69  }
70  out[ i ] /= RangeFieldType( 3 );
71  }
72  }
73 
75  void partial (const std::array<unsigned int, 3>& order,
76  const typename Traits::DomainType& in, // position
77  std::vector<typename Traits::RangeType>& out) const // return value
78  {
79  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
80  if (totalOrder == 0) {
81  evaluateFunction(in, out);
82  } else if (totalOrder == 1) {
83  out.resize(size());
84  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
85 
86  using RangeFieldType = typename Traits::RangeFieldType;
87  RangeFieldType y[3][5] = { { 1.0, 0.0, 0.0, 2*in[0], 0.0 },
88  { 0.0, 1.0, 0.0, -2*in[1], 2*in[1] },
89  { 0.0, 0.0, 1.0, 0.0, -2*in[2] } };
90 
91  for (std::size_t i = 0; i < size(); ++i) {
92  out[i] = RangeFieldType{0};
93  for (std::size_t j = 0; j < 5; ++j)
94  out[i] += coefficients[i][j+1] * y[direction][j];
95  out[i] /= RangeFieldType{3};
96  }
97  } else {
98  DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
99  }
100  }
101 
103  unsigned int order () const
104  {
105  return 2;
106  }
107  };
108 
109 
110 
111  // RannacherTurek3DLocalBasis::coefficients
112  // ----------------------------------------
113 
114  template< class D, class R >
115  const int RannacherTurek3DLocalBasis< D, R >
116  ::coefficients[ 6 ][ 6 ] = {{ 2, -7, 2, 2, 4, 2 },
117  { -1, -1, 2, 2, 4, 2 },
118  { 2, 2, -7, 2, -2, 2 },
119  { -1, 2, -1, 2, -2, 2 },
120  { 2, 2, 2, -7, -2, -4 },
121  { -1, 2, 2, -1, -2, -4 }};
122 
123 } //namespace Dune
124 
125 #endif // #ifndef DUNE_RANNACHER_TUREK_3D_LOCALBASIS_HH
Dune::RannacherTurek3DLocalBasis::partial
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: rannacherturek3dlocalbasis.hh:75
Dune::LocalBasisTraits
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:31
Dune::RannacherTurek3DLocalBasis::size
unsigned int size() const
number of shape functions
Definition: rannacherturek3dlocalbasis.hh:28
Dune
Definition: bdfmcube.hh:15
Dune::RannacherTurek3DLocalBasis::evaluateJacobian
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
evaluate jacobian of all shape functions
Definition: rannacherturek3dlocalbasis.hh:52
Dune::LocalBasisTraits::RangeFieldType
RF RangeFieldType
Export type for range field.
Definition: common/localbasis.hh:46
Dune::LocalBasisTraits::DomainType
D DomainType
domain type
Definition: common/localbasis.hh:43
Dune::RannacherTurek3DLocalBasis::Traits
LocalBasisTraits< D, 3, FieldVector< D, 3 >, R, 1, FieldVector< R, 1 >, FieldMatrix< R, 1, 3 > > Traits
Definition: rannacherturek3dlocalbasis.hh:25
Dune::RannacherTurek3DLocalBasis::evaluateFunction
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
evaluate all shape functions
Definition: rannacherturek3dlocalbasis.hh:34
localbasis.hh
Dune::RannacherTurek3DLocalBasis::order
unsigned int order() const
polynomial order of the shape functions
Definition: rannacherturek3dlocalbasis.hh:103
Dune::RannacherTurek3DLocalBasis
Definition: rannacherturek3dlocalbasis.hh:18