Coordinate Systems

CoordinateSystem

The base class for coordinate system handling on meshes.

class underworld3.coordinates.CoordinateSystem[source]

Bases: object

This class is attached to a mesh to provide programmatic access to coordinate system operations.

CoordinateSystem.R - the coordinates in the natural reference (sympy symbols) CoordinateSystem.xR - the coordinates in the natural reference (sympy symbols in terms of mesh.X) CoordinateSystem._Rot - the matrix to rotate from X to R, written symbolically CoordinateSystem._xRot - the matrix to rotate from X to R, written in terms of mesh.X

Need these:

Coordinates N -> Native, xi_1, xi_2 by default, but over-written with meaningful symbols if possible Coordinates R -> Natural XYZ or R Theta Phi or R Theta z Coordinates X -> Cartesian XYZ

nRotr - Matrix to map N to R rRotn - nRotr.T nRotx - Matrix to map N to X xRotn - nRotx.T xRotr = xRotn * nRotr rRotx = rRotn * nRotx

Either R or X will be an alias for N depending on whether the DM is Cartesian (Alias has the same vectors but different names)

__init__(mesh, system=CoordinateSystemType.CARTESIAN)[source]
Parameters:

system (CoordinateSystemType | None)

property coords

Mesh vertex coordinates as a numerical array.

Returns the mesh node coordinates in physical units. When the mesh has reference quantities set, coordinates are scaled from internal model units to physical units.

This is the primary interface for accessing coordinate data for initialization, visualization, and coordinate-based calculations.

Returns:

Coordinate array with shape (N, dim) where N is the number of mesh vertices. Returns UnitAwareArray if mesh has units configured.

Return type:

numpy.ndarray or UnitAwareArray

Examples

>>> # Get all coordinates
>>> coords = mesh.X.coords
>>> x_values = coords[:, 0]
>>> y_values = coords[:, 1]
>>> # Use for field initialization
>>> temperature.array[:, 0, 0] = 300 + 100 * mesh.X.coords[:, 0]
>>> # Pass to evaluation functions
>>> values = uw.function.evaluate(expression, mesh.X.coords)
>>> # Get bounds
>>> x_min, x_max = mesh.X.coords[:, 0].min(), mesh.X.coords[:, 0].max()

See also

units

Coordinate units (e.g., ‘kilometer’).

mesh.points

Alias for mesh.X.coords.

var.coords

Variable-specific DOF coordinates (may differ from mesh vertices).

property units

Physical units of the coordinate system.

Returns the units for coordinates accessed via .coords. This is derived from the mesh’s reference quantities and indicates the physical scale of the coordinate data.

Returns:

Unit string (e.g., ‘kilometer’, ‘meter’), or None if no reference quantities are set.

Return type:

str or None

Examples

>>> print(mesh.X.units)  # 'kilometer'
>>> coords = mesh.X.coords  # Coordinates in kilometers

See also

coords

Coordinate data array.

mesh.units

Same as mesh.X.units.

property with_units

Coordinate symbols with unit information.

DEPRECATED (2025-11-26): Following the Transparent Container Principle, this now returns raw coordinate symbols. Units are derived on demand via uw.get_units() which finds the _units attribute on coordinate atoms.

Examples:
>>> x, y = mesh.X.with_units  # Same as mesh.X[0], mesh.X[1]
>>> area = x * y  # Raw SymPy Mul; uw.get_units(area) → km**2
>>> uw.get_units(x)  # → kilometer (derived from _units attribute)
Returns:

tuple: Coordinate symbols (x, y) or (x, y, z)

property shape

Shape of the symbolic coordinate matrix.

property is_symbol

SymPy type check - CoordinateSystem contains symbols but is not itself a symbol.

property is_Matrix

SymPy type check - CoordinateSystem behaves like a Matrix.

property is_scalar

SymPy type check - CoordinateSystem is a Matrix, not a scalar.

property is_number

SymPy type check - CoordinateSystem is not a number.

property is_commutative

SymPy type check - delegate to underlying matrix.

property X: MutableDenseMatrix

Cartesian coordinate system providing symbolic and numerical access.

This is the primary interface for mesh coordinates, providing both:

  • Symbolic access: mesh.X[0], mesh.X[1] for use in equations

  • Numerical access: mesh.X.coords for coordinate data arrays

