Skip to content

API Reference

Seafloor Age Tracking

SeafloorAgeTracker

Seafloor age tracker using Lagrangian particle tracking with C++ backend.

Uses pygplates.TopologicalModel.reconstruct_geometry() for efficient point advection with built-in collision detection.

Key features: - GPlately-compatible: Matches GPlately's SeafloorGrid output - C++ backend: Fast reconstruction using pygplates internals - Icosahedral initialization: Full ocean coverage from start - Continental filtering: Via polygon queries with caching - Checkpointing: Save/restore state for restarts

Parameters:

Name Type Description Default
rotation_files list of str

Paths to rotation model files (.rot).

required
topology_files list of str

Paths to topology/plate boundary files (.gpml/.gpmlz).

required
continental_polygons str

Path to continental polygon file. If None, tracers are not removed when they enter continental regions.

None
config TracerConfig

Configuration parameters. If None, uses defaults.

None
verbose bool

Deprecated. Use GTRACK_LOGLEVEL environment variable instead. Set GTRACK_LOGLEVEL=INFO for progress messages or DEBUG for details.

True

Examples:

>>> # Initialize tracker
>>> tracker = SeafloorAgeTracker(
...     rotation_files=['rotations.rot'],
...     topology_files=['topologies.gpmlz'],
...     continental_polygons='continents.gpmlz'
... )
>>>
>>> # Initialize with icosahedral mesh (GPlately-compatible)
>>> tracker.initialize(starting_age=200)
>>>
>>> # Step forward (decreasing geological age toward present)
>>> for target_age in range(199, -1, -1):
...     cloud = tracker.step_to(target_age)
...     xyz = cloud.xyz
...     ages = cloud.get_property('age')

__init__(rotation_files: Union[str, List[str]], topology_files: Union[str, List[str]], continental_polygons: Optional[str] = None, config: Optional[TracerConfig] = None, verbose: bool = True)

initialize(starting_age: float, method: str = 'icosahedral', refinement_levels: Optional[int] = None, initial_ocean_mean_spreading_rate: Optional[float] = None, age_distance_law: Optional[Callable[[np.ndarray, float], np.ndarray]] = None) -> int

Initialize tracers for given geological age.

Parameters:

Name Type Description Default
starting_age float

Starting geological age (Ma). Tracers are placed based on ocean structure at this age.

required
method str

Initialization method: - 'icosahedral': Full ocean mesh with computed ages (GPlately-compatible) - 'ridge_only': Tracers only at ridges with age=0 (legacy gtrack)

'icosahedral'
refinement_levels int

Icosahedral mesh refinement level. If None, uses config default.

None
initial_ocean_mean_spreading_rate float

Spreading rate for age calculation (mm/yr). If None, uses config default.

None
age_distance_law callable

Custom function to convert distance to age. Signature: (distances_km, spreading_rate_mm_yr) -> ages_myr

None

Returns:

Type Description
int

Number of tracers initialized.

Examples:

>>> # GPlately-compatible initialization
>>> tracker.initialize(starting_age=200)
>>>
>>> # Higher resolution
>>> tracker.initialize(starting_age=200, refinement_levels=6)
>>>
>>> # Custom age calculation
>>> def my_age_law(distances, rate):
...     return distances / (rate / 2) * 1.1  # 10% older
>>> tracker.initialize(starting_age=200, age_distance_law=my_age_law)

initialize_from_cloud(cloud: PointCloud, current_age: float) -> int

Initialize from existing PointCloud.

Use this to restart from a checkpoint or to provide custom initial tracer positions.

Parameters:

Name Type Description Default
cloud PointCloud

Point cloud with 'age' property containing material age of each tracer (time since ridge formation).

required
current_age float

Current geological age (Ma).

required

Returns:

Type Description
int

Number of tracers initialized.

Raises:

Type Description
ValueError

If cloud does not have 'age' property.

step_to(target_age: float) -> PointCloud

Evolve tracers to target geological age using C++ backend.

Can only step forward (decreasing geological age toward 0).

Parameters:

Name Type Description Default
target_age float

Target geological age (Ma). Must be less than current_age.

required

Returns:

Type Description
PointCloud

Point cloud with 'age' property containing material ages.

Raises:

Type Description
RuntimeError

If tracker is not initialized.

ValueError

If target_age > current_age (can only go forward).

get_current_state() -> PointCloud

