gravity_field#

This module contains a set of factory functions for setting up the gravitational potential models of celestial bodies in an environment. Below a short overview of aspects of some of the gravity field models in order to aid in properly selecting an choosing a model. Unlike most other environment model options in Tudat, there are multiple options for creating either a spherical harmonic gravity field, and a point mass gravity field:

  • Point-mass gravity field: defining the gravitational parameter manually (central()) or requiring the gravitational parameter to be extracted from Spice (central_spice()).

  • Spherical harmonic gravity field: defining all the settings manually (spherical_harmonic()), loading a pre-defined model for a solar system body (from_file_spherical_harmonic()) or calculating the spherical harmonic coefficients (up to a given degree) based on an ellipsoidal homogeneous mass distribution (spherical_harmonic_triaxial_body())

Rigid body properties will always be created automatically when a body is endowed with a gravity field, as described below:

  • Point-mass gravity field: mass computed from gravitational parameter; zero inertia tensor, and center of mass at origin of body-fixed frame

  • Spherical harmonic gravity field: mass computed from gravitational parameter, center of mass computed from degree 1 gravity field coefficients, inertia tensor as described in spherical_harmonic()

  • Polyhedron gravity field: mass computed from gravitational parameter, center of mass and inertia tensor computed from homogeneous mas distribution inside body

References#

Functions#

central(gravitational_parameter)

Factory function for central gravity field settings object.

central_spice()

Factory function to create central gravity field settings from Spice settings.

spherical_harmonic(gravitational_parameter, ...)

Factory function for creating a spherical harmonics gravity field settings object.

sh_triaxial_ellipsoid_from_density(axis_a, ...)

Factory function for spherical harmonics gravity field settings object from triaxial ellipsoid parameters, using the density to define the mass distribution.

sh_triaxial_ellipsoid_from_gravitational_parameter(...)

Factory function for spherical harmonics gravity field settings object from triaxial ellipsoid parameters, using the gravitational parameter to define the mass distribution..

from_file_spherical_harmonic(file, ...[, ...])

Factory function to load a custom spherical harmonics gravity field settings from a file.

predefined_spherical_harmonic(predefined_model)

Factory function for spherical harmonics gravity field settings of a predefined model.

polyhedron_from_mu(gravitational_parameter, ...)

Factory function for creating a polyhedron gravity field settings object, using the gravitational parameter.

polyhedron_from_density(density, ...[, ...])

Factory function for creating a polyhedron gravity field settings object, using the density.

central(gravitational_parameter: float) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function for central gravity field settings object.

Factory function for settings object, defining a point-mass gravity field model with user-defined gravitational parameter \(\mu\). The gravitational potential is the defined as:

\[U(\mathbf{r})=\frac{\mu}{||\mathbf{r}||}\]

with \(\mathbf{r}\) the position vector measured from the body’s center of mass.

Parameters:

gravitational_parameter (float) – Gravitational parameter defining the point-mass gravity field.

Returns:

Instance of the GravityFieldSettings derived CentralGravityFieldSettings class

Return type:

CentralGravityFieldSettings

Examples

In this example, we create GravityFieldSettings for Earth using a simple central gravity field model:

# define parameters describing central gravity model
gravitational_parameter = 3.986e14
# create gravity field settings
body_settings.get( "Earth" ).gravity_field_settings = environment_setup.gravity_field.central( gravitational_parameter )
central_spice() tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function to create central gravity field settings from Spice settings.

Factory function for settings object, defining a point-mass gravity field model. This function provides the same model as central()), but with gravitational parameter \(\mu\) from Spice.

Returns:

Instance of the GravityFieldSettings class of gravity field type central_spice

Return type:

GravityFieldSettings

Examples

In this example, we create GravityFieldSettings for Earth using a simple central gravity field model and data from Spice:

# create gravity field settings
body_settings.get( "Earth" ).gravity_field_settings = environment_setup.gravity_field.central_spice( )
spherical_harmonic(gravitational_parameter: float, reference_radius: float, normalized_cosine_coefficients: numpy.ndarray[numpy.float64[m, n]], normalized_sine_coefficients: numpy.ndarray[numpy.float64[m, n]], associated_reference_frame: str) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function for creating a spherical harmonics gravity field settings object.

