Rodin::Geometry::PointBase class

Base class for spatial points on a discrete mesh.

This class represents a point in the computational domain, maintaining both reference and physical coordinates along with the associated polytope.

Mathematical Foundation

A Point represents the tuple $ (x, r, p, \tau) $ where:

  • $ \tau $ is a polytope in the mesh $ \mathcal{T}_h $
  • $ x : K \rightarrow \tau $ is the geometric transformation
  • $ r \in K $ are the reference coordinates
  • $ p = x(r) \in \tau $ are the physical coordinates

The transformation satisfies:

\[ p = x(r) \]

Coordinate Systems

Points maintain two coordinate systems:

  • Reference coordinates $ r $ : coordinates in the reference element $ K $
  • Physical coordinates $ p $ : coordinates in the physical element $ \tau $

Jacobian Information

The class provides access to:

  • Jacobian matrix $ \mathbf{J}_x(r) = \frac{\partial x}{\partial r} $
  • Jacobian determinant $ |\mathbf{J}_x(r)| $
  • Inverse Jacobian $ \mathbf{J}_x^{-1}(r) $
  • Distortion measure

Thread Safety

This class is not thread-safe. Each thread should use its own Point instances.

Derived classes

class Point final
Concrete implementation of a spatial point on a mesh.

Public types

enum class Coordinates { Reference, Physical }
Enumeration of coordinate types.

Constructors, destructors, conversion operators

PointBase(const Polytope& polytope) explicit
Constructs a point on a polytope.
PointBase(const Polytope& polytope, const Math::SpatialPoint& pc) explicit
Constructs a point with physical coordinates.
PointBase(const Polytope& polytope, Math::SpatialPoint&& pc) explicit
Constructs a point with physical coordinates (move).
PointBase(const PointBase& other)
Copy constructor.
PointBase(PointBase&& other)
Move constructor.

Public functions

auto getDimension(Coordinates coords = Coordinates::Physical) const -> size_t
Gets the dimension of coordinates.
auto operator()(size_t i) const -> Real
Accesses the i-th physical coordinate.
auto asVector() const -> const auto&
Gets physical coordinates as a vector.
auto x() const -> Real
Gets the x-coordinate (first component).
auto y() const -> Real
Gets the y-coordinate (second component).
auto z() const -> Real
Gets the z-coordinate (third component).
auto norm() const -> Real
Computes the Euclidean norm.
auto stableNorm() const -> Real
Computes the stable norm (avoids overflow).
auto blueNorm() const -> Real
Computes the Blue norm.
template<size_t p>
auto lpNorm() const -> Real
Computes the Lp norm.
auto squaredNorm() const -> Real
Computes the squared Euclidean norm.
auto operator<(const PointBase& p) const -> bool
Lexicographical comparison operator.
auto getPolytope() const -> const Polytope&
Gets the associated polytope.
auto getPhysicalCoordinates() const -> const Math::SpatialVector<Real>&
Gets the physical coordinates.
auto getCoordinates(Coordinates coords = Coordinates::Physical) const -> const Math::SpatialVector<Real>&
Gets coordinates of specified type.
auto getJacobian() const -> const Math::SpatialMatrix<Real>& virtual
Computes the Jacobian matrix at this point.
auto getJacobianDeterminant() const -> Real
Gets the Jacobian determinant at this point.
auto getJacobianInverse() const -> const Math::SpatialMatrix<Real>&
Computes the inverse Jacobian matrix.
auto getDistortion() const -> Real
Computes the distortion measure at this point.
auto setPolytope(const Polytope& polytope) -> PointBase&
Sets the polytope for this point.
auto getReferenceCoordinates() const -> const Math::SpatialVector<Real>& pure virtual
Gets the reference coordinates (pure virtual).

Enum documentation

enum class Rodin::Geometry::PointBase::Coordinates

Enumeration of coordinate types.

Enumerators
Reference

Reference coordinates $ r \in K $ .

Physical

Physical coordinates $ p \in \tau $ .

Function documentation

Rodin::Geometry::PointBase::PointBase(const Polytope& polytope) explicit

Constructs a point on a polytope.

Parameters
polytope in The polytope containing this point