Get current tracers as PointCloud without evolving.

Returns:

Type Description
PointCloud

Current tracer positions with 'age' property.

save_checkpoint(filepath: str) -> None

Save state to checkpoint file.

Parameters:

Name Type Description Default
filepath str

Path to save checkpoint (.npz format).

required

load_checkpoint(filepath: str) -> None

Load state from checkpoint file.

Parameters:

Name Type Description Default
filepath str

Path to checkpoint file.

required

compute_ages(target_age: float, starting_age: float, rotation_files: Union[str, List[str]], topology_files: Union[str, List[str]], continental_polygons: Optional[str] = None, config: Optional[TracerConfig] = None, verbose: bool = True) -> PointCloud classmethod

One-shot computation of seafloor ages (functional interface).

Creates a tracker, initializes at starting_age, and evolves to target_age in a single call.

Parameters:

Name Type Description Default
target_age float

Target geological age (Ma).

required
starting_age float

Starting geological age (Ma).

required
rotation_files list of str

Paths to rotation model files (.rot).

required
topology_files list of str

Paths to topology/plate boundary files (.gpml/.gpmlz).

required
continental_polygons str

Path to continental polygon file.

None
config TracerConfig

Configuration parameters.

None
verbose bool

Print progress information.

True

Returns:

Type Description
PointCloud

Point cloud with 'age' property.

Examples:

>>> cloud = SeafloorAgeTracker.compute_ages(
...     target_age=100,
...     starting_age=200,
...     rotation_files=['rotations.rot'],
...     topology_files=['topologies.gpmlz']
... )
>>> ages = cloud.get_property('age')

TracerConfig dataclass

Configuration parameters for seafloor age tracking.

Parameters match GPlately's SeafloorGrid for compatibility.

Attributes:

Name Type Description
time_step float

Time step size in Myr (default: 1.0)

earth_radius float

Earth's radius in meters (default: 6.3781e6)

Collision Detection (C++ backend)

velocity_delta_threshold : float Minimum velocity difference to trigger collision check, in km/Myr. Converted to cm/yr (divide by 10) for pygplates C++ API. Default: 7.0 km/Myr (= 0.7 cm/yr) distance_threshold_per_myr : float Base distance threshold for collision detection, in km/Myr. Default: 10.0 km/Myr

Initialization

default_refinement_levels : int Number of icosahedral mesh refinement levels. Level 5 = ~10,242 points, Level 6 = ~40,962 points. Default: 5 initial_ocean_mean_spreading_rate : float Mean spreading rate in mm/yr for initial age calculation. Used to compute age = distance / (rate / 2). Default: 75.0 mm/yr (GPlately default)

MOR Seed Generation

ridge_sampling_degrees : float Ridge tessellation resolution in degrees (~50 km at equator). Default: 0.5 degrees spreading_offset_degrees : float Angle to rotate new points off ridge, in degrees. Default: 0.01 degrees (~1 km)

Continental Handling

continental_cache_size : int Number of timesteps to cache for continental polygon queries. Default: 10

Reinitialization

reinit_k_neighbors : int Number of nearest neighbors for interpolation during reinitialization. Default: 5 reinit_max_distance : float Maximum distance (meters) for valid interpolation neighbors. Default: 500e3 (500 km)

Examples:

>>> # Use default configuration (GPlately-compatible)
>>> config = TracerConfig()
>>>
>>> # Custom configuration with higher resolution
>>> config = TracerConfig(
...     default_refinement_levels=6,  # ~40,962 points
...     ridge_sampling_degrees=0.25,  # ~25 km resolution
...     time_step=0.5                 # 0.5 Myr timesteps
... )

velocity_delta_threshold_cm_yr: float property

Velocity threshold in cm/yr for pygplates C++ API.

The C++ API expects cm/yr, while we store km/Myr for readability. Conversion: 1 km/Myr = 0.1 cm/yr

__post_init__()

Validate configuration parameters.

to_dict() -> dict

Convert configuration to dictionary.

Returns:

Type Description
dict

Configuration as dictionary

from_dict(config_dict: dict) classmethod

Create configuration from dictionary.

Parameters:

Name Type Description Default
config_dict dict

Dictionary with configuration parameters

required

Returns:

Type Description
TracerConfig

Configuration object


Point Rotation

PointCloud dataclass

Container for points with associated properties.

