KASKADE 7 development version
barycentric.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2020 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef BARYCENTRIC_HH
14#define BARYCENTRIC_HH
15
16#include <cassert>
17#include <numeric> // std::accumulate, with clang++
18
19#include "dune/common/fvector.hh"
20#include "dune/common/fmatrix.hh"
21
22namespace Kaskade
23{
31 template <class CoordType, int dim>
33 {
35 zeta[dim] = 1;
36 for (int i=0; i<dim; ++i)
37 {
38 zeta[i] = x[i];
39 zeta[dim] -= x[i];
40 }
41
42 // Shape functions vanishing on the boundary will probably rely on
43 // one barycentric coordinate being exactly zero, in order not to be
44 // affected by large penalty values from Dirichlet boundary
45 // conditions. Hence we here enforce exactly zero values for very
46 // small barcyentric coordinates.
47 constexpr CoordType eps = std::numeric_limits<CoordType>::epsilon();
48 for (int i=0; i<=dim; ++i)
49 if (std::abs(zeta[i]) < 100*eps) // well, that's probably *exactly* on the boundary
50 zeta[i] = 0;
51
52 return zeta;
53 }
54
55 // ----------------------------------------------------------------------------------------------
56
61 template <class CoordType, int dim,
62 std::enable_if_t<std::is_floating_point_v<CoordType>,int> = 0>
64 {
65 constexpr CoordType eps = std::numeric_limits<CoordType>::epsilon();
66
67 auto b = barycentric(x); // This is already modified such that b>=0.
68 b /= (b.one_norm()+eps); // Now sum(b) = 1
69
70 std::copy(b.begin(),b.begin()+dim,x.begin());
71 return x;
72 }
73
74 // ----------------------------------------------------------------------------------------------
75
89 template <class CoordType, int dim>
91 {
93 zeta[0] = 1;
94 for (int i=0; i<dim; ++i) {
95 zeta[i+1] = x[i];
96 zeta[0] -= x[i];
97 }
98
99 // Shape functions vanishing on the boundary will probably rely on
100 // one barycentric coordinate being exactly zero, in order not to be
101 // affected by large penalty values from Dirichlet boundary
102 // conditions. Hence we here enforce exactly zero values for very
103 // small barcyentric coordinates.
104 constexpr CoordType eps = std::numeric_limits<CoordType>::epsilon();
105 for (int i=0; i<=dim; ++i)
106 if (std::abs(zeta[i]) < 100*eps) // well, that's probably *exactly* on the boundary
107 zeta[i] = 0;
108
109 return zeta;
110 }
111
121 template <class CoordType, int dim>
123 {
125 for (int i=0; i<dim; ++i)
126 {
127 B[0][i] = -1.0;
128 B[i+1][i] = 1.0;
129 }
130 return B;
131 }
132
133 // ----------------------------------------------------------------------------------------------
134
135
147 template <size_t dim>
148 std::array<int,dim+1> barycentric(std::array<int,dim> const& x, int bsum)
149 {
150 std::array<int,dim+1> zeta;
151 std::copy(begin(x),end(x),begin(zeta));
152 zeta[dim] = bsum - std::accumulate(begin(x),end(x),0);
153 return zeta;
154 }
155
156} /* end of namespace Kaskade */
157
158#endif
Dune::FieldVector< CoordType, dim+1 > barycentric(Dune::FieldVector< CoordType, dim > const &x)
Computes the barycentric coordinates of a point in the unit simplex.
Definition: barycentric.hh:32
Dune::FieldMatrix< CoordType, dim+1, dim > barycentric2Derivative(Dune::FieldVector< CoordType, dim > const &x)
Computes the derivative of the barycentric coordinates of a point in the unit simplex.
Definition: barycentric.hh:122
Dune::FieldVector< CoordType, dim+1 > barycentric2(Dune::FieldVector< CoordType, dim > const &x)
Computes barycentric coordinates of a point in the unit simplex.
Definition: barycentric.hh:90
Dune::FieldVector< CoordType, dim > projectToUnitSimplex(Dune::FieldVector< CoordType, dim > x)
Projects a vector onto the unit simplex.
Definition: barycentric.hh:63