Rodin::Geometry::PointBase::PointBase(const Polytope& polytope, const Math::SpatialPoint& pc) explicit

Constructs a point with physical coordinates.

Parameters
polytope in The polytope containing this point
pc in Physical coordinates

Rodin::Geometry::PointBase::PointBase(const Polytope& polytope, Math::SpatialPoint&& pc) explicit

Constructs a point with physical coordinates (move).

Parameters
polytope in The polytope containing this point
pc in Physical coordinates

size_t Rodin::Geometry::PointBase::getDimension(Coordinates coords = Coordinates::Physical) const

Gets the dimension of coordinates.

Parameters
coords in Type of coordinates (Reference or Physical)
Returns Dimension of the specified coordinate system

Real Rodin::Geometry::PointBase::operator()(size_t i) const

Accesses the i-th physical coordinate.

Parameters
in Coordinate index
Returns i-th physical coordinate value

const auto& Rodin::Geometry::PointBase::asVector() const

Gets physical coordinates as a vector.

Returns Reference to physical coordinate vector

Real Rodin::Geometry::PointBase::x() const

Gets the x-coordinate (first component).

Returns x-coordinate value

Real Rodin::Geometry::PointBase::y() const

Gets the y-coordinate (second component).

Returns y-coordinate value

Real Rodin::Geometry::PointBase::z() const

Gets the z-coordinate (third component).

Returns z-coordinate value

Real Rodin::Geometry::PointBase::norm() const

Computes the Euclidean norm.

Returns Euclidean norm $ \|p\|_2 $

Real Rodin::Geometry::PointBase::stableNorm() const

Computes the stable norm (avoids overflow).

Returns Stable norm

Real Rodin::Geometry::PointBase::blueNorm() const

Computes the Blue norm.

Returns Blue norm

template<size_t p>
Real Rodin::Geometry::PointBase::lpNorm() const

Computes the Lp norm.

Template parameters
p Norm order
Returns Lp norm $ \|p\|_p $

Real Rodin::Geometry::PointBase::squaredNorm() const

Computes the squared Euclidean norm.

Returns Squared norm $ \|p\|_2^2 $

bool Rodin::Geometry::PointBase::operator<(const PointBase& p) const

Lexicographical comparison operator.

Parameters
in Point to compare with
Returns True if this point is lexicographically less than p

const Polytope& Rodin::Geometry::PointBase::getPolytope() const

Gets the associated polytope.

Returns Reference to the polytope containing this point

const Math::SpatialVector<Real>& Rodin::Geometry::PointBase::getPhysicalCoordinates() const

Gets the physical coordinates.

Returns Reference to physical coordinate vector

const Math::SpatialVector<Real>& Rodin::Geometry::PointBase::getCoordinates(Coordinates coords = Coordinates::Physical) const

Gets coordinates of specified type.

Parameters
coords in Type of coordinates to retrieve
Returns Reference to coordinate vector

const Math::SpatialMatrix<Real>& Rodin::Geometry::PointBase::getJacobian() const virtual

Computes the Jacobian matrix at this point.

Returns Jacobian matrix $ \mathbf{J}_x(r) $

The Jacobian is computed and cached on first access.

Real Rodin::Geometry::PointBase::getJacobianDeterminant() const

Gets the Jacobian determinant at this point.

Returns Determinant $ |\mathbf{J}_x(r)| $

The determinant is cached on first access.

const Math::SpatialMatrix<Real>& Rodin::Geometry::PointBase::getJacobianInverse() const

Computes the inverse Jacobian matrix.

Returns Inverse Jacobian $ \mathbf{J}_x^{-1}(r) $

The inverse is computed and cached on first access.

Real Rodin::Geometry::PointBase::getDistortion() const

Computes the distortion measure at this point.

Returns Distortion value

Measures how much the transformation distorts space at this point.

PointBase& Rodin::Geometry::PointBase::setPolytope(const Polytope& polytope)

Sets the polytope for this point.

Parameters
polytope in New polytope
Returns Reference to this point for method chaining

const Math::SpatialVector<Real>& Rodin::Geometry::PointBase::getReferenceCoordinates() const pure virtual

Gets the reference coordinates (pure virtual).

Returns Reference coordinate vector

Must be implemented by derived classes.