Stores points in Cartesian XYZ format internally (matches gadopt). Properties (lithospheric_depth, etc.) are stored separately from positions.

Parameters:

Name Type Description Default
xyz ndarray

Cartesian coordinates, shape (N, 3), in meters. Points should lie on Earth's surface (radius ~6.3781e6 m).

required
properties dict

Dictionary mapping property names to arrays of shape (N,). Properties are preserved during rotation operations.

dict()
plate_ids ndarray

Plate IDs for each point, shape (N,). Required for rotation.

None

Examples:

>>> xyz = np.random.randn(1000, 3)
>>> from gtrack.geometry import normalize_to_sphere
>>> xyz = normalize_to_sphere(xyz)  # Project to Earth's surface
>>> cloud = PointCloud(xyz=xyz)
>>> cloud.add_property('lithospheric_depth', np.random.rand(1000) * 100e3)

xyz: np.ndarray instance-attribute

latlon: np.ndarray property

Get lat/lon coordinates (computed from XYZ).

Returns:

Type Description
ndarray

Array of shape (N, 2) with [lat, lon] in degrees. Latitude: -90 to 90, Longitude: -180 to 180.

lonlat: np.ndarray property

Get lon/lat coordinates (computed from XYZ).

Returns:

Type Description
ndarray

Array of shape (N, 2) with [lon, lat] in degrees. Longitude: -180 to 180, Latitude: -90 to 90.

__init__(xyz: np.ndarray, properties: Dict[str, np.ndarray] = dict(), plate_ids: Optional[np.ndarray] = None) -> None

from_latlon(latlon: np.ndarray, properties: Optional[Dict[str, np.ndarray]] = None) -> PointCloud classmethod

Create PointCloud from lat/lon coordinates.

Parameters:

Name Type Description Default
latlon ndarray

Coordinates, shape (N, 2) with [lat, lon] in degrees.

required
properties dict

Properties to attach to the points.

None

Returns:

Type Description
PointCloud

New PointCloud with XYZ coordinates computed from lat/lon.

Examples:

>>> latlon = np.array([[45.0, -120.0], [30.0, 90.0]])
>>> cloud = PointCloud.from_latlon(latlon)

add_property(name: str, values: np.ndarray) -> None

Add or update a property.

Parameters:

Name Type Description Default
name str

Name of the property.

required
values ndarray

Property values, shape (N,).

required

Raises:

Type Description
ValueError

If values length doesn't match number of points.

get_property(name: str) -> np.ndarray

Get a property by name.

Parameters:

Name Type Description Default
name str

Name of the property.

required

Returns:

Type Description
ndarray

Property values.

Raises:

Type Description
KeyError

If property not found.

remove_property(name: str) -> None

Remove a property.

Parameters:

Name Type Description Default
name str

Name of the property to remove.

required

subset(mask: np.ndarray) -> PointCloud

Create subset of points using boolean mask.

Parameters:

Name Type Description Default
mask ndarray

Boolean mask, shape (N,). True values are kept.

required

Returns:

Type Description
PointCloud

New PointCloud with subset of points.

copy() -> PointCloud

Create a deep copy.

Returns:

Type Description
PointCloud

Deep copy of this PointCloud.

PointRotator

Rotate points between geological ages using plate reconstructions.

This class provides the main API for rotating user-provided points according to plate tectonic reconstructions.

Key Features: - Cartesian XYZ internal representation (matches gadopt) - Properties stored separately from positions and preserved during rotation - Batched rotation for performance (10-50x speedup by grouping by plate_id) - Handles undefined plates with warnings

Parameters:

Name Type Description Default
rotation_files list of str

Paths to rotation model files (.rot).

required
topology_files list of str

Paths to topology/plate boundary files (.gpmlz). Used for partitioning if static_polygons not provided.

None
static_polygons str

Path to static polygons for plate ID assignment. If None, uses topology_files for partitioning.

None

Examples:

>>> rotator = PointRotator(
...     rotation_files=['rotations.rot'],
...     topology_files=['topologies.gpmlz'],
...     static_polygons='static_polygons.gpmlz'
... )
>>>
>>> # Load user points
>>> cloud = PointCloud.from_latlon(my_latlon_array)
>>>
>>> # Assign plate IDs at present day
>>> cloud = rotator.assign_plate_ids(cloud, at_age=0.0)
>>>
>>> # Rotate to 50 Ma
>>> rotated = rotator.rotate(cloud, from_age=0.0, to_age=50.0)

