Python API - polyhedral_gravity#
Module Overview#
The evaluation of the polyhedral gravity model requires the following parameters:
Name |
---|
Polyhedral Mesh (either as vertices & faces or as polyhedral source files) |
Constant Density \(\rho\) |
In the Python Interface, you define these parameters as polyhedral_gravity.Polyhedron
from polyhedral_gravity import Polyhedron, GravityEvaluable, evaluate, PolyhedronIntegrity, NormalOrientation
polyhedron = Polyhedron(
polyhedral_source=(vertices, faces), # (N,3) and (M,3) array-like
density=density, # Density of the Polyhedron, Unit must match to the mesh's scale
normal_orientation=NormalOrientation.OUTWARDS, # Possible values OUTWARDS (default) or INWARDS
integrity_check=PolyhedronIntegrity.VERIFY, # Possible values AUTOMATIC (default), VERIFY, HEAL, DISABLE
)
Note
Tsoulis et al.’s formulation requires that the normals point OUTWARDS
.
The implementation can handle both cases and also can automatically determine the property if initiall set wrong.
Using AUTOMATIC
(default for first-time-user) or VERIFY
raises a ValueError
if the polyhedral_gravity.NormalOrientation
is wrong.
Using HEAL
will re-order the vertex sorting to fix errors.
Using DISABLE
will turn this check off and avoid \(O(n^2)\) runtime complexcity of this check! Highly recommened, when you “know your mesh”!
The polyhedron’s mesh’s units must match with the constant density! For example, if the mesh is in \([m]\), then the constant density should be in \([\frac{kg}{m^3}]\).
Afterwards one can use polyhedral_gravity.evaluate()
the gravity at a single point P via:
potential, acceleration, tensor = evaluate(
polyhedron=polyhedron,
computation_points=P,
parallel=True,
)
or via use the cached approach polyhedral_gravity.GravityEvaluable
(desriable for subsequent evaluations using the same polyhedral_gravity.Polyhedron
)
evaluable = GravityEvaluable(polyhedron=polyhedron)
potential, acceleration, tensor = evaluable(
computation_points=P,
parallel=True,
)
Note
If P
would be an array of points, the return value would be a List[Tuple[potential, acceleration, tensor]]
!
The calculation outputs the following parameters for every Computation Point P. The units of the respective output depend on the units of the input parameters (mesh and density)!
Hence, if e.g. your mesh is in \(km\), the density must match. Further, the output units will match the input units.
Name |
If mesh \([m]\) and density \([\frac{kg}{m^3}]\) |
Comment |
---|---|---|
\(V\) |
\(\frac{m^2}{s^2}\) or \(\frac{J}{kg}\) |
The potential or also called specific energy |
\(V_x\), \(V_y\), \(V_z\) |
\(\frac{m}{s^2}\) |
The gravitational acceleration in the three cartesian directions |
\(V_{xx}\), \(V_{yy}\), \(V_{zz}\), \(V_{xy}\), \(V_{xz}\), \(V_{yz}\) |
\(\frac{1}{s^2}\) |
The spatial rate of change of the gravitational acceleration |
This model’s output obeys to the geodesy and geophysics sign conventions. Hence, the potential \(V\) for a polyhedron with a mass \(m > 0\) is defined as positive.
The accelerations \(V_x\), \(V_y\), \(V_z\) are defined as
Accordingly, the second derivative tensor is defined as the derivative of \(\textbf{g}\).
Polyhedron#
- class polyhedral_gravity.Polyhedron#
A constant density Polyhedron stores the mesh data consisting of vertices and triangular faces. The density and the coordinate system in which vertices and faces are defined need to have the same scale/ units. Tsoulis et al.’s polyhedral gravity model requires that the plane unit normals of every face are pointing outwards of the polyhedron. Otherwise the results are negated. The class by default enforces this constraints and offers utility to (automatically) make the input data obey to this constraint.
- __getitem__(self: polyhedral_gravity.Polyhedron, index: int) Annotated[list[Annotated[list[float], FixedSize(3)]], FixedSize(3)] #
Returns the the three cooridnates of the vertices making the face at the requested index. This does not return the face as list of vertex indices, but resolved with the actual coordinates.
- Parameters:
index – The index of the face
- Returns:
The resolved face
- Return type:
\((3, 3)\)-array-like
- Raises:
IndexError if face index is out-of-bounds –
- __init__(self: polyhedral_gravity.Polyhedron, polyhedral_source: Union[tuple[list[Annotated[list[float], FixedSize(3)]], list[Annotated[list[int], FixedSize(3)]]], list[str]], density: float, normal_orientation: polyhedral_gravity.NormalOrientation = <NormalOrientation.OUTWARDS: 0>, integrity_check: polyhedral_gravity.PolyhedronIntegrity = <PolyhedronIntegrity.AUTOMATIC: 2>) None #
Creates a new Polyhedron from vertices and faces and a constant density. If the integrity_check is not set to DISABLE, the mesh integrity is checked (so that it fits the specification of the polyhedral model by Tsoulis et al.)
- Parameters:
polyhedral_source – The vertices (\((N, 3)\)-array-like) and faces (\((M, 3)\)-array-like) of the polyhedron as pair or The filenames of the files containing the vertices & faces as list of strings
density – The constant density of the polyhedron, it must match the mesh’s units, e.g. mesh in \([m]\) then density in \([kg/m^3]\)
normal_orientation – The pointing direction of the mesh’s plane unit normals, i.e., either
OUTWARDS
orINWARDS
of the polyhedron. One ofpolyhedral_gravity.NormalOrientation
. (default:OUTWARDS
)integrity_check –
Conducts an Integrity Check (degenerated faces/ vertex ordering) depending on the values. One of
polyhedral_gravity.PolyhedronIntegrity
:AUTOMATIC
(Default): Prints to stdout and throws ValueError if normal_orientation is wrong/ inconsistenVERIFY
: LikeAUTOMATIC
, but does not print to stdoutDISABLE
: Recommened, when you are familiar with the mesh to avoid \(O(n^2)\) runtime cost. Disables ALL checksHEAL
: Automatically fixes the normal_orientation and vertex ordering to the correct values
- Raises:
ValueError – If the faces array does not contain a reference to vertex 0 indicating an index start at 1
ValueError – If
integrity_check
is set toAUTOMATIC
orVERIFY
and the mesh is inconsistent
Note
The
integrity_check
is automatically enabled to avoid wrong results due to the wrong vertex ordering. The check requires \(O(n^2)\) operations. You want to turn this off, when you know you mesh!
- __repr__(self: polyhedral_gravity.Polyhedron) str #
str
: A string representation of this polyhedron
- check_normal_orientation(self: polyhedral_gravity.Polyhedron) tuple[polyhedral_gravity.NormalOrientation, set[int]] #
Returns a tuple consisting of majority plane unit normal orientation, i.e. the direction in which at least more than half of the plane unit normals point, and the indices of the faces violating this orientation, i.e. the faces whose plane unit normals point in the other direction. The set of incides vioalting the property is empty if the mesh has a clear ordering. The set contains values if the mesh is incosistent.
- Returns:
Tuple consisting consisting of majority plane unit normal orientation and the indices of the faces violating this orientation.
Note
This utility is mainly for diagnostics and debugging purposes. If the polyhedron is constrcuted with integrity_check set to
AUTOMATIC
orVERIFY
, the construction fails anyways. If set toHEAL
, this method should return an empty set (but maybe a different ordering than initially specified) Only if set toDISABLE
, then this method might actually return a set with faulty indices. Hence, if you want to know your mesh error. Construct the polyhedron withintegrity_check=DISABLE
and call this method.
- property density#
The density of the polyhedron (Read/ Write)
- Type:
float
- property faces#
The faces of the polyhedron (Read-Only)
- Type:
(M, 3)-array-like of
int
- property normal_orientation#
The orientation of the plane unit normals (Read-Only)
- property vertices#
The vertices of the polyhedron (Read-Only)
- Type:
(N, 3)-array-like of
float
Enums to specify MeshChecks#
- class polyhedral_gravity.NormalOrientation#
The orientation of the plane unit normals of the polyhedron. Tsoulis et al. equations require the normals to point outwards of the polyhedron. If the opposite hold, the result is negated. The implementation can handle both cases.
Members:
OUTWARDS : Outwards pointing plane unit normals
INWARDS : Inwards pointing plane unit normals
- class polyhedral_gravity.PolyhedronIntegrity#
The pointing direction of the normals of a Polyhedron. They can either point outwards or inwards the polyhedron.
Members:
DISABLE : All activities regarding MeshChecking are disabled. No runtime overhead!
VERIFY : Only verification of the NormalOrientation. A misalignment (e.g. specified OUTWARDS, but is not) leads to a runtime_error. Runtime Cost \(O(n^2)\)
AUTOMATIC : Like
VERIFY
, but also informs the user about the option in any case on the runtime costs. This is the implicit default option. Runtime Cost: \(O(n^2)\) and output to stdout in every case!HEAL : Verification and Autmatioc Healing of the NormalOrientation. A misalignemt does not lead to a runtime_error, but to an internal correction of vertices ordering. Runtime Cost: \(O(n^2)\)
GravityModel#
Single Function#
- polyhedral_gravity.evaluate(polyhedron: polyhedral_gravity.Polyhedron, computation_points: Annotated[list[float], FixedSize(3)] | list[Annotated[list[float], FixedSize(3)]], parallel: bool = True) tuple[float, Annotated[list[float], FixedSize(3)], Annotated[list[float], FixedSize(6)]] | list[tuple[float, Annotated[list[float], FixedSize(3)], Annotated[list[float], FixedSize(6)]]] #
Evaluates the polyhedral gravity model for a given constant density polyhedron at a given computation point.
- Parameters:
polyhedron – The polyhedron for which to evaluate the gravity model
computation_points – The computation points as tuple or list of points
parallel – If
True
, the computation is done in parallel (default:True
)
- Returns:
Either a triplet of potential \(V\), acceleration \([V_x, V_y, V_z]\) and second derivatives \([V_{xx}, V_{yy}, V_{zz}, V_{xy},V_{xz}, V_{yz}]\) at the computation points or if multiple computation points are given a list of these triplets
Cached Evaluation#
- class polyhedral_gravity.GravityEvaluable#
A class to evaluate the polyhedral gravity model for a given constant density polyhedron at a given computation point. It provides a
poylhedral_gravity.GravityEvaluable.__call__()
method to evaluate the polyhedral gravity model for computation points while also caching the polyhedron & intermediate results over the lifetime of the object.- __call__(self: polyhedral_gravity.GravityEvaluable, computation_points: Annotated[list[float], FixedSize(3)] | list[Annotated[list[float], FixedSize(3)]], parallel: bool = True) tuple[float, Annotated[list[float], FixedSize(3)], Annotated[list[float], FixedSize(6)]] | list[tuple[float, Annotated[list[float], FixedSize(3)], Annotated[list[float], FixedSize(6)]]] #
Evaluates the polyhedral gravity model for a given constant density polyhedron at a given computation point.
- Parameters:
computation_points – The computation points as tuple or list of points
parallel – If
True
, the computation is done in parallel (default:True
)
- Returns:
Either a triplet of potential \(V\), acceleration \([V_x, V_y, V_z]\) and second derivatives \([V_{xx}, V_{yy}, V_{zz}, V_{xy},V_{xz}, V_{yz}]\) at the computation points or if multiple computation points are given a list of these triplets
- __init__(self: polyhedral_gravity.GravityEvaluable, polyhedron: polyhedral_gravity.Polyhedron) None #
Creates a new GravityEvaluable for a given constant density polyhedron. It provides a
poylhedral_gravity.GravityEvaluable.__call__()
method to evaluate the polyhedral gravity model for computation points while also caching the polyhedron & intermediate results over the lifetime of the object.- Parameters:
polyhedron – The polyhedron for which to evaluate the gravity model
- __repr__(self: polyhedral_gravity.GravityEvaluable) str #
str
: A string representation of this GravityEvaluable