DUNE PDELab (git)

basis.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 © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5
6#ifndef DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH
7#define DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH
8
9#include <cstddef>
10#include <vector>
11
14
15#include <dune/localfunctions/common/localtoglobaladaptors.hh>
16#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
17#include <dune/localfunctions/whitney/edges0.5/common.hh>
18
19namespace Dune {
20
22 //
23 // Basis
24 //
25
27
35 template<class Geometry, class RF>
37 private EdgeS0_5Common<Geometry::mydimension, typename Geometry::ctype>
38 {
39 public:
41 struct Traits {
42 typedef typename Geometry::ctype DomainField;
43 static const std::size_t dimDomainLocal = Geometry::mydimension;
44 static const std::size_t dimDomainGlobal = Geometry::coorddimension;
47
48 typedef RF RangeField;
49 static const std::size_t dimRange = dimDomainLocal;
51
53 };
54
55 private:
56 typedef Dune::Impl::LagrangeSimplexLocalBasis<typename Traits::DomainField,
57 typename Traits::RangeField,
58 Traits::dimDomainLocal,
59 1 // Polynomial order
60 > P1LocalBasis;
62
63 static const P1LocalBasis& p1LocalBasis;
64 static const std::size_t dim = Traits::dimDomainLocal;
65
67 using Base::refelem;
68 using Base::s;
69
70 // global values of the Jacobians (gradients) of the p1 basis
71 std::vector<typename P1Basis::Traits::Jacobian> p1j;
72 // edge sizes and orientations
73 std::vector<typename Traits::DomainField> edgel;
74
75 public:
77
83 template<typename VertexOrder>
84 EdgeS0_5Basis(const Geometry& geo, const VertexOrder& vertexOrder) :
85 p1j(s, typename P1Basis::Traits::Jacobian(0)), edgel(s)
86 {
87 // use some arbitrary position to evaluate jacobians, they are constant
88 static const typename Traits::DomainLocal xl(0);
89
90 // precompute Jacobian (gradients) of the p1 element
91 P1Basis(p1LocalBasis, geo).evaluateJacobian(xl, p1j);
92
93 // calculate edge sizes and orientations
94 for(std::size_t i = 0; i < s; ++i) {
95 edgel[i] = (geo.corner(refelem.subEntity(i,dim-1,0,dim))-
96 geo.corner(refelem.subEntity(i,dim-1,1,dim))
97 ).two_norm();
98 const typename VertexOrder::iterator& edgeVertexOrder =
99 vertexOrder.begin(dim-1, i);
100 if(edgeVertexOrder[0] > edgeVertexOrder[1])
101 edgel[i] *= -1;
102 }
103 }
104
106 std::size_t size () const { return s; }
107
109 void evaluateFunction(const typename Traits::DomainLocal& xl,
110 std::vector<typename Traits::Range>& out) const
111 {
112 out.assign(s, typename Traits::Range(0));
113
114 // compute p1 values -- use the local basis directly for that, local and
115 // global values are identical for scalars
116 std::vector<typename P1LocalBasis::Traits::RangeType> p1v;
117 p1LocalBasis.evaluateFunction(xl, p1v);
118
119 for(std::size_t i = 0; i < s; i++) {
120 const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
121 const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
122 out[i].axpy( p1v[i0], p1j[i1][0]);
123 out[i].axpy(-p1v[i1], p1j[i0][0]);
124 out[i] *= edgel[i];
125 }
126 }
127
130 std::vector<typename Traits::Jacobian>& out) const
131 {
132 out.resize(s);
133
134 for(std::size_t i = 0; i < s; i++) {
135 const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
136 const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
137 for(std::size_t j = 0; j < dim; j++)
138 for(std::size_t k = 0; k < dim; k++)
139 out[i][j][k] = edgel[i] *
140 (p1j[i0][0][k]*p1j[i1][0][j]-p1j[i1][0][k]*p1j[i0][0][j]);
141 }
142 }
143
145 void partial (const std::array<unsigned int, dim>& order,
146 const typename Traits::DomainLocal& in, // position
147 std::vector<typename Traits::Range>& out) const // return value
148 {
149 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
150 if (totalOrder == 0) {
151 evaluateFunction(in, out);
152 } else if (totalOrder==1) {
153 auto const k = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
154 out.resize(size());
155
156 for (std::size_t i = 0; i < s; i++)
157 {
158 const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
159 const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
160 for(std::size_t j = 0; j < dim; j++)
161 out[i][j] = edgel[i] *
162 (p1j[i0][0][k]*p1j[i1][0][j] - p1j[i1][0][k]*p1j[i0][0][j]);
163 }
164 } else {
165 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
166 }
167 }
168
170 std::size_t order () const { return 1; }
171 };
172
173 template<class Geometry, class RF>
174 const typename EdgeS0_5Basis<Geometry, RF>::P1LocalBasis&
175 EdgeS0_5Basis<Geometry, RF>::p1LocalBasis = P1LocalBasis();
176
177} // namespace Dune
178
179#endif // DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH
Basis for order 0.5 (lowest order) edge elements on simplices.
Definition: basis.hh:38
EdgeS0_5Basis(const Geometry &geo, const VertexOrder &vertexOrder)
Construct an EdgeS0_5Basis.
Definition: basis.hh:84
void evaluateJacobian(const typename Traits::DomainLocal &, std::vector< typename Traits::Jacobian > &out) const
Evaluate all Jacobians.
Definition: basis.hh:129
void evaluateFunction(const typename Traits::DomainLocal &xl, std::vector< typename Traits::Range > &out) const
Evaluate all shape functions.
Definition: basis.hh:109
std::size_t size() const
number of shape functions
Definition: basis.hh:106
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainLocal &in, std::vector< typename Traits::Range > &out) const
Evaluate partial derivatives of all shape functions.
Definition: basis.hh:145
std::size_t order() const
Polynomial order of the shape functions.
Definition: basis.hh:170
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:91
Wrapper class for geometries.
Definition: geometry.hh:71
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: geometry.hh:219
static constexpr int coorddimension
dimension of embedding coordinate system
Definition: geometry.hh:97
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:100
static constexpr int mydimension
geometry dimension
Definition: geometry.hh:94
Default exception for dummy implementations.
Definition: exceptions.hh:263
Convert a simple scalar local basis into a global basis.
Definition: localtoglobaladaptors.hh:65
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:279
Dune namespace.
Definition: alignedallocator.hh:13
export type traits for function signature
Definition: basis.hh:41
Common base class for edge elements.
Definition: common.hh:23
RefElem refelem
The reference element for this edge element.
Definition: common.hh:30
std::size_t s
The number of base functions.
Definition: common.hh:38
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 23, 23:29, 2024)