Dune Core Modules (2.9.1)
affinegeometry.hh
Go to the documentation of this file.
3// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
53 static void Ax ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &ret )
64 static void ATx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &ret )
75 static void AB ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, n, p > &B, FieldMatrix< ctype, m, p > &ret )
89 static void ATBT ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, p, m > &B, FieldMatrix< ctype, n, p > &ret )
177 static void Lx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
188 static void LTx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
235 static bool cholesky_L ( const FieldMatrix< ctype, n, n > &A, FieldMatrix< ctype, n, n > &ret, const bool checkSingular = false )
344 static bool spdInvAx ( FieldMatrix< ctype, n, n > &A, FieldVector< ctype, n > &x, const bool checkSingular = false )
426 static void leftInvAx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &y )
462 static bool xTRightInvA ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &y )
519 typedef Geo::ReferenceElement< Geo::ReferenceElementImplementation< ctype, mydimension > > ReferenceElement;
532 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
548 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
635 const JacobianTransposed &jacobianTransposed ([[maybe_unused]] const LocalCoordinate &local) const
646 const JacobianInverseTransposed &jacobianInverseTransposed ([[maybe_unused]] const LocalCoordinate &local) const
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:484
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition: affinegeometry.hh:543
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:536
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition: affinegeometry.hh:497
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition: affinegeometry.hh:561
static const int mydimension
Dimension of the geometry.
Definition: affinegeometry.hh:491
AffineGeometry(const ReferenceElement &refElement, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from reference element, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:528
JacobianInverse jacobianInverse(const LocalCoordinate &local) const
Obtain the Jacobian's inverse.
Definition: affinegeometry.hh:668
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition: affinegeometry.hh:553
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: affinegeometry.hh:618
FieldMatrix< ctype, mydimension, coorddimension > JacobianInverse
Type for the inverse Jacobian matrix.
Definition: affinegeometry.hh:515
FieldMatrix< ctype, coorddimension, mydimension > Jacobian
Type for the Jacobian matrix.
Definition: affinegeometry.hh:512
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: affinegeometry.hh:646
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition: affinegeometry.hh:506
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: affinegeometry.hh:567
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: affinegeometry.hh:564
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition: affinegeometry.hh:509
static const int coorddimension
Dimension of the world space.
Definition: affinegeometry.hh:494
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:581
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition: affinegeometry.hh:573
Jacobian jacobian(const LocalCoordinate &local) const
Obtain the Jacobian.
Definition: affinegeometry.hh:657
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition: affinegeometry.hh:500
bool affine() const
Always true: this is an affine geometry.
Definition: affinegeometry.hh:558
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: affinegeometry.hh:635
Volume volume() const
Obtain the volume of the element.
Definition: affinegeometry.hh:624
FieldMatrix< K, COLS, ROWS > transposed() const
Return transposed of the matrix as FieldMatrix.
Definition: fmatrix.hh:172
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelement.hh:52
CoordinateField volume() const
obtain the volume of the reference element
Definition: referenceelement.hh:241
decltype(auto) type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelement.hh:171
int size(int c) const
number of subentities of codimension c
Definition: referenceelement.hh:94
decltype(auto) position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelement.hh:203
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:126
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.
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
unspecified-type ReferenceElement
Returns the type of reference element for the argument type T.
Definition: referenceelements.hh:497
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:170
A unique label for each type of element that can occur in a grid.
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 21, 23:30, 2024)