assign_plate_ids(cloud: PointCloud, at_age: float, use_static_polygons: bool = True, remove_undefined: bool = True) -> PointCloud

Assign plate IDs to points based on their positions.

Parameters:

Name Type Description Default
cloud PointCloud

Points to assign plate IDs to.

required
at_age float

Geological age at which to assign plate IDs (Ma). Use 0.0 for present-day positions.

required
use_static_polygons bool

If True and static_polygons are loaded, use them for assignment. Otherwise use topology_features.

True
remove_undefined bool

If True, remove points with undefined plate IDs (plate_id=0) and emit a warning. If False, keep points with plate_id=0.

True

Returns:

Type Description
PointCloud

Cloud with plate_ids assigned. May have fewer points if remove_undefined=True and some points had undefined plates.

Warns:

Type Description
UserWarning

If any points have undefined plate IDs.

rotate(cloud: PointCloud, from_age: float, to_age: float, reassign_plate_ids: bool = False) -> PointCloud

Rotate points from one geological age to another.

Uses batched rotation operations for performance (10-50x speedup by grouping points by plate_id).

Parameters:

Name Type Description Default
cloud PointCloud

Points to rotate. Must have plate_ids assigned.

required
from_age float

Source geological age (Ma).

required
to_age float

Target geological age (Ma).

required
reassign_plate_ids bool

If True, reassign plate IDs at target age. If False, keep original plate IDs.

False

Returns:

Type Description
PointCloud

Rotated points with same properties.

Raises:

Type Description
ValueError

If cloud does not have plate_ids assigned.

Notes

Direction of rotation: - from_age=0, to_age=50: Rotate present-day positions to 50 Ma - from_age=50, to_age=0: Rotate 50 Ma positions to present day

Examples:

>>> # Rotate present-day continental points to 50 Ma
>>> rotated = rotator.rotate(cloud, from_age=0.0, to_age=50.0)

rotate_incremental(cloud: PointCloud, from_age: float, to_age: float, time_step: float = 1.0, reassign_at_each_step: bool = True) -> PointCloud

Rotate points incrementally through geological time.

Useful when plate IDs may change during rotation (e.g., points crossing plate boundaries over long time spans).

Parameters:

Name Type Description Default
cloud PointCloud

Points to rotate. Must have plate_ids assigned.

required
from_age float

Source geological age (Ma).

required
to_age float

Target geological age (Ma).

required
time_step float

Time step for incremental rotation (Myr).

1.0
reassign_at_each_step bool

If True, reassign plate IDs after each time step. This handles points that cross plate boundaries.

True

Returns:

Type Description
PointCloud

Rotated points.

Examples:

>>> # Rotate with 1 Myr steps, reassigning plates
>>> rotated = rotator.rotate_incremental(
...     cloud, from_age=0.0, to_age=100.0, time_step=1.0
... )

PolygonFilter

Filter points by polygon containment.

Supports filtering by: - Continental polygons (keep only continental points) - Custom polygons (user-provided) - Exclusion zones (remove points inside polygons)

Parameters:

Name Type Description Default
polygon_files str or list of str

Path(s) to polygon files (.gpml, .gpmlz).

required
rotation_files str or list of str

Paths to rotation model files (.rot).

required

Examples:

>>> filter = PolygonFilter(
...     polygon_files='continental_polygons.gpmlz',
...     rotation_files=['rotations.rot']
... )
>>>
>>> # Keep only continental points
>>> continental_cloud = filter.filter_inside(cloud, at_age=0.0)
>>>
>>> # Remove continental points (keep oceanic)
>>> oceanic_cloud = filter.filter_outside(cloud, at_age=0.0)

get_containment_mask(cloud: PointCloud, at_age: float) -> np.ndarray

Get boolean mask of points inside polygons.

Parameters:

Name Type Description Default
cloud PointCloud

Points to check.

required
at_age float

Geological age at which to check containment (Ma). Use 0.0 for present-day polygons.

required

Returns:

Type Description
ndarray

Boolean mask, shape (N,), True for points inside polygons.

filter_inside(cloud: PointCloud, at_age: float) -> PointCloud

Keep only points inside polygons.

Parameters:

Name Type Description Default
cloud PointCloud

Points to filter.

required
at_age float

