coordinax.representations

Contents

coordinax.representations#

The coordinax.representations module defines representation descriptors for geometric data and the representation-aware conversion API.

Overview#

A representation is the triple \(R = (K, B, S)\):

  • \(K\): geometry kind (what sort of geometric object the data represents)

  • \(B\): basis kind (how components are interpreted as basis components)

  • \(S\): semantic kind (what the object means physically/mathematically)

This is separate from charts and manifolds:

  • Charts define coordinate systems.

  • Manifolds define the underlying space.

  • Representations define geometric meaning and transformation law category.

Quick Start#

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> import unxt as u

>>> # Canonical point representation.
>>> rep = cxr.point

>>> # Convert point data between charts while preserving representation.
>>> q_cart = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> q_sph = cxr.cconvert(q_cart, cxc.cart3d, rep, cxc.sph3d, rep)
>>> q_sph
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

>>> # Build a reusable conversion map.
>>> to_sph = cxr.cmap(cxc.cart3d, cxr.point, cxc.sph3d)
>>> to_sph(q_cart)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

# Change tangent components between basis conventions in the same chart.
>>> v = {"r": u.Q(1.0, "km/s"), "theta": u.Q(0.0, "rad/s"), "phi": u.Q(0.0, "rad/s")}
>>> at = {"r": u.Q(2.0, "km"), "theta": u.Q(3.0, "rad"), "phi": u.Q(4.0, "rad")}
>>> v2 = cxr.change_basis(v, cxc.sph3d, cxr.coord_basis, cxr.phys_basis, at=at)
>>> v2
{'r': Q(1., 'km / s'), 'theta': Q(0., 'km / s'), 'phi': Q(0., 'km / s')}

Functional API#

  • cconvert: representation-aware coordinate conversion API

  • change_basis: same-chart tangent basis conversion API

  • cmap: partial-function builder around cconvert

  • guess_basis_kind: infer basis kind from dimensions/data

  • guess_geometry_kind: infer geometric kind from dimensions/data

  • guess_rep: infer full representation from dimensions/data

  • guess_semantic_kind: infer semantic kind from dimensions/data

For point data, cconvert dispatches through chart-level point conversion laws.

Available Objects#

Conversion Functions#

  • cconvert: convert data across charts/representations

  • change_basis: change tangent basis without changing chart

  • cmap: build reusable conversion callables

  • guess_basis_kind: infer a basis kind

  • guess_geometry_kind: infer a geometry kind

  • guess_rep: infer a full representation

  • guess_semantic_kind: infer a semantic kind

Representation Descriptor#

  • Representation: immutable descriptor (geom_kind, basis, semantic_kind)

  • point: canonical point representation

point is equivalent to (point_geom, no_basis, loc).

Geometry Kind#

  • AbstractGeometry: base class for geometric kind descriptors

  • PointGeometry / point_geom: point geometric kind

Basis Kind#

  • AbstractBasis: base class for basis descriptors

  • NoBasis / no_basis: basis kind used for affine point data

  • CoordinateBasis / coord_basis: coordinate-basis tangent components

  • PhysicalBasis / phys_basis: physical-basis tangent components

Semantic Kind#

  • AbstractSemanticKind: base class for semantic descriptors

  • Location / loc: location semantic kind

Notes#

  • Representations are orthogonal to charts: chart choice and representation choice are independent concerns.

  • The current built-in flow is point-first, centered on (point_geom, no_basis, loc).

  • change_basis currently supports tangent-data conversions between CoordinateBasis \(\rightleftarrows\) PhysicalBasis, including Cartesian, spherical, and metric-based conversions.

Representation descriptors for geometric coordinate data.

This subpackage defines the representation layer of coordinax: static, immutable objects that describe what coordinate or component data means geometrically.

In coordinax, a representation is separate from a chart.

  • A chart describes how a manifold is coordinatized locally.

  • A representation describes what sort of geometric object the data represents, whether it is basis-dependent, and what semantic interpretation it carries.

Concretely, a representation is an ordered triple $$ R = (K, B, S), $$ where

  • $K$ is the geometric kind,

  • $B$ is the basis kind, and

  • $S$ is the semantic kind.

For a point, the canonical representation is $$ (mathrm{PointGeometry},, mathrm{NoBasis},, mathrm{Location}), $$

which indicates that the data represents a point on a manifold, that it is not expressed in a basis-dependent linear space, and that its semantic meaning is a location.

This separation is important because the same chart may be used to express different geometric objects, while the same representation may be expressed in many different charts. The role of this subpackage is therefore to provide the metadata needed to interpret geometric data correctly and to dispatch representation-aware conversions such as cconvert.

The main public objects provided here are:

  • Representation, the full representation descriptor,

  • AbstractGeometry and concrete geometric kinds such as PointGeometry,

  • AbstractBasis and concrete basis kinds such as NoBasis,

  • AbstractSemanticKind and concrete semantic kinds such as Location, and

  • cconvert, the central conversion function that converts data between charts and, more generally, between compatible representations.

Now let’s work through some examples.

>>> import coordinax.representations as cxr

Construct the canonical point representation explicitly:

>>> rep = cxr.Representation(cxr.PointGeometry(), cxr.NoBasis(), cxr.Location())
>>> rep
point

This is also provided as a convenient pre-defined instance.

>>> rep == cxr.point
True

Representations are used together with charts to convert geometric data while preserving its interpretation:

>>> import coordinax.charts as cxc
>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> q = cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d, cxr.point)
>>> q
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

Here p is interpreted as point data in Cartesian 3D coordinates, and q is the same point expressed in spherical coordinates.

The important distinction is that the chart changes, while the representation remains the same:

>>> from_rep = cxr.point
>>> to_rep = cxr.point
>>> q = cxr.cconvert(p, cxc.cart3d, from_rep, cxc.sph3d, to_rep)

This reflects the current design focus of coordinax.representations: point data first, with the representation structure already in place for later extension to tangent, cotangent, and other geometric objects.

coordinax.representations.add(*args: Any, **kwargs: Any)#

Add two coordinate data objects.

This is an abstract API definition. See the main coordinax package for concrete implementations.

coordinax.representations.add(lhs: Any, lhs_chart: AbstractChart, lhs_rep: Representation, rhs: Any, rhs_chart: AbstractChart, rhs_rep: Representation, /) Any
Parameters:
Return type:

Any

Add two coordinate data objects via Cartesian round-trip.

Both operands are converted to the ambient Cartesian chart of lhs_chart, added component-wise, then converted back to lhs_chart. If lhs_chart is already Cartesian (or has no global Cartesian), rhs is converted into lhs_chart and added directly.

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> import unxt as u
>>> p1 = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> p2 = {"x": u.Q(4, "m"), "y": u.Q(5, "m"), "z": u.Q(6, "m")}
>>> cxr.add(p1, cxc.cart3d, cxr.point, p2, cxc.cart3d, cxr.point)
{'x': Q(5, 'm'), 'y': Q(7, 'm'), 'z': Q(9, 'm')}
coordinax.representations.add(lhs: coordinax.vectors._src.point.Point, rhs: coordinax.vectors._src.point.Point, /) coordinax.vectors._src.point.Point
Parameters:
Return type:

Any

Add two points.

For non-Cartesian charts the operation converts both operands to the ambient Cartesian chart, adds there, and converts the result back to the lhs chart. For Cartesian charts the addition is direct.

The result keeps the lhs chart and representation.

>>> import coordinax.main as cx
>>> v1 = cx.Point.from_([1, 2, 3], "m")
>>> v2 = cx.Point.from_([4, 5, 6], "m")
>>> print(cxr.add(v1, v2))
<Point: chart=Cart3D (x, y, z) [m]
    [5 7 9]>
coordinax.representations.add(lhs: coordinax.vectors._src.tangent.Tangent, rhs: coordinax.vectors._src.tangent.Tangent, /) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Add two tangent vectors component-wise.

Tangent spaces are genuine vector spaces: addition is component-wise in any chart basis (no Cartesian round-trip is needed or correct). Both operands must share the same representation (chart + basis + semantic).

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v1 = cx.Tangent.from_(
...     {"x": u.Q(1.0, "m/s"), "y": u.Q(2.0, "m/s"), "z": u.Q(3.0, "m/s")},
...     cxc.cart3d, cxr.coord_vel,
... )
>>> v2 = cx.Tangent.from_(
...     {"x": u.Q(4.0, "m/s"), "y": u.Q(5.0, "m/s"), "z": u.Q(6.0, "m/s")},
...     cxc.cart3d, cxr.coord_vel,
... )
>>> result = cxr.add(v1, v2)
>>> result["x"]
Q(5., 'm / s')
Parameters:
Return type:

Any

coordinax.representations.cconvert(*args: Any, **kwargs: Any)#

Transform the current vector to the target chart.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc

Define a point in Cartesian coordinates:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}

Convert it to spherical coordinates:

>>> cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d, cxr.point)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}
coordinax.representations.cconvert(obj: NoneType, /, *fixed_args: Any, **fixed_kw: Any) Any
Parameters:
Return type:

Any

Return a partial function for vector conversion.

Convert a point from Cartesian coordinates to spherical coordinates:

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc

Define a point in Cartesian coordinates:

>>> q = {"x": 1.0, "y": 2.0, "z": 3.0}

Convert it to spherical coordinates:

>>> map = cxr.cconvert(None, cxc.cart3d, cxr.point, cxc.sph3d)
>>> map(q)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}
coordinax.representations.cconvert(x: Any, from_chart: AbstractChart, from_rep: Representation, to_chart: AbstractChart, to_rep: Representation, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) Any
Parameters:
Return type:

Any

Convert point data between charts.

Convert a point from Cartesian coordinates to spherical coordinates:

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc

Define a point in Cartesian coordinates:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}

Convert it to spherical coordinates:

>>> q = cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d, cxr.point)
>>> q
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output q represents the same geometric point but expressed in the target chart.

The representation remains unchanged; only the chart changes:

>>> cxr.cconvert(q, cxc.sph3d, cxr.point, cxc.cart3d, cxr.point)
{'x': Array(1., dtype=float64), 'y': Array(2., dtype=float64),
 'z': Array(3., dtype=float64)}

Let’s work through more examples.

Cartesian to Spherical (with units):

