Skip to content

Gravity Models

Gravity acceleration functions including point-mass and spherical harmonic models.

Note

For conceptual explanations and examples, see Gravity Models in the Learn section.

Point-Mass Gravity

accel_point_mass_gravity builtin

accel_point_mass_gravity(r_object: ndarray, r_central_body: ndarray, gm: float) -> ndarray

Compute acceleration due to point-mass gravity.

Accepts either a 3D position vector or a 6D state vector for r_object.

Parameters:

Name Type Description Default
r_object ndarray

Position (length 3) or state (length 6) of the object. Units: (m)

required
r_central_body ndarray

Position vector of the central body. Units: (m)

required
gm float

Gravitational parameter. Units: (m³/s²)

required

Returns:

Type Description
ndarray

np.ndarray: Acceleration due to gravity. Units: (m/s²)

Example
1
2
3
4
5
6
import brahe as bh
import numpy as np

r_object = np.array([bh.R_EARTH, 0.0, 0.0])
r_central_body = np.array([0.0, 0.0, 0.0])
a_grav = bh.accel_point_mass_gravity(r_object, r_central_body, bh.GM_EARTH)

Spherical Harmonic Gravity

accel_gravity_spherical_harmonics builtin

accel_gravity_spherical_harmonics(r_eci: ndarray, R_i2b: ndarray, gravity_model: GravityModel, n_max: int, m_max: int) -> ndarray

Compute acceleration due to spherical harmonic gravity model.

This function computes the gravitational acceleration on an object using a spherical harmonic expansion of Earth's gravity field. It transforms the position to body-fixed coordinates, evaluates the spherical harmonics, and transforms the acceleration back to the inertial frame.

Accepts either a 3D position vector or a 6D state vector for r_eci.

Parameters:

Name Type Description Default
r_eci ndarray

Position (length 3) or state (length 6) in ECI frame. Units: (m)

required
R_i2b ndarray

Rotation matrix from ECI to body-fixed frame (3x3)

required
gravity_model GravityModel

Gravity model to use

required
n_max int

Maximum degree to evaluate (n_max <= model.n_max)

required
m_max int

Maximum order to evaluate (m_max <= min(n_max, model.m_max))

required

Returns:

Type Description
ndarray

np.ndarray: Acceleration in ECI frame. Units: (m/s²)

Raises:

Type Description
Exception

If n_max or m_max exceed model limits or if m_max > n_max

Example
import brahe as bh
import numpy as np

# Initialize EOP for frame transformations
bh.initialize_eop()

# Load gravity model
model = bh.GravityModel.from_model_type(bh.GravityModelType.JGM3)

# Create test position
r_eci = np.array([6525.919e3, 1710.416e3, 2508.886e3])

# Get rotation matrix (or use identity for simplified case)
R = np.eye(3)

# Compute acceleration
a_grav = bh.accel_gravity_spherical_harmonics(r_eci, R, model, 20, 20)
print(f"Acceleration: {a_grav} m/s²")

Gravity Model Class

GravityModel

GravityModel()

Spherical harmonic gravity model for high-fidelity gravitational acceleration computation.

This class represents a spherical harmonic expansion of Earth's gravitational potential, allowing for accurate modeling of Earth's non-uniform gravity field.

Initialize instance.

gm property

gm: Any

Gravitational parameter. Units: (m³/s²)

m_max property

m_max: Any

Maximum order of spherical harmonic expansion.

model_errors property

model_errors: Any

Error estimation type.

model_name property

model_name: Any

Name of the gravity model.

n_max property

n_max: Any

Maximum degree of spherical harmonic expansion.

normalization property

normalization: Any

Coefficient normalization convention.

radius property

radius: Any

Reference radius. Units: (m)

tide_system property

tide_system: Any

Tide system convention.

compute_spherical_harmonics method descriptor

compute_spherical_harmonics(r_body: ndarray, n_max: int, m_max: int) -> ndarray

Compute gravitational acceleration in body-fixed frame using spherical harmonics.

