Skip to content
Open
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
99 changes: 95 additions & 4 deletions README.md
Original file line number Diff line number Diff line change
Expand Up @@ -428,6 +428,67 @@ points, cells = om.generate_multiscale_mesh([sdf1, sdf2], [el1, el2])

Global meshes are defined in WGS84 but meshed in a stereographic projection. Regional refinement can be added as additional domains.

#### Understanding Stereographic Projections for Global Meshes

Global ocean meshes require special handling due to the distortion
inherent in projecting a sphere onto a 2D plane. OceanMesh uses a
**north-polar stereographic projection** for global meshing, which:

- Preserves angles (conformal projection)
- Minimises distortion near the poles
- Allows efficient 2D mesh generation algorithms to work on projected
coordinates

**Cartopy Integration (Optional Dependency)**

For accurate stereographic projections, OceanMesh can use
[cartopy](https://scitools.org.uk/cartopy/docs/latest/)'s
``NorthPolarStereo`` coordinate reference system. Install cartopy for
global meshing via the ``global`` extra:

```bash
pip install oceanmesh[global]
```

or manually:

```bash
pip install cartopy
```

If cartopy is not installed, OceanMesh falls back to legacy analytic
formulas (less accurate but functional).

**Scale Factor (k0) Configuration**

The stereographic projection uses a **scale factor** :math:`k_0` that
controls distortion at different latitudes. The local scale factor at
latitude :math:`\phi` is

.. math::

k(\phi) = \frac{2 k_0}{1 + \sin \phi}.

Common :math:`k_0` values:

- :math:`k_0 = 1.0` (default): standard stereographic with true scale
at the pole.
- :math:`k_0 = 0.994`: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).

The scale factor affects mesh sizing: smaller :math:`k_0` reduces
element sizes near the poles, larger :math:`k_0` increases them.
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The math notation mixes inline LaTeX (:math:) with block-level reStructuredText directive (.. math::). In a Markdown file, these reStructuredText constructs won't render properly. Consider using standard Markdown math notation with $ for inline math and $$ for display math, or moving this documentation to a .rst file if Sphinx/reStructuredText formatting is required.

Suggested change
The stereographic projection uses a **scale factor** :math:`k_0` that
controls distortion at different latitudes. The local scale factor at
latitude :math:`\phi` is
.. math::
k(\phi) = \frac{2 k_0}{1 + \sin \phi}.
Common :math:`k_0` values:
- :math:`k_0 = 1.0` (default): standard stereographic with true scale
at the pole.
- :math:`k_0 = 0.994`: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller :math:`k_0` reduces
element sizes near the poles, larger :math:`k_0` increases them.
The stereographic projection uses a **scale factor** $k_0$ that
controls distortion at different latitudes. The local scale factor at
latitude $\phi$ is
$$
k(\phi) = \frac{2 k_0}{1 + \sin \phi}
$$
Common $k_0$ values:
- $k_0 = 1.0$ (default): standard stereographic with true scale
at the pole.
- $k_0 = 0.994$: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller $k_0$ reduces
element sizes near the poles, larger $k_0$ increases them.

Copilot uses AI. Check for mistakes.
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The reStructuredText-style inline math notation (:math:) won't render properly in a Markdown file. Use standard Markdown math notation with $ delimiters instead, e.g., $k_0$.

Suggested change
The stereographic projection uses a **scale factor** :math:`k_0` that
controls distortion at different latitudes. The local scale factor at
latitude :math:`\phi` is
.. math::
k(\phi) = \frac{2 k_0}{1 + \sin \phi}.
Common :math:`k_0` values:
- :math:`k_0 = 1.0` (default): standard stereographic with true scale
at the pole.
- :math:`k_0 = 0.994`: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller :math:`k_0` reduces
element sizes near the poles, larger :math:`k_0` increases them.
The stereographic projection uses a **scale factor** $k_0$ that
controls distortion at different latitudes. The local scale factor at
latitude $\\phi$ is
.. math::
k(\phi) = \frac{2 k_0}{1 + \sin \phi}.
Common $k_0$ values:
- $k_0 = 1.0$ (default): standard stereographic with true scale
at the pole.
- $k_0 = 0.994$: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller $k_0$ reduces
element sizes near the poles, larger $k_0$ increases them.