Factory function for settings object, defining a gravity field model through spherical harmonic expansion. The associated_reference_frame must be the same frame ID as the target frame of the body’s rotation model. It represents the frame in which the spherical harmonic field is defined.

The gravitational potential is the defined as:

\[U(\mathbf{r})=\sum_{l=0}^{l_{max}}\sum_{m=0}^{l}\mu\left(\frac{{R}^{l}}{r^{l+1}}\right)\bar{P}_{lm}(\sin\phi)\left(\bar{C}_{lm}\cos m\theta+\bar{S}_{lm}\sin m\theta\right)\]

with \(\mathbf{r}\) the position vector of the evaluation point, measured from the body’s center of mass. The angles \(\phi\) and \(\theta\) are the body-fixed latitude and longitude of the evaluation point, and \(\bar{P}_{lm}\) is the associated Legendre polynomial (at degree/order \(l/m\)).

For the spherical harmonic gravity field (including other spherical harmonic functions), the normalized mean moment of inertia must be set by the user, to allow an inertia tensor to be computed. This is done using the scaled_mean_moment_of_inertia attribute of the SphericalHarmonicsGravityFieldSettings class, as in the example below

# Add gravity field model settings to body of spherical harmonic type
body_settings.get( "Mars" ).gravity_field = ...

# Add setting for moment of inertia
body_settings.get( "Mars" ).gravity_field.scaled_mean_moment_of_inertia = 0.365

This code snippet will automatically create a rigid body properties for Mars, with the inertia tensor computed from this value of 0.365 and the degree 2 gravity field coefficients. Note that, if gravity field variations are used for the body, time-variability of the degree 1- and 2- coefficients will be reflected in time-variability of the body’s center of mass and inertia tensor.

Note: Spherical harmonic coefficients used for this environment model must always be fully normalized. To normalize un-normalized spherical harmonic coefficients, see normalize_spherical_harmonic_coefficients().

Parameters:
  • gravitational_parameter (float) – Gravitational parameter \(\mu\) of gravity field.

  • reference_radius (float) – Reference radius \(R\) of spherical harmonic field expansion.

  • normalized_cosine_coefficients (numpy.ndarray) – Cosine spherical harmonic coefficients (geodesy normalized). Entry (i,j) denotes coefficient \(\bar{C}_{ij}\) at degree i and order j. As such, note that entry (0,0) of cosine coefficients should be equal to 1.

  • normalized_sine_coefficients (numpy.ndarray) – Sine spherical harmonic coefficients (geodesy normalized). Entry (i,j) denotes coefficient \(\bar{S}_{ij}\) at degree i and order j.

  • associated_reference_frame (str) – Identifier for body-fixed reference frame with which the spherical harmonics coefficients are associated.

Returns:

Instance of the GravityFieldSettings derived SphericalHarmonicsGravityFieldSettings class

Return type:

SphericalHarmonicsGravityFieldSettings

Examples

In this example, we create GravityFieldSettings for Earth using a spherical harmonics gravity model:

# Define the spherical harmonics gravity model
gravitational_parameter = 3986004.415E+8
reference_radius = 6378136.3
# Normalized coefficients taken from https://cddis.nasa.gov/archive/egm96/general_info/egm96_to360.ascii
# The above file is described in https://cddis.nasa.gov/archive/egm96/general_info/readme.egm96
normalized_cosine_coefficients = [
    [1,                   0,                   0,                   0],
    [0,                   0,                   0,                   0],
    [-0.484165371736E-03, -0.186987635955E-09, 0.243914352398E-05,  0],
    [0.957254173792E-06,  0.202998882184E-05,  0.904627768605E-06,  0.721072657057E-06]
]
normalized_sine_coefficients = [
    [0,                   0,                   0,                   0],
    [0,                   0,                   0,                   0],
    [0,                   0.119528012031E-08,  -0.140016683654E-05, 0],
    [0,                   0.248513158716E-06,  -0.619025944205E-06, 0.141435626958E-05]
]
associated_reference_frame = "IAU_Earth"
# Create the gravity field settings and add them to the body "Earth"
body_settings.get( "Earth" ).gravity_field_settings = environment_setup.gravity_field.spherical_harmonic(
    gravitational_parameter,
    reference_radius,
    normalized_cosine_coefficients,
    normalized_sine_coefficients,
    associated_reference_frame )