Geological age at which to check containment (Ma).

required

Returns:

Type Description
PointCloud

Points inside polygons.

Examples:

>>> # Keep only continental points at present day
>>> continental = filter.filter_inside(cloud, at_age=0.0)

filter_outside(cloud: PointCloud, at_age: float) -> PointCloud

Keep only points outside polygons.

Parameters:

Name Type Description Default
cloud PointCloud

Points to filter.

required
at_age float

Geological age at which to check containment (Ma).

required

Returns:

Type Description
PointCloud

Points outside polygons.

Examples:

>>> # Remove continental points (keep oceanic) at present day
>>> oceanic = filter.filter_outside(cloud, at_age=0.0)

get_statistics(cloud: PointCloud, at_age: float) -> dict

Get statistics about polygon containment.

Parameters:

Name Type Description Default
cloud PointCloud

Points to analyze.

required
at_age float

Geological age at which to check containment (Ma).

required

Returns:

Type Description
dict

Statistics including: - total: Total number of points - inside: Number of points inside polygons - outside: Number of points outside polygons - inside_fraction: Fraction of points inside


I/O Functions

load_points_numpy(filepath: Union[str, Path], xyz_columns: Tuple[int, int, int] = (0, 1, 2), property_columns: Optional[Dict[str, int]] = None) -> PointCloud

Load points from numpy file.

Parameters:

Name Type Description Default
filepath str or Path

Path to .npy or .npz file.

required
xyz_columns tuple

Column indices for x, y, z coordinates (for .npy files).

(0, 1, 2)
property_columns dict

Mapping from property name to column index (for .npy files).

None

Returns:

Type Description
PointCloud

Loaded points.

Examples:

>>> # Load from .npz with xyz and properties
>>> cloud = load_points_numpy('points.npz')
>>>
>>> # Load from .npy with specific columns
>>> cloud = load_points_numpy(
...     'data.npy',
...     xyz_columns=(0, 1, 2),
...     property_columns={'depth': 3, 'temperature': 4}
... )

load_points_latlon(filepath: Union[str, Path], latlon_columns: Tuple[int, int] = (0, 1), property_columns: Optional[Dict[str, int]] = None, delimiter: str = ',', skip_header: int = 0) -> PointCloud

Load points from lat/lon text file (CSV, etc.).

Parameters:

Name Type Description Default
filepath str or Path

Path to text file.

required
latlon_columns tuple

Column indices for lat, lon (in degrees).

(0, 1)
property_columns dict

Mapping from property name to column index.

None
delimiter str

Column delimiter.

','
skip_header int

Number of header lines to skip.

0

Returns:

Type Description
PointCloud

Loaded points.

Examples:

>>> cloud = load_points_latlon(
...     'points.csv',
...     latlon_columns=(0, 1),
...     property_columns={'depth': 2}
... )

save_points_numpy(cloud: PointCloud, filepath: Union[str, Path], include_properties: bool = True) -> None

Save points to numpy format.

Parameters:

Name Type Description Default
cloud PointCloud

Points to save.

required
filepath str or Path

Output path (.npy or .npz).

required
include_properties bool

If True, save properties (requires .npz format).

True

Examples:

>>> save_points_numpy(cloud, 'output.npz')

save_points_latlon(cloud: PointCloud, filepath: Union[str, Path], delimiter: str = ',', header: Optional[str] = None, include_properties: bool = True) -> None

Save points to lat/lon text file.

Parameters:

Name Type Description Default
cloud PointCloud

Points to save.

required
filepath str or Path

Output path.

required
delimiter str

Column delimiter.

','
header str

Header line to write.

None
include_properties bool

If True, include properties as additional columns.

True

Examples:

>>> save_points_latlon(cloud, 'output.csv', header='lat,lon,depth')

PointCloudCheckpoint

Checkpoint manager for PointCloud state.

Provides save/load functionality with metadata for checkpointing during long-running simulations.

Examples:

>>> checkpoint = PointCloudCheckpoint()
>>>
>>> # Save with metadata
>>> checkpoint.save(cloud, 'checkpoint_50Ma.npz', geological_age=50.0)
>>>
>>> # Load and get metadata
>>> cloud, metadata = checkpoint.load('checkpoint_50Ma.npz')
>>> print(metadata['geological_age'])  # 50.0

