dune-localfunctions  2.9.0
rannacherturek2dlocalbasis.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_2D_LOCALBASIS_HH
6 #define DUNE_RANNACHER_TUREK_2D_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  {
23  R, 1, FieldVector< R, 1 >,
24  FieldMatrix< R, 1, 2 > > Traits;
25 
27  unsigned int size () const
28  {
29  return 4;
30  }
31 
33  inline void evaluateFunction ( const typename Traits::DomainType &in,
34  std::vector< typename Traits::RangeType > &out ) const
35  {
36  out.resize(4);
37  typename Traits::DomainFieldType qbase = in[0]*in[0]-in[1]*in[1];
38  out[0] = .75 - 2*in[0] + in[1] + qbase;
39  out[1] = -.25 + in[1] + qbase;
40  out[2] = .75 + in[0] - 2*in[1] - qbase;
41  out[3] = -.25 + in[0] - qbase;
42  }
43 
45  inline void evaluateJacobian ( const typename Traits::DomainType &in,
46  std::vector< typename Traits::JacobianType > &out ) const
47  {
48  out.resize(4);
49 
50  // see http://www.dune-project.org/doc/doxygen/html/classDune_1_1C1LocalBasisInterface.html#d6f8368f8aa43439cc7ef10419f6e2ea
51  // out[i][j][k] = d_k \phi^i_j , where \phi^i_j is the j'th component of the i'th shape function.
52 
53  out[0][0][0] = -2 + 2*in[0]; out[0][0][1] = 1 - 2*in[1];
54  out[1][0][0] = 2*in[0]; out[1][0][1] = 1 - 2*in[1];
55  out[2][0][0] = 1 - 2*in[0]; out[2][0][1] = -2 + 2*in[1];
56  out[3][0][0] = 1 - 2*in[0]; out[3][0][1] = 2*in[1];
57  }
58 
60  void partial (const std::array<unsigned int, 2>& order,
61  const typename Traits::DomainType& in, // position
62  std::vector<typename Traits::RangeType>& out) const // return value
63  {
64  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
65  if (totalOrder == 0) {
66  evaluateFunction(in, out);
67  } else if (totalOrder == 1) {
68  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
69  out.resize(size());
70 
71  switch (direction) {
72  case 0:
73  out[0] = -2 + 2*in[0];
74  out[1] = 2*in[0];
75  out[2] = 1 - 2*in[0];
76  out[3] = 1 - 2*in[0];
77  break;
78  case 1:
79  out[0] = 1 - 2*in[1];
80  out[1] = 1 - 2*in[1];
81  out[2] = -2 + 2*in[1];
82  out[3] = 2*in[1];
83  break;
84  default:
85  DUNE_THROW(RangeError, "Component out of range.");
86  }
87  } else if (totalOrder == 2) {
88  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 2));
89  out.resize(size());
90 
91  switch (direction) {
92  case 0:
93  out[0] = out[1] = 2;
94  out[2] = out[3] =-2;
95  break;
96  case 1:
97  out[0] = out[1] =-2;
98  out[2] = out[3] = 2;
99  break;
100  default:
101  out[0] = out[1] = out[2] = out[3] = 0;
102  break;
103  }
104  } else {
105  out[0] = out[1] = out[2] = out[3] = 0;
106  }
107  }
108 
110  unsigned int order () const
111  {
112  // must be 2 here since it contains x^2 and x^2
113  return 2;
114  }
115  };
116 
117 } //namespace Dune
118 
119 #endif // #ifndef DUNE_RANNACHER_TUREK_2D_LOCALBASIS_HH
Definition: bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:34
D DomainType
domain type
Definition: common/localbasis.hh:42
DF DomainFieldType
Export type for domain field.
Definition: common/localbasis.hh:36
Definition: rannacherturek2dlocalbasis.hh:21
LocalBasisTraits< D, 2, FieldVector< D, 2 >, R, 1, FieldVector< R, 1 >, FieldMatrix< R, 1, 2 > > Traits
Definition: rannacherturek2dlocalbasis.hh:24
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: rannacherturek2dlocalbasis.hh:60
unsigned int size() const
number of shape functions
Definition: rannacherturek2dlocalbasis.hh:27
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
evaluate all shape functions
Definition: rannacherturek2dlocalbasis.hh:33
unsigned int order() const
polynomial order of the shape functions
Definition: rannacherturek2dlocalbasis.hh:110
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
evaluate jacobian of all shape functions
Definition: rannacherturek2dlocalbasis.hh:45