sh_triaxial_ellipsoid_from_density(axis_a: float, axis_b: float, axis_c: float, density: float, maximum_degree: int, maximum_order: int, associated_reference_frame: str, gravitational_constant: float = 6.67259e-11) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.SphericalHarmonicsGravityFieldSettings#

Factory function for spherical harmonics gravity field settings object from triaxial ellipsoid parameters, using the density to define the mass distribution.

Factory function for settings object, defining a gravity field model through spherical harmonic expansion of a homogeneous triaxial ellipsoid, same as spherical_harmonic The constant mass distribution in the specified ellipsoid shape is expanded to obtain a spherical harmonic coefficient representation. Gravity fields from this setting object are expressed in normalized spherical harmonic coefficients. The constant mass distribution is defined by the density and gravitational constant (optional). The body-fixed x-, y- and z- axes are assumed to be along the A-, B- and C- axes. This function implements the models of (see Balmino [1]).

Parameters:
  • axis_a (float) – Dimension of largest axis of triaxial ellipsoid.

  • axis_b (float) – Dimension of intermediate axis of triaxial ellipsoid.

  • axis_c (float) – Dimension of smallest axis of triaxial ellipsoid.

  • density (float) – Density of ellipsoid.

  • maximum_degree (int) – Maximum degree of spherical harmonics expansion.

  • maximum_order (int) – Maximum order of spherical harmonics expansion.

  • associated_reference_frame (str) – Identifier for body-fixed reference frame with which the spherical harmonics coefficients are associated.

  • gravitational_constant (float, default=physical_constants::GRAVITATIONAL_CONSTANT) – Gravitational constant G of the gravity field.

Returns:

Instance of the GravityFieldSettings derived SphericalHarmonicsGravityFieldSettings class

Return type:

SphericalHarmonicsGravityFieldSettings

Examples

In this example, we create GravityFieldSettings for Earth using the expansion of a homogeneous triaxial ellipsoid into a spherical harmonics gravity model:

# Create the gravity field settings for Earth with Spherical Harmonics from a triaxial ellipsoid
body_settings.get( "Earth" ).gravity_field_settings = environment_setup.gravity_field.spherical_harmonic_triaxial_ellipsoid_from_density(
    axis_a=6378171.88,
    axis_b=6378102.03,
    axis_c=6356752.24,
    density=5520,
    maximum_degree=50,
    maximum_order=50,
    associated_reference_frame="IAU_Earth" )
sh_triaxial_ellipsoid_from_gravitational_parameter(axis_a: float, axis_b: float, axis_c: float, maximum_degree: int, maximum_order: int, associated_reference_frame: str, gravitational_parameter: float) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.SphericalHarmonicsGravityFieldSettings#

Factory function for spherical harmonics gravity field settings object from triaxial ellipsoid parameters, using the gravitational parameter to define the mass distribution..

Factory function for settings object, defining a gravity field model through spherical harmonic expansion of a homogeneous triaxial ellipsoid, same as spherical_harmonic The constant mass distribution in the specified ellipsoid shape is expanded to obtain a spherical harmonic coefficient representation. Gravity fields from this setting object are expressed in normalized spherical harmonic coefficients. The constant mass distribution is defined by the gravitational parameter. The body-fixed x-, y- and z- axes are assumed to be along the A-, B- and C- axes. This function implements the models of (see Balmino [1]).

Parameters:
  • axis_a (float) – Dimension of largest axis of triaxial ellipsoid.

  • axis_b (float) – Dimension of intermediate axis of triaxial ellipsoid.

  • axis_c (float) – Dimension of smallest axis of triaxial ellipsoid.

  • maximum_degree (int) – Maximum degree of spherical harmonics expansion.

  • maximum_order (int) – Maximum order of spherical harmonics expansion.

  • associated_reference_frame (str) – Identifier for body-fixed reference frame with which the spherical harmonics coefficients are associated.

  • gravitational_parameter (float) – Gravitational parameter \(\mu\) of gravity field.

Returns:

Instance of the GravityFieldSettings derived SphericalHarmonicsGravityFieldSettings class