Parameters:

Name Type Description Default
r_body ndarray

Position vector in body-fixed frame. Units: (m)

required
n_max int

Maximum degree to evaluate (n_max <= model.n_max)

required
m_max int

Maximum order to evaluate (m_max <= min(n_max, model.m_max))

required

Returns:

Type Description
ndarray

np.ndarray: Acceleration in body-fixed frame. Units: (m/s²)

Raises:

Type Description
Exception

If n_max or m_max exceed model limits or if m_max > n_max

Example
1
2
3
4
5
6
7
8
import brahe as bh
import numpy as np

model = bh.GravityModel.from_model_type(bh.GravityModelType.JGM3)
r_body = np.array([6525.919e3, 1710.416e3, 2508.886e3])

# Compute using 20x20 expansion
a_body = model.compute_spherical_harmonics(r_body, 20, 20)

from_file builtin

from_file(filepath: str) -> GravityModel

Load gravity model from a .gfc file.

Parameters:

Name Type Description Default
filepath str

Path to the .gfc gravity model file

required

Returns:

Name Type Description
GravityModel GravityModel

Loaded gravity model

Raises:

Type Description
Exception

If file cannot be loaded or parsed

Example
1
2
3
4
5
import brahe as bh

# Load custom gravity model from file
model = bh.GravityModel.from_file("path/to/model.gfc")
print(f"Model: {model.model_name}, {model.n_max}x{model.m_max}")

from_model_type builtin

from_model_type(model_type: GravityModelType) -> GravityModel

Load a gravity model from a GravityModelType.

Parameters:

Name Type Description Default
model_type GravityModelType

Which model to load (packaged or from file)

required

Returns:

Name Type Description
GravityModel GravityModel

Loaded gravity model

Raises:

Type Description
Exception

If file loading fails (for FromFile variant)

Example
import brahe as bh

# Load JGM3 70x70 model
model = bh.GravityModel.from_model_type(bh.GravityModelType.JGM3)
print(f"Loaded {model.model_name}")

# Load EGM2008 360x360 model
model_hifi = bh.GravityModel.from_model_type(bh.GravityModelType.EGM2008_360)

# Load from custom file
custom_type = bh.GravityModelType.from_file("/path/to/model.gfc")
model_custom = bh.GravityModel.from_model_type(custom_type)

get method descriptor

get(n: int, m: int) -> Tuple

Get spherical harmonic coefficients for a specific degree and order.

Parameters:

Name Type Description Default
n int

Degree (0 <= n <= n_max)

required
m int

Order (0 <= m <= min(n, m_max))

required

Returns:

Name Type Description
tuple Tuple

(C_nm, S_nm) coefficients

Raises:

Type Description
Exception

If n or m are out of bounds

Example
1
2
3
4
5
6
7
import brahe as bh

model = bh.GravityModel.from_model_type(bh.GravityModelType.JGM3)

# Get J2 coefficient (C20)
c20, s20 = model.get(2, 0)
print(f"J2 = {-c20}")

set_max_degree_order method descriptor

set_max_degree_order(n: int, m: int) -> Any

Truncate the gravity model to a smaller degree and order to save memory.

This method resizes the internal coefficient matrix in-place, discarding coefficients beyond the specified limits. This is useful when loading a high-fidelity model but only needing a lower-degree expansion.

Parameters:

Name Type Description Default
n int

New maximum degree (must be <= current n_max)

required
m int

New maximum order (must be <= min(n, current m_max))

required

Raises:

Type Description
ValueError

If m > n or if n/m exceed current model limits

Example
1
2
3
4
5
6
7
8
9
import brahe as bh

# Load full EGM2008 360x360 model
model = bh.GravityModel.from_model_type(bh.GravityModelType.EGM2008_360)
print(f"Before: {model.n_max}x{model.m_max}")  # 360x360

# Truncate to 70x70 to save memory
model.set_max_degree_order(70, 70)
print(f"After: {model.n_max}x{model.m_max}")  # 70x70

See Also