Copilot uses AI. Check for mistakes.
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The reStructuredText-style inline math notation (:math:) won't render properly in a Markdown file. Use standard Markdown math notation with $ delimiters instead, e.g., $k_0 = 1.0$ for inline math.

Suggested change
The stereographic projection uses a **scale factor** :math:`k_0` that
controls distortion at different latitudes. The local scale factor at
latitude :math:`\phi` is
.. math::
k(\phi) = \frac{2 k_0}{1 + \sin \phi}.
Common :math:`k_0` values:
- :math:`k_0 = 1.0` (default): standard stereographic with true scale
at the pole.
- :math:`k_0 = 0.994`: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller :math:`k_0` reduces
element sizes near the poles, larger :math:`k_0` increases them.
The stereographic projection uses a **scale factor** $k_0$ that
controls distortion at different latitudes. The local scale factor at
latitude $\phi$ is
$$k(\phi) = \frac{2 k_0}{1 + \sin \phi}.$$
Common $k_0$ values:
- $k_0 = 1.0$ (default): standard stereographic with true scale
at the pole.
- $k_0 = 0.994$: used in EPSG:3413 (Arctic) and EPSG:3031
(Antarctic) to minimise distortion at standard parallels (~70°N/S).
The scale factor affects mesh sizing: smaller $k_0$ reduces
element sizes near the poles, larger $k_0$ increases them.

Copilot uses AI. Check for mistakes.

References
~~~~~~~~~~