>>> import unxt as u
>>> p = {"x": u.Q(1.0, "m"), "y": u.Q(0.0, "m"), "z": u.Q(0.0, "m")}
>>> cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d, cxr.point)
{'r': Q(1., 'm'), 'theta': Q(1.57079633, 'rad'), 'phi': Q(0., 'rad')}

Cylindrical to Cartesian (without units):

>>> p = {"rho": 3.0, "phi": 0, "z": 4.0}
>>> cxr.cconvert(p, cxc.cyl3d, cxr.point, cxc.cart3d, cxr.point)
{'x': Array(3., dtype=float64, ...), 'y': Array(0., dtype=float64, ...),
 'z': 4.0}

Polar to Cartesian (2D):

>>> p = {"r": u.Q(5.0, "m"), "theta": u.Q(90, "deg")}
>>> cxr.cconvert(p, cxc.polar2d, cxr.point, cxc.cart2d, cxr.point)
{'x': Q(3.061617e-16, 'm'), 'y': Q(5., 'm')}

Between Spherical variants (Spherical to LonLatSpherical):

>>> p = {"r": u.Q(1.0, "m"), "theta": u.Q(45, "deg"), "phi": u.Q(0, "deg")}
>>> cxr.cconvert(p, cxc.sph3d, cxr.point, cxc.lonlat_sph3d, cxr.point)
{'lon': Q(0, 'deg'), 'lat': Q(45., 'deg'), 'distance': Q(1., 'm')}

Identity conversion (same chart):

>>> p = {"x": u.Q(2.0, "m"), "y": u.Q(3.0, "m")}
>>> cxr.cconvert(p, cxc.cart2d, cxr.point, cxc.cart2d, cxr.point) is p
True
coordinax.representations.cconvert(x: Any, from_chart: AbstractChart, from_rep: Representation, to_chart: AbstractChart, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) Any
Parameters:
Return type:

Any

Convert point data between charts.

Convert a point from Cartesian coordinates to spherical coordinates:

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc

Define a point in Cartesian coordinates:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}

Convert it to spherical coordinates:

>>> q = cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d)
>>> q
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output q represents the same geometric point but expressed in the target chart.

The representation remains unchanged; only the chart changes:

>>> cxr.cconvert(q, cxc.sph3d, cxr.point, cxc.cart3d)
{'x': Array(1., dtype=float64), 'y': Array(2., dtype=float64),
 'z': Array(3., dtype=float64)}
coordinax.representations.cconvert(x: Any, from_chart: AbstractChart, from_geom: PointGeometry, from_rep: Representation, to_chart: AbstractChart, to_geom: PointGeometry, to_rep: Representation, /, *, usys: AbstractUnitSystem | None = None) Any
Parameters:
Return type:

Any

Convert point data between charts.

This function delegates to coordinax.charts.pt_map. The representation arguments are checked to ensure they correspond to canonical point data:

$$(mathrm{PointGeometry},, mathrm{NoBasis},, mathrm{Location}).$$

Convert a point from Cartesian coordinates to spherical coordinates:

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc

Define a point in Cartesian coordinates:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}

Convert it to spherical coordinates:

>>> cxr.cconvert(p, cxc.cart3d, cxr.point_geom, cxr.point,
...                 cxc.sph3d, cxr.point_geom, cxr.point)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}
coordinax.representations.cconvert(x: Any, from_chart: AbstractChart, from_geom: TangentGeometry, from_rep: Representation, to_chart: AbstractChart, to_geom: TangentGeometry, to_rep: Representation, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) Any
Parameters:
Return type:

Any

Convert tangent data between charts via Jacobian pushforward.

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"r": jnp.array(5.0), "theta": jnp.array(1.0), "phi": jnp.array(2.0)}
>>> at = {"r": jnp.array(3.0), "theta": jnp.array(0.5), "phi": jnp.array(0.0)}
>>> cxr.cconvert(v, cxc.sph3d, cxr.tangent_geom, cxr.coord_disp,
...              cxc.sph3d, cxr.tangent_geom, cxr.phys_disp, at=at)
{'r': Array(5., dtype=float64, ...),
 'theta': Array(3., dtype=float64, ...),
 'phi': Array(..., dtype=float64, ...)}