Return type:

SphericalHarmonicsGravityFieldSettings

from_file_spherical_harmonic(file: str, maximum_degree: int, maximum_order: int, associated_reference_frame: str = '', gravitational_parameter_index: int = 0, reference_radius_index: int = 1) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function to load a custom spherical harmonics gravity field settings from a file.

Factory function to load a custom spherical harmonics gravity field settings from a file. The file should contain fully normalized spherical harmonic coefficients. The associated gravitational paramerer and reference radius should be given in m^3/s^2 and m, respectively. The file format should be the same as that used for the files in the directories here. Specifically, the file should contain

  • The first line should be a series of text blocks (typically numerical data). Two of these blocks (by default the first and second one) should be the gravitational parameter and reference radius, respectively. The text block should be separated by spaces, tabs and/or commas

  • Each subsequent line should contain a set of spherical harmonic coefficients (first ordered in ascending order by degree, then in ascending order by order), where the first, second, third and fourth value of the line should be: degree \(l\), order \(m\), normalized cosine coefficient \(\bar{C}_{lm}\), normalized sine coefficient \(\bar{S}_{lm}\). Additional entries (for instance with coefficient uncertainties) are ignored.

Parameters:
  • file (str) – Full file path and name where th gravity field file is located

  • maximum_degree (int) – Maximum degree of the coefficients that are to be loaded

  • maximum_order (int) – Maximum order of the coefficients that are to be loaded

  • associated_reference_frame (str, default = "") – Name of the body-fixed reference frame to which the gravity field is to be fixed. If left empty, this reference frame will automatically be set to the body-fixed frame defined by this body’s rotation (see rotation_model for specifying rotation models).

  • gravitational_parameter_index (int, default = 0) – Index of the values in the file header (first line of file) that contains the gravitational parameter

  • reference_radius_index (int, default = 1) – Index of the values in the file header (first line of file) that contains the reference radius

Returns:

Instance of the GravityFieldSettings derived SphericalHarmonicsGravityFieldSettings class

Return type:

SphericalHarmonicsGravityFieldSettings

Examples

In this example, we create GravityFieldSettings for the Moon using GGGRX spherical harmonics gravity model, up to degree and order 300:

# Create the gravity field settings for Moon with Spherical Harmonics loaded in from a Spherical Harmonics file
body_settings.get("Moon").gravity_field_settings = environment_setup.gravity_field.from_file_spherical_harmonic(
    r"...\.tudat\resource\gravity_models\Moon\gggrx_1200l_sha.tab", 300, 300)
predefined_spherical_harmonic(predefined_model: tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.PredefinedSphericalHarmonicsModel, maximum_degree: int = -1) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function for spherical harmonics gravity field settings of a predefined model.

Factory function for spherical harmonics gravity field settings of a predefined model

Parameters:
  • predefined_model (PredefinedSphericalHarmonicsModel) – Identified for gravity field model that is to be loaded

  • maximum_degree (int, default = -1) – Maximum degree and order to which the coefficients are to be loaded. If value is negative, all coefficients for the specified gravity field are loaded

Returns:

Instance of the GravityFieldSettings derived SphericalHarmonicsGravityFieldSettings class

Return type:

SphericalHarmonicsGravityFieldSettings

Examples

In this example, we create GravityFieldSettings for Earth using EGM96 spherical harmonics gravity model, up to degree and order 32:

# Create the gravity field settings for Earth with Spherical Harmonics from a triaxial ellipsoid
body_settings.get( "Earth" ).gravity_field_settings = environment_setup.gravity_field.predefined_spherical_harmonic(
    environment_setup.gravity_field.egm96, 32 )
polyhedron_from_mu(gravitational_parameter: float, vertices_coordinates: numpy.ndarray[numpy.float64[m, n]], vertices_defining_each_facet: numpy.ndarray[numpy.int32[m, n]], associated_reference_frame: str, gravitational_constant: float = 6.67259e-11) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function for creating a polyhedron gravity field settings object, using the gravitational parameter.

Factory function for settings object, defining a gravity field model through a polyhedron. The associated_reference_frame must be the same frame ID as the target frame of the body’s rotation model. It represents the frame in which the polyhedron field is defined.

