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