DUNE PDELab (git)

mappedgeometry.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#ifndef DUNE_GEOMETRY_MAPPEDGEOMETRY_HH
6#define DUNE_GEOMETRY_MAPPEDGEOMETRY_HH
7
8#include <cassert>
9#include <limits>
10#include <optional>
11#include <stdexcept>
12#include <type_traits>
13
14#include <dune/common/copyableoptional.hh>
18#include <dune/common/math.hh>
19#include <dune/common/transpose.hh>
20#include <dune/geometry/affinegeometry.hh> // for FieldMatrixHelper
22#include <dune/geometry/referenceelements.hh>
23#include <dune/geometry/type.hh>
24#include <dune/geometry/utility/algorithms.hh>
25#include <dune/geometry/utility/convergence.hh>
26
27namespace Dune {
28
64template <class Map, class Geo>
66{
67public:
69 using LocalCoordinate = typename Geo::LocalCoordinate;
70
72 using GlobalCoordinate = std::remove_reference_t<decltype(std::declval<Map>()(std::declval<typename Geo::GlobalCoordinate>()))>;
73
75 using ctype = typename Geo::ctype;
76
78 static constexpr int mydimension = LocalCoordinate::size();
79
81 static constexpr int coorddimension = GlobalCoordinate::size();
82
84 using Volume = std::remove_reference_t<decltype(Dune::power(std::declval<ctype>(),mydimension))>;
85
88
91
94
97
98private:
100 using ReferenceElement = typename ReferenceElements::ReferenceElement;
101
102protected:
103 using MatrixHelper = Impl::FieldMatrixHelper<ctype>;
104
105 // type of the mapping representation the geometry parametrization
106 using Mapping = Map;
107
108 // type of the geometry that is wrapped
109 using Geometry = Geo;
110
111 // type of a mapping representing the derivative of `Map` w.r.t. `GlobalCoordinate`
112 using DerivativeMapping = std::remove_reference_t<decltype(derivative(std::declval<Map>()))>;
113
114public:
123 template <class Geo_, class Map_,
124 std::enable_if_t<Dune::IsInteroperable<Map, Map_>::value, int> = 0,
125 std::enable_if_t<Dune::IsInteroperable<Geo, Geo_>::value, int> = 0>
126 MappedGeometry (Map_&& mapping, Geo_&& geometry, bool affine = false)
127 : mapping_(std::forward<Map_>(mapping))
128 , dMapping_(derivative(*mapping_))
129 , geometry_(std::forward<Geo_>(geometry))
130 , affine_(affine)
131 {}
132
137 bool affine () const
138 {
139 return affine_;
140 }
141
144 {
145 return geometry_.type();
146 }
147
149 int corners () const
150 {
151 return geometry_.corners();
152 }
153
156 {
157 assert( (i >= 0) && (i < corners()) );
158 return mapping()(geometry_.corner(i));
159 }
160
163 {
164 return mapping()(geometry_.center());
165 }
166
177 {
178 return mapping()(geometry_.global(local));
179 }
180
197 LocalCoordinate local (const GlobalCoordinate& y, Impl::GaussNewtonOptions<ctype> opts = {}) const
198 {
199 LocalCoordinate x = refElement().position(0,0);
200 Impl::GaussNewtonErrorCode err = Impl::gaussNewton(
201 [&](const LocalCoordinate& local) { return this->global(local); },
202 [&](const LocalCoordinate& local) { return this->jacobianTransposed(local); },
203 y, x, opts
204 );
205
206 if (err != Impl::GaussNewtonErrorCode::OK)
208 "Local coordinate can not be recovered from global coordinate, error code = " << int(err) << "\n"
209 << " (global(x) - y).two_norm() = " << (global(x) - y).two_norm()
210 << " > tol = " << opts.absTol);
211
212 return x;
213 }
214
226 {
227 return MatrixHelper::sqrtDetAAT(jacobianTransposed(local));
228 }
229
238 Volume volume (Impl::ConvergenceOptions<ctype> opts = {}) const
239 {
241 if (affine())
242 return vol0;
243
244 using std::abs;
245 for (int p = 2; p < opts.maxIt; ++p) {
247 if (abs(vol1 - vol0) < opts.absTol)
248 return vol1;
249
250 vol0 = vol1;
251 }
252 return vol0;
253 }
254
257 {
258 Volume vol(0);
259 for (const auto& qp : quadRule)
260 vol += integrationElement(qp.position()) * qp.weight();
261 return vol;
262 }
263
269 Jacobian jacobian (const LocalCoordinate& local) const
270 {
271 auto&& jLocal = geometry_.jacobian(local);
272 auto&& jMapping = (*dMapping_)(geometry_.global(local));
273 return jMapping * jLocal;
274 }
275
282 {
283 return transpose(jacobian(local));
284 }
285
294 {
295 JacobianInverse out;
296 MatrixHelper::leftInvA(jacobian(local), out);
297 return out;
298 }
299
308 {
309 return transpose(jacobianInverse(local));
310 }
311
313 friend ReferenceElement referenceElement (const MappedGeometry& geometry)
314 {
315 return geometry.refElement();
316 }
317
318protected:
319 // the internal stored reference element
320 ReferenceElement refElement () const
321 {
322 return referenceElement(geometry_);
323 }
324
325private:
326 // internal reference to the stored mapping
327 const Mapping& mapping () const
328 {
329 return *mapping_;
330 }
331
332 // internal reference to the stored localgeometry
333 const Geometry& geometry () const
334 {
335 return geometry_;
336 }
337
338private:
340 CopyableOptional<Mapping> mapping_;
341
343 CopyableOptional<DerivativeMapping> dMapping_;
344
346 Geometry geometry_;
347
349 bool affine_;
350};
351
352// deduction guides
353template <class Map, class Geo>
354MappedGeometry (const Map&, const Geo&)
355 -> MappedGeometry<Map,Geo>;
356
357template <class Map, class Geo>
358MappedGeometry (const Map&, const Geo&, bool)
359 -> MappedGeometry<Map,Geo>;
360
361} // end namespace Dune
362
363#endif // DUNE_GEOMETRY_MAPPEDGEOMETRY_HH
An implementation of the Geometry interface for affine geometries.
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Geometry parametrized by a LocalFunction and a LocalGeometry.
Definition: mappedgeometry.hh:66
std::remove_reference_t< decltype(std::declval< Map >()(std::declval< typename Geo::GlobalCoordinate >()))> GlobalCoordinate
type of global coordinates
Definition: mappedgeometry.hh:72
GlobalCoordinate center() const
Map the center of the wrapped geometry.
Definition: mappedgeometry.hh:162
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: mappedgeometry.hh:307
JacobianInverse jacobianInverse(const LocalCoordinate &local) const
Obtain the Jacobian's inverse.
Definition: mappedgeometry.hh:293
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: mappedgeometry.hh:281
MappedGeometry(Map_ &&mapping, Geo_ &&geometry, bool affine=false)
Constructor from mapping to parametrize the geometry.
Definition: mappedgeometry.hh:126
static constexpr int mydimension
geometry dimension
Definition: mappedgeometry.hh:78
GeometryType type() const
Obtain the element type from the reference element.
Definition: mappedgeometry.hh:143
typename Geo::LocalCoordinate LocalCoordinate
type of local coordinates
Definition: mappedgeometry.hh:69
typename Geo::ctype ctype
coordinate type
Definition: mappedgeometry.hh:75
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: mappedgeometry.hh:155
static constexpr int coorddimension
coordinate dimension
Definition: mappedgeometry.hh:81
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the coordinate mapping.
Definition: mappedgeometry.hh:176
Volume volume(const QuadratureRule< ctype, mydimension > &quadRule) const
Obtain the volume of the mapping's image by given quadrature rules.
Definition: mappedgeometry.hh:256
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: mappedgeometry.hh:149
friend ReferenceElement referenceElement(const MappedGeometry &geometry)
Obtain the reference-element related to this geometry.
Definition: mappedgeometry.hh:313
Volume volume(Impl::ConvergenceOptions< ctype > opts={}) const
Obtain the volume of the mapping's image.
Definition: mappedgeometry.hh:238
std::remove_reference_t< decltype(Dune::power(std::declval< ctype >(), mydimension))> Volume
type of volume
Definition: mappedgeometry.hh:84
bool affine() const
Is this mapping affine? Not in general, since we don't know anything about the mapping.
Definition: mappedgeometry.hh:137
Jacobian jacobian(const LocalCoordinate &local) const
Obtain the Jacobian.
Definition: mappedgeometry.hh:269
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: mappedgeometry.hh:225
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:214
A container for all quadrature rules of dimension dim
Definition: quadraturerules.hh:260
A few common exception classes.
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
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:39
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
auto transpose(const Matrix &matrix)
Return the transposed of the given matrix.
Definition: transpose.hh:183
STL namespace.
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:128
typename Container::ReferenceElement ReferenceElement
The reference element type.
Definition: referenceelements.hh:146
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)