The gravitational potential, acceleration, Laplacian of potential and Hessian of potential are computed according to Werner and Scheeres [2].

This function uses the gravitational parameter to define the gravity field. To instead use the density constant see polyhedron_from_density(). Since both models tend to be computationally intensive, it is recommended to use polyhedra with the lowest number of facets that allows meeting the desired accuracy. The number of facets of a polyhedron model can be reduced using any mesh processing software, for example PyMeshLab. Additionally, different functions to process a polyhedron are available in Polyhedron utilities.

Parameters:
  • gravitational_parameter (float) – Gravitational parameter \(\mu\) of gravity field.

  • vertices_coordinates (numpy.ndarray) – Cartesian coordinates of each polyhedron vertex. Entry (i,j) denotes vertex i, coordinate j (one row per vertex, 3 columns).

  • vertices_defining_each_facet (numpy.ndarray) – Index (0 based) of the vertices constituting each facet. Entry (i,j) denotes facet i, and the jth vertex of the facet (one row per facet, 3 columns). In each row, the vertices’ indices should be ordered counterclockwise when seen from the outside of the polyhedron.

  • associated_reference_frame (str) – Identifier for body-fixed reference frame with which the spherical harmonics coefficients are associated.

  • gravitational_constant (float, default=GRAVITATIONAL_CONSTANT) – Newton’s gravitational constant G, used to computed the density

Returns:

Instance of the GravityFieldSettings derived PolyhedronGravityFieldSettings class

Return type:

PolyhedronGravityFieldSettings

polyhedron_from_density(density: float, vertices_coordinates: numpy.ndarray[numpy.float64[m, n]], vertices_defining_each_facet: numpy.ndarray[numpy.int32[m, n]], associated_reference_frame: str, gravitational_constant: float = 6.67259e-11) tudatpy.kernel.numerical_simulation.environment_setup.gravity_field.GravityFieldSettings#

Factory function for creating a polyhedron gravity field settings object, using the density.

Factory function for settings object, defining a gravity field model through a polyhedron. The associated_reference_frame must be the same frame ID as the target frame of the body’s rotation model. It represents the frame in which the polyhedron field is defined.

The gravitational potential, acceleration, Laplacian of potential and Hessian of potential are computed according to Werner and Scheeres [2].

This function uses the density to define the gravity field. To instead use the gravitational parameter see polyhedron_from_mu().

Parameters:
  • density (float, default=TUDAT_NAN) – Density of the polyhedron.

  • vertices_coordinates (numpy.ndarray) – Cartesian coordinates of each polyhedron vertex. Entry (i,j) denotes vertex i, coordinate j (one row per vertex, 3 columns).

  • vertices_defining_each_facet (numpy.ndarray) – Index (0 based) of the vertices constituting each facet. Entry (i,j) denotes facet i, and the jth vertex of the facet (one row per facet, 3 columns). In each row, the vertices’ indices should be ordered counterclockwise when seen from the outside of the polyhedron.

  • associated_reference_frame (str) – Identifier for body-fixed reference frame with which the spherical harmonics coefficients are associated.

  • gravitational_constant (float, default=GRAVITATIONAL_CONSTANT) – Newton’s gravitational constant G, used to computed the gravitational parameter

Returns:

Instance of the GravityFieldSettings derived PolyhedronGravityFieldSettings class

Return type:

PolyhedronGravityFieldSettings

Enumerations#

GravityFieldType

Enumeration of gravity field types.

PredefinedSphericalHarmonicsModel

Enumeration of predefined spherical harmonics models.

class GravityFieldType#

Enumeration of gravity field types.

Enumeration of gravity field types supported by tudat.

Members:

central_gravity :

central_spice_gravity :

spherical_harmonic_gravity :

polyhedron_gravity :

ring_gravity : No documentation found.

property name#
class PredefinedSphericalHarmonicsModel#

Enumeration of predefined spherical harmonics models.

Enumeration of predefined spherical harmonics models supported by tudat, for which thee coefficient files are automatically available (downloaded from here). The directory where these files are stored can be extracted using the get_gravity_models_path() function.

Members:

egm96 :

Coefficients for EGM96 Earth gravity field up to degree and order 200, (see link )

ggm02c :