Returns:

Symbolic coordinate matrix that also supports .coords for data.

Return type:

sympy.Matrix

Examples

>>> # Symbolic access for equations
>>> x, y = mesh.X
>>> temperature_expr = 300 + 100 * x
>>> # Component access
>>> x_coord = mesh.X[0]  # Symbolic x-coordinate
>>> y_coord = mesh.X[1]  # Symbolic y-coordinate
>>> # Numerical coordinate data
>>> coords = mesh.X.coords  # (N, dim) array of vertex positions
>>> x_values = mesh.X.coords[:, 0]
>>> # Coordinate units (if set)
>>> units = mesh.X.units  # e.g., 'kilometer'

See also

R

Natural (curvilinear) coordinates for spherical/cylindrical meshes.

geo

Geographic coordinates for regional spherical meshes.

N

Raw BaseScalar coordinates for JIT compilation.

property x: MutableDenseMatrix

Alias for natural coordinates (lowercase convention).

property N: MutableDenseMatrix

Raw BaseScalar coordinates for derivatives and JIT compilation.

property R: MutableDenseMatrix

Natural (curvilinear) coordinates \((r, \theta, \phi)\) or geographic.

property r: MutableDenseMatrix

Symbolic natural coordinate variables (for mathematical expressions).

property xR: MutableDenseMatrix

Alias for R (backward compatibility).

property geo

Geographic coordinates for GEOGRAPHIC meshes.

Provides access to longitude, latitude, depth coordinates and geographic basis vectors on ellipsoidal (WGS84) meshes.

Returns:

Object with .lon, .lat, .depth, .coords, unit vectors, etc. Use .view() for a complete summary of available properties.

Return type:

GeographicCoordinateAccessor

Raises:

AttributeError – If coordinate system is not GEOGRAPHIC. The error message indicates what coordinate system IS available for this mesh.

Examples

>>> lon = mesh.X.geo.lon         # Longitude data array
>>> geo_coords = mesh.X.geo.coords  # (N, 3) array [lon, lat, depth]
>>> mesh.X.geo.view()            # Show all available properties
property spherical

Spherical/polar coordinates for SPHERICAL and CYLINDRICAL2D meshes.

Provides access to radius and angle coordinates, plus basis vectors: - 3D (SPHERICAL): r, θ (colatitude), φ (longitude) - 2D (CYLINDRICAL2D/Annulus): r, θ (polar angle)

Returns:

Object with .r, .theta, .coords, unit vectors, etc. For 3D also .phi. Use .view() for a complete summary.

Return type:

SphericalCoordinateAccessor

Raises:

AttributeError – If coordinate system is not SPHERICAL or CYLINDRICAL2D. The error message indicates what IS available for this mesh.

Examples

>>> r = mesh.X.spherical.r           # Radius data array
>>> sph_coords = mesh.X.spherical.coords  # (N, 2) or (N, 3) array
>>> mesh.X.spherical.view()          # Show all available properties
property rRotN: MutableDenseMatrix

Rotation matrix from Cartesian (N) to natural (R) coordinates.

Transforms vectors from Cartesian basis to the local curvilinear basis. For spherical: rows are \((\hat{r}, \hat{\theta}, \hat{\phi})\).

property xRotN: MutableDenseMatrix

Rotation matrix from Cartesian (N) to Cartesian (X) - typically identity.

property geoRotN: MutableDenseMatrix

Ellipsoidal rotation matrix for GEOGRAPHIC coordinate systems.

Transforms Cartesian vectors to the local geographic frame: - Row 0: geodetic up (perpendicular to ellipsoid surface) - Row 1: north (meridional, along ellipsoid surface) - Row 2: east (azimuthal, along ellipsoid surface)

This is the ellipsoidal equivalent of rRotN for spherical coordinates. For a sphere (a=b), this reduces to the spherical rotation matrix.

For an ellipsoid, the geodetic normal differs from the radial direction by up to ~10 arcminutes at mid-latitudes, which is significant for regional models at scales of 10-100 km.

Returns:

3×3 rotation matrix, or None if not GEOGRAPHIC coordinate system.

Return type:

sympy.Matrix

Examples

Transform a Cartesian velocity to geographic components:

