Dune Core Modules (2.9.0)

dualq1localbasis.hh
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_DUAL_Q1_LOCALBASIS_HH
6#define DUNE_DUAL_Q1_LOCALBASIS_HH
7
8#include <array>
9#include <numeric>
10
13
14#include <dune/localfunctions/common/localbasis.hh>
15
16namespace Dune
17{
27 template<class D, class R, int dim>
29 {
30 public:
33
34 void setCoefficients(const std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)>& coefficients)
35 {
36 coefficients_ = coefficients;
37 }
38
40 unsigned int size () const
41 {
42 return 1<<dim;
43 }
44
46 inline void evaluateFunction (const typename Traits::DomainType& in,
47 std::vector<typename Traits::RangeType>& out) const
48 {
49 // compute q1 values
50 std::vector<typename Traits::RangeType> q1Values(size());
51
52 for (size_t i=0; i<size(); i++) {
53
54 q1Values[i] = 1;
55
56 for (int j=0; j<dim; j++)
57 // if j-th bit of i is set multiply with in[j], else with 1-in[j]
58 q1Values[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
59
60 }
61
62 // compute the dual values by using that they are linear combinations of q1 functions
63 out.resize(size());
64 for (size_t i=0; i<size(); i++)
65 out[i] = 0;
66
67 for (size_t i=0; i<size(); i++)
68 for (size_t j=0; j<size(); j++)
69 out[i] += coefficients_[i][j]*q1Values[j];
70
71
72 }
73
75 inline void
76 evaluateJacobian (const typename Traits::DomainType& in, // position
77 std::vector<typename Traits::JacobianType>& out) const // return value
78 {
79 // compute q1 jacobians
80 std::vector<typename Traits::JacobianType> q1Jacs(size());
81
82 // Loop over all shape functions
83 for (size_t i=0; i<size(); i++) {
84
85 // Loop over all coordinate directions
86 for (int j=0; j<dim; j++) {
87
88 // Initialize: the overall expression is a product
89 // if j-th bit of i is set to -1, else 1
90 q1Jacs[i][0][j] = (i & (1<<j)) ? 1 : -1;
91
92 for (int k=0; k<dim; k++) {
93
94 if (j!=k)
95 // if k-th bit of i is set multiply with in[j], else with 1-in[j]
96 q1Jacs[i][0][j] *= (i & (1<<k)) ? in[k] : 1-in[k];
97
98 }
99
100 }
101
102 }
103
104 // compute the dual jacobians by using that they are linear combinations of q1 functions
105 out.resize(size());
106 for (size_t i=0; i<size(); i++)
107 out[i] = 0;
108
109 for (size_t i=0; i<size(); i++)
110 for (size_t j=0; j<size(); j++)
111 out[i].axpy(coefficients_[i][j],q1Jacs[j]);
112
113 }
114
116 void partial (const std::array<unsigned int, dim>& order,
117 const typename Traits::DomainType& in, // position
118 std::vector<typename Traits::RangeType>& out) const // return value
119 {
120 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
121 if (totalOrder == 0) {
122 evaluateFunction(in, out);
123 } else {
124 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
125 }
126 }
127
129 unsigned int order () const
130 {
131 return 1;
132 }
133
134 private:
135 std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)> coefficients_;
136 };
137}
138#endif
Dual Lagrange shape functions of order 1 on the reference cube.
Definition: dualq1localbasis.hh:29
unsigned int size() const
number of shape functions
Definition: dualq1localbasis.hh:40
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualq1localbasis.hh:129
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualq1localbasis.hh:46
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualq1localbasis.hh:76
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: dualq1localbasis.hh:116
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Default exception for dummy implementations.
Definition: exceptions.hh:263
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Dune namespace.
Definition: alignedallocator.hh:13
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:42
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)