- [GitHub PR #87](https://github.com/CHLNDDEV/oceanmesh/pull/87):
discussion of scale-factor improvements.
- [Seamsh stereographic example](https://git.immc.ucl.ac.be/jlambrechts/seamsh/-/blob/ca380b59fe4d2ea57ffbf08a7b0c70bdf7df1afb/examples/6-stereographics.py#L55):
similar implementation in seamsh.
- Snyder, J.P. (1987). *Map Projections - A Working Manual*. USGS
Professional Paper 1395.

<!--pytest-codeblocks:skip-->

```python
Expand All @@ -442,7 +503,17 @@ fname_global_latlon = "tests/global/global_latlon.shp"
fname_global_stereo = "tests/global/global_stereo.shp"

global_region = om.Region(extent=(-180.0, 180.0, -89.0, 90.0), crs=4326)
shoreline_global_latlon = om.Shoreline(fname_global_latlon, global_region, 1.0)

# Use k0 = 0.994 for reduced polar distortion (similar to
# EPSG:3413/3031 in polar stereographic grids).
k0 = 0.994

shoreline_global_latlon = om.Shoreline(
fname_global_latlon,
global_region,
1.0,
scale_factor=k0,
)
sdf_global_latlon = om.signed_distance_function(shoreline_global_latlon)
edge_global = om.enforce_mesh_gradation(
om.compute_minimum([
Expand All @@ -451,6 +522,7 @@ edge_global = om.enforce_mesh_gradation(
]),
gradation=0.15,
stereo=True,
scale_factor=k0,
)

aus_region = om.Region(extent=(110.0, 160.0, -45.0, -10.0), crs=4326)
Expand All @@ -464,7 +536,13 @@ edge_regional = om.enforce_mesh_gradation(
gradation=0.12,
)

shoreline_global_stereo = om.Shoreline(fname_global_stereo, global_region, 1.0, stereo=True)
shoreline_global_stereo = om.Shoreline(
fname_global_stereo,
global_region,
1.0,
stereo=True,
scale_factor=k0,
)
sdf_global_stereo = om.signed_distance_function(shoreline_global_stereo)

points, cells = om.generate_multiscale_mesh(
Expand All @@ -485,13 +563,26 @@ points, cells = om.generate_multiscale_mesh(
# Global mesh generation only (stereographic meshing)
import oceanmesh as om
from oceanmesh.region import to_lat_lon

fname = "tests/global/global_latlon.shp"
fname2 = "tests/global/global_stereo.shp"
region = om.Region(extent=(-180.00, 180.00, -89.00, 90.00), crs=4326)
shore = om.Shoreline(fname, region.bbox, 0.5)

# Use k0 = 0.994 to mirror common polar stereographic practice.
k0 = 0.994

shore = om.Shoreline(fname, region.bbox, 0.5, scale_factor=k0)
edge = om.distance_sizing_function(shore, rate=0.11)
domain = om.signed_distance_function(om.Shoreline(fname2, region.bbox, 0.5, stereo=True))
edge = om.enforce_mesh_gradation(edge, gradation=0.15, stereo=True, scale_factor=k0)

domain = om.signed_distance_function(
om.Shoreline(fname2, region.bbox, 0.5, stereo=True, scale_factor=k0)
)
points, cells = om.generate_mesh(domain, edge, stereo=True, max_iter=100)

# Points are in stereographic coordinates; convert back to lon/lat for
# post-processing or export.
lon, lat = to_lat_lon(points[:, 0], points[:, 1])
```

[Back to top](#table-of-contents)
Expand Down
9 changes: 7 additions & 2 deletions oceanmesh/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -126,6 +126,11 @@

__version__ = _version.get_versions()["version"]

from . import _version
try: # Optional global-stereo helpers (import may fail generically)
from .projections import CARTOPY_AVAILABLE
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

StereoProjection is added to __all__ but never imported into the module namespace. This will cause an AttributeError when users try to access oceanmesh.StereoProjection. Add from .projections import StereoProjection before extending __all__.

Suggested change
from .projections import CARTOPY_AVAILABLE
from .projections import CARTOPY_AVAILABLE, StereoProjection

Copilot uses AI. Check for mistakes.

__version__ = _version.get_versions()["version"]
__all__.extend(["StereoProjection", "CARTOPY_AVAILABLE"])
except ImportError:
# Projections module not importable; expose CARTOPY_AVAILABLE flag as False
CARTOPY_AVAILABLE = False
__all__.append("CARTOPY_AVAILABLE")
88 changes: 72 additions & 16 deletions oceanmesh/edgefx.py
Original file line number Diff line number Diff line change
Expand Up @@ -88,7 +88,9 @@ def enforce_mesh_size_bounds_elevation(grid, dem, bounds):
return grid


def enforce_mesh_gradation(grid, gradation=0.15, crs="EPSG:4326", stereo=False):
def enforce_mesh_gradation(
grid, gradation=0.15, crs="EPSG:4326", stereo=False, scale_factor=1.0
):
"""Enforce a mesh size gradation bound `gradation` on a :class:`grid`

Parameters
Expand Down Expand Up @@ -124,28 +126,78 @@ def enforce_mesh_gradation(grid, gradation=0.15, crs="EPSG:4326", stereo=False):
tmp = gradient_limit([*sz], elen, gradation, 10000, cell_size)
tmp = np.reshape(tmp, (sz[0], sz[1]), "F")
if stereo:
logger.info("Global mesh: fixing gradient on the north pole...")
# max distortion at the pole: 2 / 180 * PI / (1 - cos(lat))**2
dx_stereo = grid.dx * 1 / 180 * np.pi / 2
# in stereo projection, all north hemisphere is contained in the unit sphere
# we want to fix the gradient close to the north pole,
# so we extract all the coordinates between -1 and 1 in stereographic projection
from .projections import CARTOPY_AVAILABLE
from .region import _get_stereo_projection

backend = "cartopy" if CARTOPY_AVAILABLE else "hardcoded"
logger.info(
"Global mesh: fixing gradient on the north pole (backend=%s, k0=%s, distortion_correction=%s)...",
backend,
scale_factor,
"enabled" if CARTOPY_AVAILABLE else "disabled",
)

proj = _get_stereo_projection(scale_factor=scale_factor)

# Choose a reasonable physical radius (in metres) around the
# north pole within which to regularise the gradient. Here we
# use ~2000 km.
r_max = 2_000_000.0
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The polar radius value 2_000_000.0 (~2000 km) is hardcoded. Consider defining this as a configurable parameter or a named constant (e.g., POLAR_REGULARIZATION_RADIUS_M = 2_000_000.0) to improve maintainability and allow users to adjust the polar regularization extent if needed for different use cases.

Copilot uses AI. Check for mistakes.

# Resolution in stereographic metres, based on underlying lon/lat step.
dx_deg = grid.dx
dx_stereo = dx_deg * np.pi / 180.0 * 6_371_000.0
Copy link

Copilot AI Dec 8, 2025

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The Earth radius value 6_371_000.0 is hardcoded here and also appears on line 170. Consider defining this as a module-level or function-level constant (e.g., EARTH_RADIUS_M = 6_371_000.0) to avoid magic numbers and ensure consistency. This would improve maintainability if the value needs to be adjusted or if different ellipsoid models are used in the future.

Copilot uses AI. Check for mistakes.

us, vs = np.meshgrid(
np.arange(-1, 1, dx_stereo), np.arange(-1, 1, dx_stereo), indexing="ij"
np.arange(-r_max, r_max, dx_stereo),
np.arange(-r_max, r_max, dx_stereo),
indexing="ij",
)
ulon, vlat = to_lat_lon(us.ravel(), vs.ravel())

if proj is not None:
ulon, vlat = proj.to_lat_lon(us.ravel(), vs.ravel())
# Spatially-varying stereographic scale factor k(phi) on
# the intermediate polar grid.
k_polar = proj.get_scale_factor(vlat)
else:
# Fall back to module-level transform if projection could
# not be constructed. Normalise metre-based coordinates by
# Earth radius so that to_lat_lon operates on the
# non-dimensional unit-sphere coordinates expected by the
# legacy stereographic formulas used when CARTOPY_AVAILABLE
# is False. This keeps the polar patch consistent with
# to_stereo/to_lat_lon.
R_earth = 6_371_000.0
us_unit = us / R_earth
vs_unit = vs / R_earth
ulon, vlat = to_lat_lon(us_unit.ravel(), vs_unit.ravel())
k_polar = None

utmp = grid.eval((ulon, vlat))
utmp = np.reshape(utmp, us.shape)
szs = utmp.shape

# Distortion-aware gradient limiting: operate on mesh sizes
# scaled by the local stereographic factor so the limiter
# acts in approximately physical space.
if k_polar is not None:
k_polar_grid = np.reshape(k_polar, utmp.shape)
utmp_work = utmp * k_polar_grid
else:
utmp_work = utmp

szs = utmp_work.shape
szs = (szs[0], szs[1], 1)
# we choose an excessively large number for the gradiation = 10
# this is merely to fix the north pole gradient
vtmp = gradient_limit([*szs], dx_stereo, 10, 10000, utmp.flatten("F"))
vtmp = gradient_limit([*szs], dx_stereo, 10, 10000, utmp_work.flatten("F"))
vtmp = np.reshape(vtmp, (szs[0], szs[1]), "F")
# construct stereo interpolating function

# Inverse correction: return to geographic sizing by dividing
# out the stereographic distortion.
if k_polar is not None:
vtmp /= k_polar_grid

# construct stereo interpolating function in projection metres
grid_stereo = Grid(
bbox=(-1, 1, -1, 1),
bbox=(-r_max, r_max, -r_max, r_max),
dx=dx_stereo,
values=vtmp,
hmin=grid.hmin,
Expand All @@ -155,7 +207,11 @@ def enforce_mesh_gradation(grid, gradation=0.15, crs="EPSG:4326", stereo=False):
grid_stereo.build_interpolant()
# reinject back into the original grid and redo the gradient computation
xg, yg = grid.create_grid()
tmp[yg > 0] = grid_stereo.eval(to_stereo(xg[yg > 0], yg[yg > 0]))
if proj is not None:
xs, ys = proj.to_stereo(xg[yg > 0], yg[yg > 0])
else:
xs, ys = to_stereo(xg[yg > 0], yg[yg > 0])
tmp[yg > 0] = grid_stereo.eval((xs, ys))
logger.info(
"Global mesh: reinject back stereographic gradient and recomputing gradient..."
)
Expand Down
Loading
Loading