>>> v_cartesian = sympy.Matrix([[vx, vy, vz]])
>>> v_geo = mesh.CoordinateSystem.geoRotN * v_cartesian.T
>>> v_up, v_north, v_east = v_geo[0], v_geo[1], v_geo[2]
property unit_e_0: MutableDenseMatrix

First natural basis vector (radial \(\hat{r}\) for curvilinear).

property unit_e_1: MutableDenseMatrix

Second natural basis vector (\(\hat{\theta}\) for curvilinear).

property unit_e_2: MutableDenseMatrix

Third natural basis vector (\(\hat{\phi}\) for spherical, None for 2D).

property unit_i: MutableDenseMatrix

Cartesian unit vector \(\hat{i}\) (x-direction).

property unit_j: MutableDenseMatrix

Cartesian unit vector \(\hat{j}\) (y-direction).

property unit_k: MutableDenseMatrix

Cartesian unit vector \(\hat{k}\) (z-direction, None for 2D).

unit_ijk(dirn)[source]

Return Cartesian unit vector for direction index (0=i, 1=j, 2=k).

Return type:

MutableDenseMatrix

property unit_vertical: MutableDenseMatrix

Primary vertical direction for this coordinate system

property unit_horizontal: MutableDenseMatrix

Primary horizontal direction for this coordinate system

property unit_horizontal_0: MutableDenseMatrix

First horizontal direction (alias for unit_horizontal)

property unit_horizontal_1: MutableDenseMatrix

Second horizontal direction (for 3D systems)

property unit_radial: MutableDenseMatrix

Radial direction (for cylindrical/spherical coordinate systems)

property unit_tangential: MutableDenseMatrix

Tangential direction (for cylindrical coordinate systems)

property unit_meridional: MutableDenseMatrix

Meridional direction (for spherical coordinate systems)

property unit_azimuthal: MutableDenseMatrix

Azimuthal direction (for spherical coordinate systems)

property geometric_dimension_names: list

Names of geometric dimensions for this coordinate system

property primary_directions: dict

Dictionary of all available geometric directions for this mesh type

create_line_sample(start_point, direction_vector, length, num_points=50)[source]

Create sample points along a line defined by sympy expressions.

Parameters:
  • start_point (list or numpy.ndarray) – Starting point coordinates in Cartesian space

  • direction_vector (sympy.Matrix) – Direction vector (should be unit vector for accurate length)

  • length (float) – Length of the line to sample

  • num_points (int, optional) – Number of sample points to generate

Returns:

Dictionary containing: - ‘cartesian_coords’: numpy array of Cartesian coordinates for global_evaluate() - ‘natural_coords’: numpy array of natural coordinates for plotting - ‘parameters’: numpy array of parameter values along the line (0 to length)

Return type:

dict

create_profile_sample(profile_type, **params)[source]

Create sample points for common profile types in this coordinate system.

Parameters:
  • profile_type (str) – Type of profile to create. Options depend on coordinate system: - Cartesian: ‘horizontal’, ‘vertical’, ‘diagonal’ - Cylindrical: ‘radial’, ‘tangential’, ‘vertical’ - Spherical: ‘radial’, ‘meridional’, ‘azimuthal’

  • **params – Profile-specific parameters (see individual profile documentation)

Returns:

Dictionary containing: - ‘cartesian_coords’: numpy array of Cartesian coordinates for global_evaluate() - ‘natural_coords’: numpy array of natural coordinates for plotting - ‘parameters’: numpy array of parameter values along the profile

Return type:

dict

zero_matrix(shape)[source]

Matrix of spatial coordinates equivalent to zeros (but still dependent on X) - Add this when you have a matrix with a mix of constants and functions - sympy / numpy can become upset if the constants are not specific functions too.

CoordinateSystemType

Enumeration of supported coordinate system types.

class underworld3.coordinates.CoordinateSystemType[source]

Bases: Enum

Coordinate system types for mesh geometry.

Meshes can have natural coordinate systems that overlay the Cartesian coordinate system used internally for solver assembly. The coordinate system type determines how vector calculus operators (gradient, divergence, curl) are computed.

See also

underworld3.maths.vector_calculus

Operators for each coordinate system.

CARTESIAN = 0
CYLINDRICAL2D = 10
POLAR = 10
CYLINDRICAL3D = 100
SPHERICAL = 200
GEOGRAPHIC = 300