save(cloud: PointCloud, filepath: Union[str, Path], geological_age: Optional[float] = None, metadata: Optional[Dict] = None) -> None

Save checkpoint with metadata.

Parameters:

Name Type Description Default
cloud PointCloud

Points to save.

required
filepath str or Path

Output path (.npz format).

required
geological_age float

Current geological age for reference.

None
metadata dict

Additional metadata to save.

None

Examples:

>>> checkpoint.save(
...     cloud, 'state.npz',
...     geological_age=50.0,
...     metadata={'simulation_step': 100}
... )

load(filepath: Union[str, Path]) -> Tuple[PointCloud, Dict]

Load checkpoint.

Parameters:

Name Type Description Default
filepath str or Path

Path to checkpoint file.

required

Returns:

Name Type Description
cloud PointCloud

Loaded point cloud.

metadata dict

Associated metadata.

Examples:

>>> cloud, metadata = checkpoint.load('state.npz')
>>> print(f"Loaded cloud at {metadata['geological_age']} Ma")

list_checkpoints(directory: Union[str, Path], pattern: str = '*.npz') -> list

List checkpoint files in a directory.

Parameters:

Name Type Description Default
directory str or Path

Directory to search.

required
pattern str

Glob pattern for checkpoint files.

"*.npz"

Returns:

Type Description
list

Sorted list of checkpoint file paths.


Mesh Generation

create_icosahedral_mesh(refinement_levels: int = 5) -> pygplates.MultiPointOnSphere

Create an icosahedral mesh as pygplates.MultiPointOnSphere.

This matches GPlately's mesh generation for compatibility. The mesh starts as a 12-vertex icosahedron and is recursively refined by bisecting each face.

Parameters:

Name Type Description Default
refinement_levels int

Number of refinement levels. Higher values create denser meshes. Level 5 produces ~10,242 points. Level 6 produces ~40,962 points.

5

Returns:

Type Description
MultiPointOnSphere

The mesh points as a MultiPointOnSphere object.

Examples:

>>> mesh = create_icosahedral_mesh(refinement_levels=5)
>>> len(list(mesh.get_points()))
10242

create_icosahedral_mesh_latlon(refinement_levels: int = 5) -> Tuple[np.ndarray, np.ndarray]

Create an icosahedral mesh returning latitude and longitude arrays.

Parameters:

Name Type Description Default
refinement_levels int

Number of refinement levels.

5

Returns:

Name Type Description
lats ndarray

Latitudes in degrees, shape (N,)

lons ndarray

Longitudes in degrees, shape (N,)

Examples:

>>> lats, lons = create_icosahedral_mesh_latlon(refinement_levels=5)
>>> len(lats)
10242

create_icosahedral_mesh_xyz(refinement_levels: int = 5, radius: float = 1.0) -> np.ndarray

Create an icosahedral mesh returning XYZ coordinates.

Parameters:

Name Type Description Default
refinement_levels int

Number of refinement levels.

5
radius float

Radius of the sphere (1.0 for unit sphere, or Earth radius in meters).

1.0

Returns:

Name Type Description
xyz ndarray

XYZ coordinates, shape (N, 3)

Examples:

>>> xyz = create_icosahedral_mesh_xyz(refinement_levels=5)
>>> xyz.shape
(10242, 3)

mesh_point_count(refinement_levels: int) -> int

Calculate the number of points in an icosahedral mesh.

The formula is: 10 * 4^levels + 2

Parameters:

Name Type Description Default
refinement_levels int

Number of refinement levels.

required

Returns:

Type Description
int

Number of points in the mesh.

Examples:

>>> mesh_point_count(5)
10242
>>> mesh_point_count(6)
40962

Logging

enable_verbose() -> None

Enable verbose output (INFO level).

Convenience function equivalent to set_log_level(logging.INFO).

enable_debug() -> None

Enable debug output (DEBUG level).

Convenience function equivalent to set_log_level(logging.DEBUG).

disable_logging() -> None

Disable all gtrack logging output.

Convenience function equivalent to set_log_level(logging.CRITICAL + 1).

set_log_level(level: int) -> None

Set the log level for all gtrack loggers.

Parameters:

Name Type Description Default
level int

Logging level (e.g., logging.DEBUG, logging.INFO).

required

Examples:

>>> import logging
>>> from gtrack.logging import set_log_level
>>> set_log_level(logging.DEBUG)  # Enable debug output