>>> v = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> at = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> cxr.cconvert(v, cxc.cart2d, cxr.coord_disp, cxc.polar2d, cxr.coord_disp, at=at)
{'r': Array(1., ...), 'theta': Array(0., ...)}
coordinax.representations.cconvert(from_vec: coordinax.vectors._src.point.Point, to_chart: AbstractChart, /, *, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.point.Point
Parameters:
Return type:

Any

Convert a point from one chart to another.

>>> import unxt as u
>>> import coordinax.main as cx
>>> vec = cx.Point.from_([1, 1, 1], "m")
>>> print(vec)
<Point: chart=Cart3D (x, y, z) [m]
    [1 1 1]>
>>> sph_vec = cx.cconvert(vec, cx.sph3d)
>>> print(sph_vec)
<Point: chart=Spherical3D (r[m], theta[rad], phi[rad])
    [1.732 0.955 0.785]>
coordinax.representations.cconvert(from_vec: coordinax.vectors._src.point.Point, from_chart: AbstractChart, to_chart: AbstractChart, /, *, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.point.Point
Parameters:
Return type:

Any

Convert a vector from one chart to another.

>>> import unxt as u
>>> import coordinax.main as cx
>>> vec = cx.Point.from_([1, 1, 1], "m")
>>> sph_vec = cx.cconvert(vec, cx.cart3d, cx.sph3d)
>>> print(sph_vec)
<Point: chart=Spherical3D (r[m], theta[rad], phi[rad])
    [1.732 0.955 0.785]>
coordinax.representations.cconvert(from_vec: coordinax.vectors._src.tangent.Tangent, to_chart: AbstractChart, /, *, at: Any = None, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Convert a tangent Tangent from one chart to another.

The at parameter provides the base point at which the tangent map (Jacobian pushforward) is evaluated. It may be a Point instance (whose .data is used) or a raw CDict.

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = cx.Tangent.from_(
...     {"x": u.Q(1.0, "m/s"), "y": u.Q(0.0, "m/s"), "z": u.Q(0.0, "m/s")},
...     cxc.cart3d, cxr.coord_basis, cxr.vel,
... )
>>> pt = cx.Point.from_([1.0, 0.0, 0.0], "m")
>>> v_sph = cx.cconvert(v, cxc.sph3d, at=pt)
>>> v_sph.chart
Spherical3D(M=Rn(3))
coordinax.representations.cconvert(from_vec: coordinax.vectors._src.tangent.Tangent, from_chart: AbstractChart, to_chart: AbstractChart, /, *, at: Any = None, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Convert a tangent Tangent from one chart to another (explicit from-chart).

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = cx.Tangent.from_(
...     {"x": u.Q(1.0, "m/s"), "y": u.Q(0.0, "m/s"), "z": u.Q(0.0, "m/s")},
...     cxc.cart3d, cxr.coord_basis, cxr.vel,
... )
>>> pt = cx.Point.from_([1.0, 0.0, 0.0], "m")
>>> v_sph = cx.cconvert(v, cxc.cart3d, cxc.sph3d, at=pt)
>>> v_sph.chart
Spherical3D(M=Rn(3))
coordinax.representations.cconvert(pv: coordinax.vectors._src.bundle.Coordinate, to_chart: AbstractChart, /, *, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.bundle.Coordinate
Parameters:
Return type:

Any

Convert a Coordinate to a new chart.

Delegates to {meth}`Coordinate.cconvert`.

>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> pt = cx.Point.from_([1.0, 0.0, 0.0], "m")
>>> pv = cx.Coordinate(point=pt)
>>> pv_sph = cx.cconvert(pv, cxc.sph3d)
>>> pv_sph.point.chart
Spherical3D(M=Rn(3))
Parameters:
Return type:

Any

coordinax.representations.cmap(*fixed_args: Any, **fixed_kw: Any)#

Return a partial function for vector conversion.

Convert a point from Cartesian coordinates to spherical coordinates:

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> import unxt as u

Define a map to convert a point from Cartesian coordinates to spherical coordinates:

>>> map = cxr.cmap(cxc.cart3d, cxr.point, cxc.sph3d)

Apply the map to a point:

>>> q = {"x": 1, "y": 2, "z": 3}
>>> map(q)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}
>>> q = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> map(q)
{'r': Q(3.74165739, 'm'), 'theta': Q(0.64052231, 'rad'),
 'phi': Q(1.10714872, 'rad')}

Apply to a Point:

>>> import coordinax.vectors as cxv
>>> vec = cxv.Point.from_(q, cxc.cart3d)
>>> map(vec)
Point(
  {'r': Q(3.74165739, 'm'), 'theta': Q(0.64052231, 'rad'),
   'phi': Q(1.10714872, 'rad')},
  chart=Spherical3D(M=Rn(3))
)
Parameters:
  • fixed_args (Any)

  • fixed_kw (Any)

Return type:

Any

coordinax.representations.guess_basis_kind(*args: Any, **kwargs: Any)#

Guess the basis kind of the given data.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> data = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> cxr.guess_basis_kind(data)
no_basis
coordinax.representations.guess_basis_kind(obj: AbstractBasis, /) AbstractBasis
Parameters:
Return type:

Any

Infer basis kind from an AbstractBasis object.

>>> import coordinax.representations as cxr
>>> basis = cxr.NoBasis()
>>> cxr.guess_basis_kind(basis) is basis
True
coordinax.representations.guess_basis_kind(dim: PhysicalType | tuple[PhysicalType, ...], /) AbstractBasis
Parameters:
Return type:

Any

Infer basis kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_basis_kind(u.dimension("length"))
no_basis
>>> cxr.guess_basis_kind((u.dimension("angle"), u.dimension("length")))
no_basis
coordinax.representations.guess_basis_kind(obj: AbstractQuantity, /) AbstractBasis
Parameters:
Return type:

Any

Infer basis kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_basis_kind(u.Q(1.0, "m"))
no_basis
coordinax.representations.guess_basis_kind(obj: dict[str, Any], /) AbstractBasis
Parameters:
Return type:

Any

Infer basis kind from the physical dimensions of a component dictionary.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> d1 = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> cxr.guess_basis_kind(d1)
no_basis
>>> d2 = {"lon": u.Q(1.0, "deg"), "lat": u.Q(2.0, "deg")}
>>> cxr.guess_basis_kind(d2)
no_basis
Parameters:
Return type:

Any

coordinax.representations.guess_geometry_kind(*args: Any, **kwargs: Any)#

Guess the geometry kind of the given data.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> data = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> cxr.guess_geometry_kind(data)
point_geom
coordinax.representations.guess_geometry_kind(obj: AbstractGeometry, /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from an AbstractGeometry object.

>>> import coordinax.representations as cxr
>>> geom = cxr.PointGeometry()
>>> cxr.guess_geometry_kind(geom) is geom
True
coordinax.representations.guess_geometry_kind(dim: PhysicalType | tuple[PhysicalType, ...], /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_geometry_kind(u.dimension("length"))
point_geom
>>> cxr.guess_geometry_kind((u.dimension("angle"), u.dimension("length")))
point_geom
>>> cxr.guess_geometry_kind(u.dimension("speed"))
tangent_geom
>>> cxr.guess_geometry_kind(u.dimension("angular speed"))
tangent_geom
>>> cxr.guess_geometry_kind(u.dimension("acceleration"))
tangent_geom
>>> cxr.guess_geometry_kind(u.dimension("angular acceleration"))
tangent_geom
coordinax.representations.guess_geometry_kind(obj: AbstractQuantity, /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_geometry_kind(u.Q(1.0, "m"))
point_geom
>>> cxr.guess_geometry_kind(u.Q(2.0, "m / s"))
tangent_geom
>>> cxr.guess_geometry_kind(u.Q(3.0, "m / s ** 2"))
tangent_geom
coordinax.representations.guess_geometry_kind(obj: dict[str, Any], /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from the physical dimensions of a component dictionary.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> d1 = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> cxr.guess_geometry_kind(d1)
point_geom
>>> d2 = {"lon": u.Q(1.0, "deg"), "lat": u.Q(2.0, "deg")}
>>> cxr.guess_geometry_kind(d2)
point_geom
>>> d3 = {"x": u.Q(1.0, "m / s"), "y": u.Q(2.0, "m / s")}
>>> cxr.guess_geometry_kind(d3)
tangent_geom
coordinax.representations.guess_geometry_kind(obj: dict[str, Any], chart: AbstractChart, /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from the physical dimensions of a component dictionary and chart.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> d = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> cxr.guess_geometry_kind(d, cxc.cart2d)
point_geom
coordinax.representations.guess_geometry_kind(obj: dict[str, Any], chart: ProlateSpheroidal3D, /) AbstractGeometry
Parameters:
Return type:

Any

Infer geometry kind from the physical dimensions of a component dictionary and chart.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> d = {"mu": u.Q(1, "km2"), "nu": u.Q(0.5, "km2"), "phi": u.Q(1, "deg")}
>>> chart = cxc.ProlateSpheroidal3D(Delta=u.StaticQuantity(1.0, "km"))
>>> cxr.guess_geometry_kind(d, chart)
point_geom
Parameters:
Return type:

Any

coordinax.representations.guess_rep(*args: Any, **kwargs: Any)#

Guess the representation of the given data.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> data = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> cxr.guess_rep(data)
point
coordinax.representations.guess_rep(obj: Representation, /) Representation
Parameters:
Return type:

Any

Infer representation from a Representation object.

>>> import coordinax.representations as cxr
>>> rep = cxr.point
>>> cxr.guess_rep(rep) is rep
True
coordinax.representations.guess_rep(obj: PointGeometry, /) Representation
Parameters:
Return type:

Any

Infer representation from a PointGeometry object.

>>> import coordinax.representations as cxr
>>> geom = cxr.PointGeometry()
>>> cxr.guess_rep(geom)
point
coordinax.representations.guess_rep(obj: PhysicalType | tuple[PhysicalType, ...] | AbstractQuantity | dict[str, Any], geom: PointGeometry, /) Representation
Parameters:
Return type:

Any

Infer point representation from data and an already-inferred PointGeometry.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> rep = cxr.guess_rep(u.dimension("length"), cxr.point_geom)
>>> rep
point
coordinax.representations.guess_rep(obj: PhysicalType | tuple[PhysicalType, ...] | AbstractQuantity | dict[str, Any], geom: TangentGeometry, /) Representation
Parameters:
Return type:

Any

Infer tangent representation from data and an already-inferred TangentGeometry.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> rep = cxr.guess_rep(u.dimension("speed"), cxr.tangent_geom)
>>> rep.geom_kind
tangent_geom
>>> rep.semantic_kind
vel
>>> rep = cxr.guess_rep(u.dimension("acceleration"), cxr.tangent_geom)
>>> rep.geom_kind
tangent_geom
>>> rep.semantic_kind
acc
coordinax.representations.guess_rep(obj: PhysicalType | tuple[PhysicalType, ...] | AbstractQuantity | dict[str, Any], /) Representation
Parameters:
Return type:

Any

Infer representation from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> rep = cxr.guess_rep(u.dimension("length"))
>>> rep
point
>>> rep = cxr.guess_rep(u.Q(1.0, "m"))
>>> rep
point
>>> rep = cxr.guess_rep(u.dimension("speed"))
>>> rep.geom_kind
tangent_geom
>>> rep.semantic_kind
vel
>>> rep = cxr.guess_rep(u.Q(1.0, "m / s"))
>>> rep.geom_kind
tangent_geom
>>> rep.semantic_kind
vel
>>> rep = cxr.guess_rep(u.dimension("acceleration"))
>>> rep.geom_kind
tangent_geom
>>> rep.semantic_kind
acc
coordinax.representations.guess_rep(obj: Any, chart: AbstractChart, /) Representation
Parameters:
Return type:

Any

Infer representation from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> rep = cxr.guess_rep(u.dimension("length"))
>>> rep
point
>>> rep = cxr.guess_rep(u.Q(1.0, "m"))
>>> rep
point
coordinax.representations.guess_rep(obj: Any, chart: AbstractChart, geom: PointGeometry, /) Representation
Parameters:
Return type:

Any

Infer representation.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> rep = cxr.guess_rep(u.Q(1.0, "m"), cxc.cart2d, point_geom)
>>> rep
point
coordinax.representations.guess_rep(obj: Any, chart: AbstractChart, geom: TangentGeometry, /) Representation
Parameters:
Return type:

Any

Infer representation.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

Speed dimensions infer Velocity:

>>> rep = cxr.guess_rep(u.Q([1.0, 2.0], "m / s"), cxc.cart2d, cxr.TangentGeometry())
>>> rep
coord_vel

Length dimensions would infer Location in general, but TangentGeometry requires a tangent semantic kind, so Displacement is returned instead:

>>> rep = cxr.guess_rep(u.Q([1.0, 2.0], "m"), cxc.cart2d, cxr.TangentGeometry())
>>> rep
coord_disp
Parameters:
Return type:

Any

coordinax.representations.guess_semantic_kind(*args: Any, **kwargs: Any)#

Guess the semantic kind of the given data.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> data = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> cxr.guess_semantic_kind(data)
loc
coordinax.representations.guess_semantic_kind(obj: AbstractSemanticKind, /) AbstractSemanticKind
Parameters:
Return type:

Any

Infer semantic kind from an AbstractSemanticKind object.

>>> import coordinax.representations as cxr
>>> sem = cxr.Location()
>>> cxr.guess_semantic_kind(sem) is sem
True
coordinax.representations.guess_semantic_kind(dim: PhysicalType | tuple[PhysicalType, ...], /) AbstractSemanticKind
Parameters:
Return type:

Any

Infer semantic kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_semantic_kind(u.dimension("length"))
loc
>>> cxr.guess_semantic_kind((u.dimension("angle"), u.dimension("length")))
loc
>>> cxr.guess_semantic_kind(u.dimension("speed"))
vel
>>> cxr.guess_semantic_kind(u.dimension("angular speed"))
vel
>>> cxr.guess_semantic_kind(u.dimension("acceleration"))
acc
>>> cxr.guess_semantic_kind(u.dimension("angular acceleration"))
acc
coordinax.representations.guess_semantic_kind(obj: AbstractQuantity, /) AbstractSemanticKind
Parameters:
Return type:

Any

Infer semantic kind from the physical dimension of a Quantity.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> cxr.guess_semantic_kind(u.Q(1.0, "m"))
loc
coordinax.representations.guess_semantic_kind(obj: dict[str, Any], /) AbstractSemanticKind
Parameters:
Return type:

Any

Infer semantic kind from the physical dimensions of a component dictionary.

>>> import unxt as u
>>> import coordinax.representations as cxr
>>> d1 = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> cxr.guess_semantic_kind(d1)
loc
>>> d2 = {"lon": u.Q(1.0, "deg"), "lat": u.Q(2.0, "deg")}
>>> cxr.guess_semantic_kind(d2)
loc
Parameters:
Return type:

Any

coordinax.representations.subtract(*args: Any, **kwargs: Any)#

Subtract two coordinate data objects.

This is an abstract API definition. See the main coordinax package for concrete implementations.

coordinax.representations.subtract(lhs: Any, lhs_chart: AbstractChart, lhs_rep: Representation, rhs: Any, rhs_chart: AbstractChart, rhs_rep: Representation, /) Any
Parameters:
Return type:

Any

Subtract two coordinate data objects via Cartesian round-trip.

Both operands are converted to the ambient Cartesian chart of lhs_chart, subtracted component-wise, then converted back to lhs_chart. If lhs_chart is already Cartesian (or has no global Cartesian), rhs is converted into lhs_chart and subtracted directly.

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> import unxt as u
>>> p1 = {"x": u.Q(4, "m"), "y": u.Q(5, "m"), "z": u.Q(6, "m")}
>>> p2 = {"x": u.Q(1, "m"), "y": u.Q(2, "m"), "z": u.Q(3, "m")}
>>> cxr.subtract(p1, cxc.cart3d, cxr.point, p2, cxc.cart3d, cxr.point)
{'x': Q(3, 'm'), 'y': Q(3, 'm'), 'z': Q(3, 'm')}
coordinax.representations.subtract(lhs: coordinax.vectors._src.point.Point, rhs: coordinax.vectors._src.point.Point, /) coordinax.vectors._src.point.Point
Parameters:
Return type:

Any

Subtract two vectors.

For non-Cartesian charts the operation converts both operands to the ambient Cartesian chart, subtracts there, and converts the result back to the lhs chart. For Cartesian charts the subtraction is direct.

The result keeps the lhs chart and representation.

>>> import coordinax.main as cx
>>> v1 = cx.Point.from_([4, 5, 6], "m")
>>> v2 = cx.Point.from_([1, 2, 3], "m")
>>> print(cxr.subtract(v1, v2))
<Point: chart=Cart3D (x, y, z) [m]
    [3 3 3]>
coordinax.representations.subtract(lhs: coordinax.vectors._src.tangent.Tangent, rhs: coordinax.vectors._src.tangent.Tangent, /) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Subtract two tangent vectors component-wise.

Tangent spaces are genuine vector spaces: subtraction is component-wise in any chart basis (no Cartesian round-trip is needed or correct). Both operands must share the same representation (chart + basis + semantic).

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v1 = cx.Tangent.from_(
...     {"x": u.Q(4.0, "m/s"), "y": u.Q(5.0, "m/s"), "z": u.Q(6.0, "m/s")},
...     cxc.cart3d, cxr.coord_vel,
... )
>>> v2 = cx.Tangent.from_(
...     {"x": u.Q(1.0, "m/s"), "y": u.Q(2.0, "m/s"), "z": u.Q(3.0, "m/s")},
...     cxc.cart3d, cxr.coord_vel,
... )
>>> result = cxr.subtract(v1, v2)
>>> result["x"]
Q(3., 'm / s')
Parameters:
Return type:

Any

coordinax.representations.tangent_map(*args: Any, **kwargs: Any)#

Compute the tangent map (Jacobian) of a chart transition.

This is an abstract API definition. See the main coordinax package for concrete implementations.

Examples

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
coordinax.representations.tangent_map(v: Any, from_chart: AbstractChart, basis: CoordinateBasis, to_chart: AbstractChart, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Push a tangent vector forward from one chart to another.

Applies the Jacobian of the chart transition map to the tangent vector components v, evaluated at the base point at.

Convert a tangent vector from Cartesian to polar 2D at the point (1, 0):

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> at = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> cxr.tangent_map(v, cxc.cart2d, cxr.coord_disp, cxc.polar2d, at=at)
{'r': Array(1., dtype=float64), 'theta': Array(0., dtype=float64)}
coordinax.representations.tangent_map(v: dict[str, Any], from_chart: AbstractChart, basis: PhysicalBasis, to_chart: AbstractChart, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Push a tangent vector forward in physical-basis components.

This dispatch applies the tangent-map pushforward while preserving the physical-basis convention by composing three steps:

  1. convert source components from physical basis to coordinate basis,

  2. apply the chart Jacobian pushforward,

  3. convert target components back to physical basis.

Convert a physical-basis tangent vector from Cartesian to spherical 3D:

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"x": u.Q(1, "m/s"), "y": u.Q(0, "m/s"), "z": u.Q(0, "m/s")}
>>> at = {"x": u.Q(1, "m"), "y": u.Q(0, "m"), "z": u.Q(0, "m")}
>>> cxr.tangent_map(v, cxc.cart3d, cxr.phys_basis, cxc.sph3d, at=at)
{'r': Q(1., 'm / s'), 'theta': Q(-0., 'm / s'), 'phi': Q(0., 'm / s')}

The same call can be made using a physical representation:

>>> v = {"x": 1, "y": 0, "z": 0}
>>> at = {"x": 1, "y": 0, "z": 0}
>>> usys = u.unitsystems.si
>>> cxr.tangent_map(v, cxc.cart3d, cxr.phys_disp, cxc.sph3d, at=at, usys=usys)
{'r': Array(1., dtype=float64), 'theta': Array(0., dtype=float64),
 'phi': Array(0., dtype=float64)}
coordinax.representations.tangent_map(v: Any, from_chart: AbstractChart, from_rep: Representation, to_chart: AbstractChart, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Push a tangent vector forward from one chart to another.

Applies the Jacobian of the chart transition map to the tangent vector components v, evaluated at the base point at.

Convert a tangent vector from Cartesian to polar 2D at the point (1, 0):

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> at = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> cxr.tangent_map(v, cxc.cart2d, cxr.coord_disp, cxc.polar2d, at=at)
{'r': Array(1., dtype=float64), 'theta': Array(0., dtype=float64)}
coordinax.representations.tangent_map(v: Any, from_chart: AbstractChart, from_rep: Representation, to_chart: AbstractChart, to_rep: Representation, /, *, at: dict[str, Any] | None = None, usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Push a tangent vector forward from one chart to another.

Applies the Jacobian of the chart transition map to the tangent vector components v, evaluated at the base point at.

Convert a tangent vector from Cartesian to polar 2D at the point (1, 0):

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> at = {"x": jnp.array(1.0), "y": jnp.array(0.0)}
>>> cxr.tangent_map(v, cxc.cart2d, cxr.coord_disp,
...                 cxc.polar2d, cxr.coord_disp, at=at)
{'r': Array(1., dtype=float64), 'theta': Array(0., dtype=float64)}
Parameters:
Return type:

Any

coordinax.representations.change_basis(*args: Any, **kwargs: Any)#

Change the basis of a tangent vector’s components.

Examples

>>> import coordinax.representations as cxr
>>> import coordinax.charts as cxc
>>> v = {"x": 1.0, "y": 0.0}
>>> at = {"x": 1.0, "y": 0.0}
>>> cxr.change_basis(v, cxc.cart2d, cxr.coord_basis, cxr.phys_basis, at=at)
{'x': 1.0, 'y': 0.0}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, M: AbstractManifold, from_basis: CoordinateBasis, to_basis: PhysicalBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from coordinate basis to physical basis using a manifold.

Retrieves the manifold’s metric and applies the appropriate transformation. For diagonal metrics (e.g. coordinax.manifolds.FlatMetric in orthogonal charts) the fast scale-factor path is taken; for general metrics (e.g. coordinax.manifolds.PullbackMetric) the Cholesky vielbein $E = L^top$ is used.

>>> import jax.numpy as jnp
>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.manifolds as cxm
>>> import coordinax.representations as cxr

Euclidean 3-D manifold in spherical coordinates (diagonal metric):

>>> M3 = cxm.R3
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(1, "rad/s"), "phi": u.Q(2, "rad/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, M3, cxr.coord_basis, cxr.phys_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(3, 'm / s'), 'phi': Q(2.87655323, 'm / s')}

Embedded two-sphere manifold — non-diagonal PullbackMetric:

>>> M = cxm.EmbeddedManifold(
...     intrinsic=cxm.S2, ambient=cxm.R3,
...     embed_map=cxm.TwoSphereIn3D(radius=u.Q(1.0, "km")),
... )
>>> v = {"theta": u.Q(1.0, "rad/s"), "phi": u.Q(2.0, "rad/s")}
>>> at = {"theta": u.Q(jnp.pi / 3, "rad"), "phi": u.Q(0.0, "rad")}
>>> cxr.change_basis(v, cxc.sph2, M, cxr.coord_basis, cxr.phys_basis, at=at)
{'theta': Q(..., 'km / s'), 'phi': Q(..., 'km / s')}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, M: AbstractManifold, from_basis: PhysicalBasis, to_basis: CoordinateBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from physical basis to coordinate basis using a manifold.

Retrieves the manifold’s metric and applies the inverse transformation. For diagonal metrics the fast scale-factor path is taken; for general metrics the Cholesky vielbein $E = L^top$ is solved as a triangular system $v = E^{-1}hat{v}$.

>>> import jax.numpy as jnp
>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.manifolds as cxm
>>> import coordinax.representations as cxr

Euclidean 3-D manifold in spherical coordinates (diagonal metric):

>>> M3 = cxm.R3
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(3, "m/s"), "phi": u.Q(2.876553, "m/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0.5, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, M3, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(1., 'rad / s'), 'phi': Q(1.99999..., 'rad / s')}

Embedded two-sphere manifold — non-diagonal PullbackMetric:

>>> M = cxm.EmbeddedManifold(
...     intrinsic=cxm.S2, ambient=cxm.R3,
...     embed_map=cxm.TwoSphereIn3D(radius=u.Q(1.0, "km")),
... )
>>> v = {"theta": u.Q(1.0, "km/s"), "phi": u.Q(2.0, "km/s")}
>>> at = {"theta": u.Q(jnp.pi / 3, "rad"), "phi": u.Q(0.0, "rad")}
>>> cxr.change_basis(v, cxc.sph2, M, cxr.phys_basis, cxr.coord_basis, at=at)
{'theta': Q(..., 'rad / s'), 'phi': Q(..., 'rad / s')}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, from_basis: PhysicalBasis, to_basis: CoordinateBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from physical basis to coordinate basis using the chart’s manifold.

Falls back to chart.M when no explicit manifold is supplied.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(3, "m/s"), "phi": u.Q(2.876553, "m/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0.5, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(1., 'rad / s'), 'phi': Q(1.99999984, 'rad / s')}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, from_basis: NoBasis | CoordinateBasis, to_basis: CoordinateBasis, /, **kw: Any) dict[str, Any]
Parameters:
Return type:

Any

Reinterpret unknown basis components as coordinate-basis components.

This is an identity on component values: no numeric basis transform is possible from an unknown basis, so values are preserved and only the representation-level basis label changes.

>>> import jax.numpy as jnp
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

The conversion is an identity map on values:

>>> v = {"x": jnp.array(1.0), "y": jnp.array(-2.0)}
>>> at = {"x": jnp.array(0.0), "y": jnp.array(0.0)}
>>> cxr.change_basis(v, cxc.cart2d, cxr.no_basis, cxr.coord_basis, at=at)
{'x': Array(1., dtype=float64, ...), 'y': Array(-2., dtype=float64, ...)}

“at” is accepted and ignored for this basis-only reinterpretation:

>>> cxr.change_basis(v, cxc.cart2d, cxr.no_basis, cxr.coord_basis)
{'x': Array(1., dtype=float64, ...), 'y': Array(-2., dtype=float64, ...)}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, from_basis: NoBasis | PhysicalBasis, to_basis: PhysicalBasis, /, **kw: Any) dict[str, Any]
Parameters:
Return type:

Any

Reinterpret unknown basis components as physical-basis components.

This conversion is only well-defined when all components share the same physical dimension. No numeric transform is applied; values are preserved.

>>> import jax.numpy as jnp
>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

Same-dimension components (even with different units) are accepted:

>>> v = {"x": u.Q(1.0, "m / s"), "y": u.Q(2.0, "km / s")}
>>> at = {"x": jnp.array(0.0), "y": jnp.array(0.0)}
>>> cxr.change_basis(v, cxc.cart2d, cxr.no_basis, cxr.phys_basis, at=at)
{'x': Q(1., 'm / s'), 'y': Q(2., 'km / s')}

Mixed dimensions are rejected:

>>> v_bad = {"x": u.Q(1.0, "m / s"), "y": u.Q(2.0, "m")}
>>> cxr.change_basis(v_bad, cxc.cart2d, cxr.no_basis, cxr.phys_basis)
Traceback (most recent call last):
...
ValueError: change_basis from NoBasis to PhysicalBasis requires all
components to have the same dimension, got ...
coordinax.representations.change_basis(v: dict[str, Any], chart: Cart0D | Cart1D | Cart2D | Cart3D | CartND, from_basis: CoordinateBasis | PhysicalBasis, to_basis: CoordinateBasis | PhysicalBasis, /, **kw: Any) dict[str, Any]
Parameters:
Return type:

Any

Change the basis used to interpret tangent components.

In a Cartesian chart every scale factor equals one,

$$h_x = h_y = h_z = 1,$$

so the coordinate basis vectors are already unit vectors ($hat{e}_i = partial_i$) and the transformation matrix is the identity ($H = I$). The coordinate basis is the physical basis, so this conversion is always the identity map and v is returned unchanged regardless of the direction of the conversion.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

Coordinate basis to physical basis in a 2-D Cartesian chart — identity:

>>> v = {"x": u.Q(3.0, "m/s"), "y": u.Q(4.0, "m/s")}
>>> at = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> cxr.change_basis(v, cxc.cart2d, cxr.coord_basis, cxr.phys_basis, at=at)
{'x': Q(3., 'm / s'), 'y': Q(4., 'm / s')}

The reverse direction is equally a no-op:

>>> cxr.change_basis(v, cxc.cart2d, cxr.phys_basis, cxr.coord_basis, at=at)
{'x': Q(3., 'm / s'), 'y': Q(4., 'm / s')}

Works for any Cartesian dimensionality; at is optional:

>>> v3 = {"x": u.Q(1.0, "m/s"), "y": u.Q(2.0, "m/s"), "z": u.Q(3.0, "m/s")}
>>> cxr.change_basis(v3, cxc.cart3d, cxr.coord_basis, cxr.phys_basis)
{'x': Q(1., 'm / s'), 'y': Q(2., 'm / s'), 'z': Q(3., 'm / s')}
coordinax.representations.change_basis(v: dict[str, Any], chart: Cart0D | Cart1D | Cart2D | Cart3D | CartND, M: AbstractManifold, from_basis: AbstractLinearBasis, to_basis: AbstractLinearBasis, /, **kw: Any) dict[str, Any]
Parameters:
Return type:

Any

Change the basis used to interpret tangent components.

In a Cartesian chart every scale factor equals one,

$$h_x = h_y = h_z = 1,$$

so the coordinate basis vectors are already unit vectors ($hat{e}_i = partial_i$) and the transformation matrix is the identity ($H = I$). The coordinate basis is the physical basis, so this conversion is always the identity map and v is returned unchanged regardless of the direction of the conversion.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.manifolds as cxm
>>> import coordinax.representations as cxr

Coordinate basis to physical basis in a 2-D Cartesian chart — identity:

>>> v = {"x": u.Q(3.0, "m/s"), "y": u.Q(4.0, "m/s")}
>>> at = {"x": u.Q(1.0, "m"), "y": u.Q(2.0, "m")}
>>> M2 = cxm.R2
>>> cxr.change_basis(v, cxc.cart2d, M2, cxr.coord_basis, cxr.phys_basis, at=at)
{'x': Q(3., 'm / s'), 'y': Q(4., 'm / s')}

The reverse direction is equally a no-op:

>>> cxr.change_basis(v, cxc.cart2d, M2, cxr.phys_basis, cxr.coord_basis, at=at)
{'x': Q(3., 'm / s'), 'y': Q(4., 'm / s')}

Works for any Cartesian dimensionality; at is optional:

>>> v3 = {"x": u.Q(1.0, "m/s"), "y": u.Q(2.0, "m/s"), "z": u.Q(3.0, "m/s")}
>>> cxr.change_basis(v3, cxc.cart3d, cxm.R3, cxr.coord_basis, cxr.phys_basis)
{'x': Q(1., 'm / s'), 'y': Q(2., 'm / s'), 'z': Q(3., 'm / s')}
coordinax.representations.change_basis(v: dict[str, Any], chart: Spherical3D, from_basis: CoordinateBasis, to_basis: PhysicalBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from coordinate basis to physical basis in a 3-D spherical chart.

In spherical coordinates $(r, theta, phi)$ the scale factors are

$$h_r = 1, quad h_theta = r, quad h_phi = r sintheta,$$

so the transformation matrix is

$$ H = begin{pmatrix}

1 & 0 & 0 \ 0 & r & 0 \ 0 & 0 & rsintheta

end{pmatrix}. $$

Given coordinate-basis components $(v^r, v^theta, v^phi)$, the physical-basis components are

$$ hat{v} = H v

implies hat{v}^r = v^r, quad hat{v}^theta = r, v^theta, quad hat{v}^phi = rsintheta, v^phi.

$$

Examples

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(1, "rad/s"), "phi": u.Q(2, "rad/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxr.coord_basis, cxr.phys_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(3, 'm / s'), 'phi': Q(2.87655323, 'm / s')}
>>> v = {"r": 5, "theta": 1, "phi": 2}  # unitless
>>> at = {"r": 3, "theta": 0.5, "phi": 0}  # unitless
>>> cxr.change_basis(v, cxc.sph3d, cxr.coord_basis, cxr.phys_basis, at=at)
{'r': 5, 'theta': 3, 'phi': Array(2.87655323, dtype=float64, ...)}
coordinax.representations.change_basis(v: dict[str, Any], chart: Spherical3D, M: EuclideanManifold, from_basis: CoordinateBasis, to_basis: PhysicalBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from coordinate basis to physical basis in a 3-D spherical chart.

Delegates to the chart-specific implementation for coordinax.manifolds.EuclideanManifold.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.manifolds as cxm
>>> import coordinax.representations as cxr
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(1, "rad/s"), "phi": u.Q(2, "rad/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxm.R3, cxr.coord_basis, cxr.phys_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(3, 'm / s'), 'phi': Q(2.87655323, 'm / s')}
>>> v = {"r": 5, "theta": 1, "phi": 2}  # unitless
>>> at = {"r": 3, "theta": 0.5, "phi": 0}  # unitless
>>> cxr.change_basis(v, cxc.sph3d, cxm.R3, cxr.coord_basis, cxr.phys_basis, at=at)
{'r': 5, 'theta': 3, 'phi': Array(2.87655323, dtype=float64, ...)}
coordinax.representations.change_basis(v: dict[str, Any], chart: Spherical3D, from_basis: PhysicalBasis, to_basis: CoordinateBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from physical basis to coordinate basis in a 3-D spherical chart.

In spherical coordinates $(r, theta, phi)$ the inverse transformation matrix is

$$ H^{-1} = begin{pmatrix}

1 & 0 & 0 \ 0 & 1/r & 0 \ 0 & 0 & 1/(rsintheta)

end{pmatrix}. $$

Given physical-basis components $(hat{v}^r, hat{v}^theta, hat{v}^phi)$, the coordinate-basis components are

$$ v = H^{-1} hat{v}

implies v^r = hat{v}^r, quad v^theta = hat{v}^theta / r, quad v^phi = hat{v}^phi / (rsintheta).

$$

Examples

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(3, "m/s"), "phi": u.Q(2.876553, "m/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0.5, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(1., 'rad / s'), 'phi': Q(1.99999984, 'rad / s')}
>>> v = {"r": 5, "theta": 3, "phi": 2.876553}  # unitless
>>> at = {"r": 3, "theta": 0.5, "phi": 0.5}  # unitless
>>> cxr.change_basis(v, cxc.sph3d, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': 5, 'theta': 1.0, 'phi': Array(1.99999984, dtype=float64, ...)}
coordinax.representations.change_basis(v: dict[str, Any], chart: Spherical3D, M: EuclideanManifold, from_basis: PhysicalBasis, to_basis: CoordinateBasis, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change from physical basis to coordinate basis in a 3-D spherical chart.

Delegates to the chart-specific implementation for EuclideanManifold.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.manifolds as cxm
>>> import coordinax.representations as cxr
>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(3, "m/s"), "phi": u.Q(2.876553, "m/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0.5, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxm.R3, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(1., 'rad / s'), 'phi': Q(1.99999984, 'rad / s')}
>>> v = {"r": 5, "theta": 3, "phi": 2.876553}  # unitless
>>> at = {"r": 3, "theta": 0.5, "phi": 0.5}  # unitless
>>> cxr.change_basis(v, cxc.sph3d, cxm.R3, cxr.phys_basis, cxr.coord_basis, at=at)
{'r': 5, 'theta': 1.0, 'phi': Array(1.99999984, dtype=float64, ...)}
coordinax.representations.change_basis(v: dict[str, Any], chart: AbstractChart, from_rep: Representation, to_rep: Representation, /, *, at: dict[str, Any], usys: AbstractUnitSystem | None = None) dict[str, Any]
Parameters:
Return type:

Any

Change basis using source and/or target Representation objects.

This is a convenience overload: the caller may pass full Representation objects for from_rep/to_rep instead of bare AbstractBasis instances. The basis is extracted from each argument and the appropriate change_basis() overload is called.

>>> import unxt as u
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

Coordinate-basis displacement to physical-basis displacement in a spherical chart, passing full Representation objects:

>>> v = {"r": u.Q(5, "m/s"), "theta": u.Q(1, "rad/s"), "phi": u.Q(2, "rad/s")}
>>> at = {"r": u.Q(3, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0, "rad")}
>>> cxr.change_basis(v, cxc.sph3d, cxr.coord_disp, cxr.phys_disp, at=at)
{'r': Q(5, 'm / s'), 'theta': Q(3, 'm / s'), 'phi': Q(2.87655323, 'm / s')}
>>> v = {"r": 5, "theta": 1, "phi": 2}  # unitless
>>> at = {"r": 3, "theta": 0.5, "phi": 0}  # unitless
>>> cxr.change_basis(v, cxc.sph3d, cxr.coord_disp, cxr.phys_disp, at=at)
{'r': 5, 'theta': 3, 'phi': Array(2.87655323, dtype=float64, ...)}
coordinax.representations.change_basis(v: coordinax.vectors._src.tangent.Tangent, to_basis: AbstractLinearBasis, /, *, at: Any = None, usys: AbstractUnitSystem | None = None) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Change the basis of a Tangent vector.

Converts the component data from the current basis to to_basis using the registered change_basis overload for dicts, then returns a new Tangent with the updated data and basis.

The at parameter provides the base point at which the scale factors are evaluated. It may be a Point instance (whose .data is used) or a raw CDict.

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr

Convert a coordinate-basis spherical velocity to the physical basis:

>>> point = cx.Point.from_(
...     {"r": u.Q(1.0, "m"), "theta": u.Q(0.5, "rad"), "phi": u.Q(0.0, "rad")},
...     cxc.sph3d,
... )
>>> vel = cx.Tangent.from_({"r": u.Q(1.0, "m/s"), "theta": u.Q(0.0, "rad/s"),
...     "phi": u.Q(0.0, "rad/s")}, cxc.sph3d, cxr.coord_vel)
>>> vel_phys = cxr.change_basis(vel, cxr.phys_basis, at=point)
>>> vel_phys.basis
phys_basis
>>> vel_phys.rep
phys_vel
coordinax.representations.change_basis(v: coordinax.vectors._src.point.Point, to_basis: AbstractLinearBasis, /, *, at: Any = None, usys: Any = None) coordinax.vectors._src.tangent.Tangent
Parameters:
Return type:

Any

Promote a Point to a Tangent with Displacement semantics.

The component data are unchanged; only the geometric interpretation is recast from a manifold point (PointGeometry) to a tangent-space displacement vector (TangentGeometry, Displacement). The resulting Tangent carries the same chart and frame as the input Point, and its basis is to_basis.

The at and usys parameters are accepted for API consistency but are not used.

>>> import unxt as u
>>> import coordinax.main as cx
>>> import coordinax.representations as cxr
>>> pt = cx.Point.from_([1.0, 2.0, 3.0], "m")
>>> disp = cxr.change_basis(pt, cxr.coord_basis)
>>> disp.semantic
dpl
>>> disp.basis
coord_basis
>>> disp.chart == pt.chart
True
Parameters:
Return type:

Any

final class coordinax.representations.Representation(geom_kind: GeomT, basis: BasisT, semantic_kind: SemanticT)#

Bases: Generic[GeomT, BasisT, SemanticT]

Representation of geometric component data.

A representation specifies what kind of geometric object component data is meant to represent, independently of the chart used to write down the coordinates or components.

In coordinax, a representation is the ordered triple

$$ R = (K, B, S), $$

where:

  • $K$ is the geometric kind (coordinax.representations.AbstractGeometry),

  • $B$ is the basis kind (coordinax.representations.AbstractBasis), and

  • $S$ is the semantic kind (coordinax.representations.AbstractSemanticKind).

Thus a Representation answers three distinct questions:

  1. What sort of geometric object is this? For example, a point.

  2. In what basis are its components written? For example, no basis for affine point data.

  3. What does the object mean? For example, a location on a manifold.

A representation is therefore not the same thing as a chart.

  • A chart specifies how a manifold is coordinatized locally: component names, ordering, dimensionalities, and the coordinate map into $mathbb{R}^n$.

  • A representation specifies how the data should be interpreted geometrically.

Equivalently: the chart determines the coordinate system, while the representation determines the geometric role of the data written in that coordinate system.

For the current point-focused design, the canonical representation is

$$ (mathrm{PointGeometry},, mathrm{NoBasis},, mathrm{Location}). $$

This indicates that the data represents a point on a manifold, that it does not live in a basis-dependent linear space, and that its semantic meaning is a location.

Parameters:

Examples

Construct the canonical point representation directly:

>>> import coordinax.representations as cxr
>>> rep = cxr.Representation(cxr.PointGeometry(), cxr.NoBasis(), cxr.Location())

coordinax also provides the predefined point representation:

>>> rep == cxr.point
True

Use the representation with cconvert to convert point data between charts while preserving its geometric interpretation:

>>> import coordinax.charts as cxc
>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> cxr.cconvert(p, cxc.cart3d, cxr.point, cxc.sph3d, cxr.point)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output represents the same point, but in the target chart.

Notes

Representation is a static, immutable descriptor object. It carries no runtime numerical data itself; it only describes how such data should be interpreted.

geom_kind: GeomT#

Geometric kind of the represented object.

basis: BasisT#

Basis kind in which components are expressed.

semantic_kind: SemanticT#

Semantic interpretation attached to the represented object.

class coordinax.representations.AbstractGeometry#

Bases: object

Abstract base class for geometric kind.

A geometric kind specifies the underlying geometric type of the data, independent of any coordinate chart or basis in which components may be written.

In the representation model used by coordinax, the full representation is determined by three orthogonal pieces of information:

  1. the geometric kind,

  2. the basis, and

  3. the semantic kind.

This class provides the first of these pieces: it answers the question “what sort of geometric object is this?” Examples include a point.

Mathematical Role:

The geometric kind determines the abstract transformation behavior of the object under coordinate changes.

  • For a point, coordinates transform by the chart transition map.

  • For a tangent vector, components transform by the pushforward (Jacobian) of the coordinate map.

  • For a cotangent vector, components transform by the pullback.

Thus geometric kind is distinct from:

  • a chart, which specifies how local coordinates are assigned, and

  • a semantic kind, which specifies the interpretation of an object of a given geometric type (for example, location, displacement, velocity, or acceleration).

Examples

>>> import coordinax.main as cx

Construct the point geometry object directly:

>>> geom = cx.PointGeometry()

With the geometry object, we construct a full representation for point data:

>>> rep = cx.Representation(geom, cx.no_basis, cx.loc)

The representation can then be used with cconvert to convert point data between charts while preserving the represented geometric object:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> cx.cconvert(p, cx.cart3d, rep, cx.sph3d, rep)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output is still point data, but expressed in the target chart.

Notes

This is a static dispatch object and carries no runtime numerical data. Concrete subclasses should represent immutable geometric categories.

canonical_name: ClassVar[str | None] = None#

Canonical name for the geometric kind.

final class coordinax.representations.PointGeometry#

Bases: AbstractGeometry

Point geometric kind.

A point geometry indicates that component data should be interpreted as the coordinates of a point on a manifold, rather than as components of a vector- or covector-like object.

Mathematical Definition:

Let $M$ be a smooth manifold. A point is an element $p in M$.

Point data is therefore affine, not linear: points do not in general form a vector space, so operations such as adding two points are not geometrically defined. Under a change of chart, point coordinates transform by the ordinary chart transition map.

Examples

Construct the point geometry object directly:

>>> import coordinax.representations as cxr
>>> geom = cxr.PointGeometry()

Use it inside a full representation for point data:

>>> rep = cxr.Representation(geom, cxr.no_basis, cxr.loc)

The representation can then be used with cconvert to convert point data between charts while preserving the represented geometric object:

>>> import coordinax.charts as cxc
>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> cxr.cconvert(p, cxc.cart3d, rep, cxc.sph3d, rep)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output is still point data, but expressed in the target chart.

Notes

PointGeometry describes only the geometric kind. The chart still determines the coordinate system, component names, and coordinate domains.

canonical_name: ClassVar = 'point_geom'#

Canonical name for the point geometry kind.

final class coordinax.representations.TangentGeometry#

Bases: AbstractGeometry

Tangent-vector geometric kind.

A tangent geometry indicates that component data should be interpreted as the components of a tangent vector at a point on a manifold, rather than as coordinates of a point.

Mathematical Definition:

Let $M$ be a smooth manifold and $p in M$ a point. A tangent vector is an element $v in T_p M$ of the tangent space at $p$.

Tangent data is linear: tangent vectors form a vector space, so adding two tangent vectors at the same point is geometrically well-defined. Under a change of chart, tangent vector components transform by the Jacobian (pushforward) of the chart transition map.

Examples

Construct the tangent geometry object directly:

>>> import coordinax.representations as cxr
>>> geom = cxr.TangentGeometry()

Use it inside a full representation for tangent data:

>>> rep = cxr.Representation(geom, cxr.coord_basis, cxr.dpl)

Notes

This is a static dispatch object and carries no runtime numerical data. TangentGeometry describes only the geometric kind. The chart still determines the coordinate system, component names, and coordinate domains.

canonical_name: ClassVar = 'tangent_geom'#

Canonical name for the tangent geometry kind.

class coordinax.representations.AbstractBasis#

Bases: object

Abstract base class for basis kind.

A basis kind specifies the component basis in which data is expressed, when such a choice is meaningful, independent of the underlying geometric type and independent of which coordinate chart is used.

In the representation model used by coordinax, the full representation is determined by three orthogonal pieces of information:

  1. the geometric kind,

  2. the basis, and

  3. the semantic kind.

This class provides the second of these pieces: it answers the question “in what basis are the components written?” Examples include no basis for affine point data, and coordinate and physical bases for tangent or cotangent objects.

Mathematical Role:

The basis determines how component values are interpreted within a fixed geometric type.

  • For a point, there is typically no basis choice: points are affine objects rather than vectors in a linear space.

  • For a tangent vector, one may distinguish between components written in a coordinate basis and components written in an orthonormal physical basis.

  • For a cotangent vector, an analogous distinction may be made between the dual coordinate basis and other dual frames.

Thus basis is distinct from:

  • a chart, which specifies how local coordinates are assigned,

  • a geometry kind, which specifies what sort of geometric object the data represents, and

  • a semantic kind, which specifies the interpretation of an object of a given geometric type (for example, location, displacement, velocity, or acceleration).

Examples

>>> import coordinax.main as cx

Construct the no-basis object directly:

>>> basis = cx.NoBasis()

With the basis object, we construct a full representation for point data:

>>> rep = cx.Representation(cx.point_geom, basis, cx.loc)

The representation can then be used with coordinax.representations.cconvert to convert point data between charts while preserving both the represented geometric object and the fact that point data has no basis-dependent components:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> q = cx.cconvert(p, cx.cart3d, rep, cx.sph3d)

The output q is still point data, but expressed in the target chart.

Notes

This is a static dispatch object and carries no runtime numerical data. Concrete subclasses should represent immutable basis categories.

canonical_name: ClassVar[str | None] = None#

Canonical name for the basis kind.

final class coordinax.representations.NoBasis#

Bases: AbstractBasis

No-basis kind.

A no-basis kind indicates that the represented data does not live in a basis-dependent linear space, so there is no meaningful choice of components with respect to a basis.

Mathematical Definition:

NoBasis is used for geometric objects whose representation is not given by expansion in a vector-space basis. The canonical example is a point on a manifold.

Let $M$ be a smooth manifold. A point is an element $p in M$, not a vector in a tangent space $T_p M$. Accordingly, point data has no associated basis: it is represented by chart coordinates, and under a change of chart those coordinates transform by the ordinary chart transition map, not by a change of basis.

Examples

Construct the no-basis object directly:

>>> import coordinax.representations as cxr
>>> basis = cxr.NoBasis()

Use it inside a full representation for point data:

>>> rep = cxr.Representation(cxr.point_geom, basis, cxr.loc)

The representation can then be used with coordinax.representations.cconvert to convert point data between charts while preserving the fact that point data has no basis-dependent components:

>>> import coordinax.charts as cxc
>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> q = cxr.cconvert(p, cxc.cart3d, rep, cxc.sph3d)

The output q is still point data, but expressed in the target chart.

>>> import wadler_lindig as wl
>>> wl.pprint(basis, canonical=False)
NoBasis()
>>> wl.pprint(basis, canonical=True)
no_basis

Notes

NoBasis does not mean “no coordinates”. It means that the represented object is not described by components in a basis-dependent linear space.

canonical_name: ClassVar = 'no_basis'#

Canonical name for the no-basis kind.

class coordinax.representations.AbstractLinearBasis#

Bases: AbstractBasis

Abstract base class for linear (basis-dependent) basis kinds.

A linear basis kind indicates that component data is expressed as components with respect to a specific choice of basis vectors in a tangent space or other linear space.

In the representation model used by coordinax, AbstractLinearBasis refines AbstractBasis by indicating that the object lives in a basis-dependent linear space, in contrast to NoBasis which is used for affine point data.

Examples

>>> import coordinax.representations as cxr
>>> isinstance(cxr.CoordinateBasis(), cxr.AbstractLinearBasis)
True
canonical_name: ClassVar[str | None] = None#

Canonical name for the basis kind.

final class coordinax.representations.CoordinateBasis#

Bases: AbstractLinearBasis

Coordinate basis kind.

A coordinate basis kind indicates that tangent-vector components are expressed in the coordinate basis induced by the current chart. In a coordinate basis, the basis vectors are the partial derivative operators $partial/partial x^i$ associated with the chart coordinates $x^i$.

Mathematical Definition:

Let $(U, varphi)$ be a chart on a smooth manifold $M$. The coordinate basis at a point $p in U$ consists of the tangent vectors $partial_i = partial/partial x^i|_p$. A tangent vector $v in T_p M$ expressed in the coordinate basis has components $v^i$ such that $v = v^i partial_i$.

Under a change of chart, coordinate-basis components transform by the Jacobian matrix of the chart transition map.

Examples

Construct the coordinate basis object directly:

>>> import coordinax.representations as cxr
>>> basis = cxr.CoordinateBasis()

Use it inside a full representation for tangent data:

>>> rep = cxr.Representation(cxr.tangent_geom, basis, cxr.dpl)

Notes

CoordinateBasis does not carry unit-length information. Components in a coordinate basis are not necessarily dimensionless: for example, the $partial_r$ component of a vector in spherical coordinates carries units of inverse length relative to the physical basis.

canonical_name: ClassVar = 'coord_basis'#

Canonical name for the coordinate basis kind.

final class coordinax.representations.PhysicalBasis#

Bases: AbstractLinearBasis

Physical (orthonormal) basis kind.

A physical basis kind indicates that tangent-vector components are expressed in an orthonormal physical basis. In a physical basis, the basis vectors are unit-length and mutually orthogonal with respect to the metric.

Mathematical Definition:

Let $M$ carry a Riemannian metric $g$ and let $(U, varphi)$ be a chart on $M$. The physical basis at a point $p in U$ consists of unit vectors $hat{e}_i = partial_i / sqrt{g_{ii}}$ (for orthogonal charts). A tangent vector $v in T_p M$ expressed in the physical basis has components $hat{v}^i$ such that $v = hat{v}^i hat{e}_i$.

Under a change of chart, physical-basis components transform by both the Jacobian and the normalization factors of the chart transition.

Examples

Construct the physical basis object directly:

>>> import coordinax.representations as cxr
>>> basis = cxr.PhysicalBasis()

Use it inside a full representation for tangent data:

>>> rep = cxr.Representation(cxr.tangent_geom, basis, cxr.vel)

Notes

PhysicalBasis components have consistent physical dimensions across all charts, unlike coordinate-basis components.

canonical_name: ClassVar = 'phys_basis'#

Canonical name for the physical basis kind.

class coordinax.representations.AbstractSemanticKind#

Bases: object

Abstract base class for semantic kind.

A semantic kind specifies the meaning attached to a represented geometric object, independent of the underlying geometric type, basis, and chart used to express its components.

In the representation model used by coordinax, the full representation is determined by three orthogonal pieces of information:

  1. the geometric kind,

  2. the basis, and

  3. the semantic kind.

This class provides the third of these pieces: it answers the question “what does this geometric object represent?” Examples include location for points, and later may include displacement, velocity, or acceleration for tangent-like objects.

Mathematical Role:

The semantic kind refines the interpretation of data within a fixed geometric type and basis.

  • For a point with NoBasis, the semantic kind Location indicates that the data represents where a point lies on a manifold.

  • For a tangent vector, different semantic kinds may distinguish displacement from velocity or acceleration, even when the underlying transformation law is the same.

  • For a cotangent object, semantic kinds may distinguish different dual interpretations without changing the underlying covector character.

Thus semantic kind is distinct from:

  • a chart, which specifies how local coordinates are assigned,

  • a geometry kind, which specifies what sort of geometric object the data represents, and

  • a basis, which specifies in what basis the components are written when such a choice is meaningful.

Examples

>>> import coordinax.main as cx

Construct the location semantic object directly:

>>> semantic = cx.Location()

With the semantic object, we construct a full representation for point data:

>>> rep = cx.Representation(cx.point_geom, cx.no_basis, semantic)

The representation can then be used with cconvert to convert point data between charts while preserving the fact that the data represents a location:

>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> cx.cconvert(p, cx.cart3d, rep, cx.sph3d, rep)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output is still point data representing the same location, but expressed in the target chart.

Notes

This is a static dispatch object and carries no runtime numerical data. Concrete subclasses should represent immutable semantic categories.

canonical_name: ClassVar[str | None] = None#

Canonical name for the geometric kind.

abstractmethod classmethod coord_dimensions(chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of the components for this semantic kind.

Parameters:

chart (AbstractChart[Any, Any, Any]) – The chart defining component names and base coordinate dimensions.

Returns:

One entry per component: the physical dimension of that component under this semantic kind, or None if the component is dimensionless.

Return type:

tuple[PhysicalType | None, ...]

final class coordinax.representations.Location#

Bases: AbstractSemanticKind

Location semantic kind.

A location semantic kind indicates that the represented data specifies where a geometric object is, rather than how it is displaced, how fast it moves, or how it accelerates.

Mathematical Definition:

Location is the canonical semantic kind for point data. Let $M$ be a smooth manifold. A point is an element $p in M$, and Location indicates that the represented data should be interpreted as the coordinates of that point in some chart.

The semantic kind Location therefore does not change the underlying geometric transformation law: for point data, coordinates still transform by the ordinary chart transition map. Instead, it records the interpretation of the point-like data as an actual position on the manifold.

Examples

Construct the location semantic object directly:

>>> import coordinax.representations as cxr
>>> semantic = cxr.Location()

Use it inside a full representation for point data:

>>> rep = cxr.Representation(cxr.point_geom, cxr.no_basis, semantic)

The representation can then be used with cconvert to convert point data between charts while preserving the fact that the data represents a location:

>>> import coordinax.charts as cxc
>>> p = {"x": 1.0, "y": 2.0, "z": 3.0}
>>> cxr.cconvert(p, cxc.cart3d, rep, cxc.sph3d, rep)
{'r': Array(3.74165739, dtype=float64, ...),
 'theta': Array(0.64052231, dtype=float64),
 'phi': Array(1.10714872, dtype=float64, ...)}

The output is still point data representing the same location, but expressed in the target chart.

Notes

Location does not by itself imply that the represented object is a point, but in the current coordinax design it is primarily used as the semantic kind paired with PointGeometry.

canonical_name: ClassVar = 'loc'#

Canonical name for the location semantic kind.

classmethod coord_dimensions(cls, chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of each component for a location.

For a location semantic kind, the dimensions are determined solely by the chart: each component’s dimension is taken directly from chart.coord_dimensions.

Examples

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> cxr.Location.coord_dimensions(cxc.cart3d)
(PhysicalType('length'), PhysicalType('length'), PhysicalType('length'))
Parameters:

chart (AbstractChart[Any, Any, Any])

Return type:

tuple[PhysicalType | None, ...]

class coordinax.representations.AbstractTangentSemanticKind#

Bases: AbstractSemanticKind

Abstract base class for tangent-vector semantic kinds.

A tangent semantic kind specifies the meaning of a tangent-vector object, distinguishing objects that are geometrically identical (both elements of a tangent space $T_p M$) but semantically different in how they are used or interpreted.

Examples include:

  • Displacement: a finite difference between two nearby points,

  • Velocity: the rate of change of position with respect to time,

  • Acceleration: the second derivative of position with respect to time.

All share the same coordinate-transformation law (Jacobian pushforward) but differ in physical dimension and physical interpretation.

Examples

>>> import coordinax.representations as cxr
>>> isinstance(cxr.Displacement(), cxr.AbstractTangentSemanticKind)
True
>>> isinstance(cxr.Velocity(), cxr.AbstractTangentSemanticKind)
True
order: ClassVar[int]#

Time-derivative order of the role (e.g. 0=pos, 1=vel, 2=acc, …).

classmethod coord_dimensions(cls, chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of each component for this tangent kind.

Each component’s base dimension is taken from chart.coord_dimensions and scaled by the appropriate power of time according to cls.order: dimension / time^order.

Examples

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.cart3d)]
['length', 'length', 'length']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.cart3d)]
['speed/...', 'speed/...', 'speed/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.cart3d)]
['acceleration', 'acceleration', 'acceleration']
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.sph3d)]
['length', 'angle', 'angle']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.sph3d)]
['speed/...', 'angular frequency/...', 'angular frequency/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.sph3d)]
['acceleration', 'angular acceleration', 'angular acceleration']
Parameters:

chart (AbstractChart[Any, Any, Any])

Return type:

tuple[PhysicalType | None, ...]

derivative()#

Return the semantic kind one step up the time-derivative ladder.

Looks up self.order + 1 in the internal order registry and returns a fresh instance of the registered class. Raises ValueError if no class is registered at that order.

This design is open for extension: registering a new AbstractTangentSemanticKind subclass at order N automatically makes kind_at_N_minus_1.derivative() return an instance of that class.

Raises:

ValueError – If no tangent semantic kind is registered at self.order + 1.

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Displacement().derivative()
vel
>>> cxr.Velocity().derivative()
acc
>>> try:
...     cxr.Acceleration().derivative()   # no Jerk registered yet
... except ValueError as e:
...     print(type(e).__name__)
ValueError
antiderivative()#

Return the semantic kind one step down the time-derivative ladder.

Looks up self.order - 1 in the internal order registry and returns a fresh instance of the registered class. Raises ValueError if no class is registered at that order.

This design is open for extension: registering a new AbstractTangentSemanticKind subclass at order N automatically makes kind_at_N_plus_1.antiderivative() return an instance of that class. For example, once an Absement class is registered at order -1, Displacement().antiderivative() will return Absement().

Raises:

ValueError – If no tangent semantic kind is registered at self.order - 1.

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Acceleration().antiderivative()
vel
>>> cxr.Velocity().antiderivative()
dpl
>>> try:
...     cxr.Displacement().antiderivative()   # no Absement registered yet
... except ValueError as e:
...     print(type(e).__name__)
ValueError
canonical_name: ClassVar[str | None] = None#

Canonical name for the geometric kind.

final class coordinax.representations.Displacement#

Bases: AbstractTangentSemanticKind

Displacement semantic kind.

A displacement semantic kind indicates that the represented tangent data is a spatial displacement — a finite difference between two nearby points on a manifold, expressed as a tangent vector.

Mathematical Definition:

Let $M$ be a smooth manifold and $p, q in M$ two nearby points. In a chart $(U, varphi)$, the displacement from $p$ to $q$ is the vector $Delta x = varphi(q) - varphi(p) in mathbb{R}^n$. As a tangent vector, this element of $T_p M$ transforms by the Jacobian of the chart transition map.

Examples

>>> import coordinax.representations as cxr
>>> semantic = cxr.Displacement()
>>> semantic.canonical_name
'dpl'
canonical_name: ClassVar = 'dpl'#

Canonical name for the displacement semantic kind.

order: ClassVar[int] = 0#

Time-derivative order of the role (0 for displacement).

derivative()#

Return the Velocity semantic kind.

Displacement has time-order 0; its time derivative is always Velocity (order 1). This override avoids a dict lookup and makes the intent explicit.

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Displacement().derivative()
vel
antiderivative()#

Return the semantic kind one step down the time-derivative ladder.

Looks up self.order - 1 in the internal order registry and returns a fresh instance of the registered class. Raises ValueError if no class is registered at that order.

This design is open for extension: registering a new AbstractTangentSemanticKind subclass at order N automatically makes kind_at_N_plus_1.antiderivative() return an instance of that class. For example, once an Absement class is registered at order -1, Displacement().antiderivative() will return Absement().

Raises:

ValueError – If no tangent semantic kind is registered at self.order - 1.

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Acceleration().antiderivative()
vel
>>> cxr.Velocity().antiderivative()
dpl
>>> try:
...     cxr.Displacement().antiderivative()   # no Absement registered yet
... except ValueError as e:
...     print(type(e).__name__)
ValueError
classmethod coord_dimensions(cls, chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of each component for this tangent kind.

Each component’s base dimension is taken from chart.coord_dimensions and scaled by the appropriate power of time according to cls.order: dimension / time^order.

Examples

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.cart3d)]
['length', 'length', 'length']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.cart3d)]
['speed/...', 'speed/...', 'speed/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.cart3d)]
['acceleration', 'acceleration', 'acceleration']
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.sph3d)]
['length', 'angle', 'angle']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.sph3d)]
['speed/...', 'angular frequency/...', 'angular frequency/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.sph3d)]
['acceleration', 'angular acceleration', 'angular acceleration']
Parameters:

chart (AbstractChart[Any, Any, Any])

Return type:

tuple[PhysicalType | None, ...]

final class coordinax.representations.Velocity#

Bases: AbstractTangentSemanticKind

Velocity semantic kind.

A velocity semantic kind indicates that the represented tangent data is a velocity vector — the first time derivative of position on a manifold.

Mathematical Definition:

Let $gamma: mathbb{R} to M$ be a smooth curve on manifold $M$. The velocity at time $t$ is the tangent vector $dot{gamma}(t) in T_{gamma(t)} M$. In a chart, velocity components $dot{x}^i = d(x^i circ gamma)/dt$ transform by the Jacobian of the chart transition map.

Examples

>>> import coordinax.representations as cxr
>>> semantic = cxr.Velocity()
>>> semantic.canonical_name
'vel'
canonical_name: ClassVar = 'vel'#

Canonical name for the velocity semantic kind.

classmethod coord_dimensions(cls, chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of each component for this tangent kind.

Each component’s base dimension is taken from chart.coord_dimensions and scaled by the appropriate power of time according to cls.order: dimension / time^order.

Examples

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.cart3d)]
['length', 'length', 'length']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.cart3d)]
['speed/...', 'speed/...', 'speed/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.cart3d)]
['acceleration', 'acceleration', 'acceleration']
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.sph3d)]
['length', 'angle', 'angle']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.sph3d)]
['speed/...', 'angular frequency/...', 'angular frequency/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.sph3d)]
['acceleration', 'angular acceleration', 'angular acceleration']
Parameters:

chart (AbstractChart[Any, Any, Any])

Return type:

tuple[PhysicalType | None, ...]

order: ClassVar[int] = 1#

Time-derivative order of the role (1 for velocity).

derivative()#

Return the semantic kind for the time derivative of this velocity.

Velocity has time-order 1. Its time derivative is acceleration (time-order 2).

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Velocity().derivative()
acc
antiderivative()#

Return the semantic kind for the time antiderivative of this velocity.

Velocity has time-order 1. Its time antiderivative is displacement (time-order 0).

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Velocity().antiderivative()
dpl
final class coordinax.representations.Acceleration#

Bases: AbstractTangentSemanticKind

Acceleration semantic kind.

An acceleration semantic kind indicates that the represented tangent data is an acceleration vector — the second time derivative of position on a manifold.

Mathematical Definition:

Let $gamma: mathbb{R} to M$ be a smooth curve on manifold $M$. The acceleration at time $t$ is the tangent vector $ddot{gamma}(t) in T_{gamma(t)} M$ (or more precisely, $nabla_{dot{gamma}}dot{gamma}$ in the covariant sense, but in flat space simply the second derivative).

Examples

>>> import coordinax.representations as cxr
>>> semantic = cxr.Acceleration()
>>> semantic.canonical_name
'acc'
classmethod coord_dimensions(cls, chart: AbstractChart[Any, Any, Any], /)#

Return the physical dimensions of each component for this tangent kind.

Each component’s base dimension is taken from chart.coord_dimensions and scaled by the appropriate power of time according to cls.order: dimension / time^order.

Examples

>>> import coordinax.charts as cxc
>>> import coordinax.representations as cxr
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.cart3d)]
['length', 'length', 'length']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.cart3d)]
['speed/...', 'speed/...', 'speed/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.cart3d)]
['acceleration', 'acceleration', 'acceleration']
>>> [str(x) for x in cxr.Displacement.coord_dimensions(cxc.sph3d)]
['length', 'angle', 'angle']
>>> [str(x) for x in cxr.Velocity.coord_dimensions(cxc.sph3d)]
['speed/...', 'angular frequency/...', 'angular frequency/...']
>>> [str(x) for x in cxr.Acceleration.coord_dimensions(cxc.sph3d)]
['acceleration', 'angular acceleration', 'angular acceleration']
Parameters:

chart (AbstractChart[Any, Any, Any])

Return type:

tuple[PhysicalType | None, ...]

derivative()#

Return the semantic kind one step up the time-derivative ladder.

Looks up self.order + 1 in the internal order registry and returns a fresh instance of the registered class. Raises ValueError if no class is registered at that order.

This design is open for extension: registering a new AbstractTangentSemanticKind subclass at order N automatically makes kind_at_N_minus_1.derivative() return an instance of that class.

Raises:

ValueError – If no tangent semantic kind is registered at self.order + 1.

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Displacement().derivative()
vel
>>> cxr.Velocity().derivative()
acc
>>> try:
...     cxr.Acceleration().derivative()   # no Jerk registered yet
... except ValueError as e:
...     print(type(e).__name__)
ValueError
canonical_name: ClassVar = 'acc'#

Canonical name for the acceleration semantic kind.

order: ClassVar[int] = 2#

Time-derivative order of the role (2 for acceleration).

antiderivative()#

Return the semantic kind for the time antiderivative of this acceleration.

Acceleration has time-order 2. Its time antiderivative is velocity (time-order 1).

Return type:

AbstractTangentSemanticKind

Examples

>>> import coordinax.representations as cxr
>>> cxr.Acceleration().antiderivative()
vel