Coefficients for the combined GGM02 Earth gravity field up to degree and order 200, (see link )

ggm02s :

Coefficients for the GRACE-only GGM02 Earth gravity field up to degree and order 160, (see link )

goco05c :

Coefficients for the GOCO05c combined Earth gravity field up to degree and order 719, (see link )

glgm3150 :

Coefficients for the GLGM3150 Moon gravity field up to degree and order 150, (see link )

lpe200 :

Coefficients for the LPE200 Moon gravity field up to degree and order 200, (see link )

gggrx1200 : No documentation found.

jgmro120d :

Coefficients for the MRO120D Moon gravity field up to degree and order 120, (see link )

jgmess160a :

Coefficients for the MESS160A Moon gravity field up to degree and order 160, (see link )

shgj180u :

Coefficients for the SHGJ180U Moon gravity field up to degree and order 180, (see link )

property name#

Classes#

GravityFieldSettings

Base class for providing settings for automatic gravity field model creation.

CentralGravityFieldSettings

GravityFieldSettings derived class defining settings of point mass gravity field.

SphericalHarmonicsGravityFieldSettings

GravityFieldSettings derived class defining settings of spherical harmonic gravity field representation.

PolyhedronGravityFieldSettings

GravityFieldSettings derived class defining settings of a polyhedron gravity field representation.

class GravityFieldSettings#

Base class for providing settings for automatic gravity field model creation.

This class is a functional base class for settings of gravity field models that require no information in addition to their type. Gravity field model classes requiring additional information must be created using an object derived from this class.

property gravity_field_type#

read-only

Type of gravity field model that is to be created.

Type:

GravityFieldType

class CentralGravityFieldSettings#

GravityFieldSettings derived class defining settings of point mass gravity field.

Derived class of GravityFieldSettings for central gravity fields, which are defined by a single gravitational parameter.

property gravitational_parameter#

Gravitational parameter of central gravity field.

Type:

float

class SphericalHarmonicsGravityFieldSettings#

GravityFieldSettings derived class defining settings of spherical harmonic gravity field representation.

Derived class of GravityFieldSettings for gravity fields, which are defined by a spherical harmonic gravity field representation.

property associated_reference_frame#

Identifier for body-fixed reference frame with which the coefficients are associated.

Type:

str

property create_time_dependent_field#

Boolean that denotes whether the field should be created as time-dependent (even if no variations are imposed initially).

Type:

bool

property gravitational_parameter#

Gravitational parameter of gravity field.

Type:

float

property normalized_cosine_coefficients#

Cosine spherical harmonic coefficients (geodesy normalized). Entry (i,j) denotes coefficient at degree i and order j.

Type:

numpy.ndarray

property normalized_sine_coefficients#

Sine spherical harmonic coefficients (geodesy normalized). Entry (i,j) denotes coefficient at degree i and order j.

Type:

numpy.ndarray

property reference_radius#

read-only

Reference radius of spherical harmonic field expansion.

Type:

float

property scaled_mean_moment_of_inertia#

Value of the scaled mean moment of inertia \(I_{xx}+I_{yy}+I_{zz}/(MR^{2})\). This value does not influence the gravity field itself, but together with the degree 2 gravity field coefficients defines the body’s inertia tensor.

Type:

float

class PolyhedronGravityFieldSettings#

GravityFieldSettings derived class defining settings of a polyhedron gravity field representation.

Derived class of GravityFieldSettings for gravity fields, which are defined by a polyhedron gravity field representation.

property associated_reference_frame#

Identifier for body-fixed reference frame with which the vertices coordinates are associated.

Type:

str

property density#

Density of the polyhedron.

Type:

float

property gravitational_parameter#

Gravitational parameter of gravity field.

Type:

float

property vertices_coordinates#

Cartesian coordinates of each polyhedron vertex. Entry (i,j) denotes vertex i, coordinate j (one row per vertex, 3 columns).

Type:

numpy.ndarray

property vertices_defining_each_facet#

Index (0 based) of the vertices constituting each facet. Entry (i,j) denotes facet i, and the jth vertex of the facet (one row per facet, 3 columns). In each row, the vertices’ indices should be ordered counterclockwise when seen from the outside of the polyhedron.

Type:

numpy.ndarray