Small Body Geophysical Analysis Tool (SBGAT)
|
Computes volume, area, shape index, center of mass, inertia tensor and principal axes of a polyhedral mesh of constant density. More...
#include <SBGATMassProperties.hpp>
Public Member Functions | |
vtkTypeMacro (SBGATMassProperties, vtkPolyDataAlgorithm) | |
void | PrintSelf (std::ostream &os, vtkIndent indent) override |
void | PrintHeader (std::ostream &os, vtkIndent indent) override |
void | PrintTrailer (std::ostream &os, vtkIndent indent) override |
double | GetVolume () const |
double | GetVolumeProjected () const |
double | GetVolumeX () |
double | GetVolumeY () |
double | GetVolumeZ () |
double | GetKx () const |
double | GetKy () const |
double | GetKz () const |
double | GetSurfaceArea () const |
double | GetMinCellArea () const |
double | GetMaxCellArea () const |
bool | CheckClosed () const |
double | GetNormalizedShapeIndex () const |
const arma::vec::fixed< 3 > & | GetCenterOfMass () const |
void | GetCenterOfMass (double *com) const |
arma::mat::fixed< 3, 3 > | GetNormalizedInertiaTensor () const |
arma::mat::fixed< 3, 3 > | GetUnitDensityInertiaTensor () const |
arma::mat::fixed< 3, 3 > | GetPrincipalAxes () const |
arma::vec::fixed< 3 > | GetPrincipalDimensions () const |
arma::vec::fixed< 3 > | GetNormalizedInertiaMoments () const |
arma::vec::fixed< 3 > | GetUnitDensityInertiaMoments () const |
double | GetAverageRadius () const |
void | SaveMassProperties (std::string path) const |
double | GetDeltaV (const int &f) const |
arma::vec::fixed< 3 > | GetDeltaCM (const int &f) const |
arma::mat::fixed< 3, 3 > | GetDeltaIOverDeltaV (const int &f) const |
arma::vec::fixed< 6 > | GetDeltaIf (const int &f) const |
![]() | |
vtkTypeMacro (SBGATFilter, vtkPolyDataAlgorithm) | |
void | PrintSelf (std::ostream &os, vtkIndent indent) override |
void | PrintHeader (std::ostream &os, vtkIndent indent) override |
void | PrintTrailer (std::ostream &os, vtkIndent indent) override |
void | SetScaleMeters () |
void | SetScaleKiloMeters () |
double | GetScaleFactor () const |
void | SetDensity (const double density) |
double | GetDensity () const |
void | GetVerticesInFacet (const int &f, double *r0, double *r1, double *r2) const |
arma::vec::fixed< 3 > | GetFacetCenter (const int &f) const |
void | GetVerticesOnEdge (const int &e, double *r0, double *r1) const |
void | GetFacetNormal (const int &f, double *n) const |
void | GetIndicesOfAdjacentFacets (const int &e, int &f0, int &f1) const |
void | GetIndicesVerticesOnEdge (const int &e, int &v0, int &v1) const |
void | GetIndicesVerticesInFacet (const int &f, int &v0, int &v1, int &v2) const |
arma::vec::fixed< 3 > | GetNonNormalizedFacetNormal (const int &f) const |
arma::vec::fixed< 3 > | GetRe (const arma::vec::fixed< 3 > &pos, const int &e) const |
arma::vec::fixed< 3 > | GetRe (const double *pos, const int &e) const |
arma::vec::fixed< 3 > | GetRf (const arma::vec::fixed< 3 > &pos, const int &f) const |
arma::vec::fixed< 3 > | GetRf (const double *pos, const int &f) const |
void | GetBoundingBox (double &xmin, double &xmax, double &ymin, double &ymax, double &zmin, double &zmax) const |
double | GetEdgeLength (const int &e) const |
int | GetN_facets () const |
int | GetN_edges () const |
int | GetN_vertices () const |
bool | GetIsInMeters () const |
bool | GetIsInKiloMeters () const |
Static Public Member Functions | |
static SBGATMassProperties * | New () |
static void | ComputeAndSaveMassProperties (vtkSmartPointer< vtkPolyData > shape, std::string path) |
![]() | |
static SBGATFilter * | New () |
Protected Member Functions | |
SBGATMassProperties () | |
~SBGATMassProperties () override | |
int | RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override |
![]() | |
SBGATFilter () | |
~SBGATFilter () override | |
void | Clear () |
int | RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override |
Protected Attributes | |
arma::vec::fixed< 3 > | center_of_mass |
arma::mat::fixed< 3, 3 > | inertia_tensor |
arma::mat::fixed< 3, 3 > | principal_axes |
arma::vec::fixed< 3 > | normalized_principal_moments |
arma::vec::fixed< 3 > | unit_density_principal_moments |
arma::mat::fixed< 3, 3 > | unit_density_inertia_tensor |
arma::vec::fixed< 3 > | principal_dimensions |
double | SurfaceArea |
double | MinCellArea |
double | MaxCellArea |
double | Volume |
double | VolumeProjected |
double | VolumeX |
double | VolumeY |
double | VolumeZ |
double | Kx |
double | Ky |
double | Kz |
double | NormalizedShapeIndex |
double | r_avg |
bool | IsClosed |
![]() | |
double ** | facet_normals |
double ** | vertices |
double | scaleFactor = 1 |
double | density |
int ** | edges |
int ** | facets |
int ** | edge_facets_ids |
bool | is_in_meters = true |
int | N_facets |
int | N_edges |
int | N_vertices |
double | bounds [6] |
Private Member Functions | |
SBGATMassProperties (const SBGATMassProperties &)=delete | |
void | operator= (const SBGATMassProperties &)=delete |
Computes volume, area, shape index, center of mass, inertia tensor and principal axes of a polyhedral mesh of constant density.
Computes the volume, the surface area, and the normalized shape index, center of mass and inertia tensor of a topologically-closed, constant-density polyhedron. This class will always use results expressed in meters
as their distance unit (e.g center-of-mass coordinates in meters, volume in m^3,...) . Unit consistency is enforced through the use of the SetScaleMeters() and SetScaleKiloMeters() method.
See "Inertia of Any Polyhedron" by Anthony R. Dobrovolskis, Icarus 124, 698–704 (1996) Article No. 0243 for further details. Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
Definition at line 42 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 56 of file SBGATMassProperties.cpp.
|
overrideprotected |
Definition at line 75 of file SBGATMassProperties.cpp.
|
privatedelete |
|
inline |
Checks whether the polydata is topologically closed or open If closed, the sum of the oriented surface area should be equal to zero
Definition at line 107 of file SBGATMassProperties.hpp.
|
static |
Computes the mass properties of the provided shape and saves the results to a JSON file
shape | pointer to considered shape |
path | savepath (ex: "mass_properties.json") |
Definition at line 401 of file SBGATMassProperties.cpp.
|
inline |
Return the average radius of the shape (that is, the radius of a sphere occupying the same volume) (m)
Definition at line 199 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the coordinates of the center of mass (m) evaluated in the frame of origin assuming a constant density distribution across the shape
Definition at line 124 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the coordinates of the center of mass (m) evaluated in the frame of origin assuming a constant density distribution across the shape
Definition at line 133 of file SBGATMassProperties.hpp.
arma::vec::fixed< 3 > SBGATMassProperties::GetDeltaCM | ( | const int & | f | ) | const |
Return coordinates of the tetrahedron's center-of-mass (m)
f | facet index |
Definition at line 543 of file SBGATMassProperties.cpp.
arma::vec::fixed< 6 > SBGATMassProperties::GetDeltaIf | ( | const int & | f | ) | const |
Return the parametrization of the the unit-density tetrahedron's inertia tensor
f | facet index |
Definition at line 583 of file SBGATMassProperties.cpp.
arma::mat::fixed< 3, 3 > SBGATMassProperties::GetDeltaIOverDeltaV | ( | const int & | f | ) | const |
Return the unit-density tetrahedron's inertia tensor divided by tetrahedron's signed volume (m^2)
f | facet index |
Definition at line 559 of file SBGATMassProperties.cpp.
double SBGATMassProperties::GetDeltaV | ( | const int & | f | ) | const |
Return signed contribution to total volume of tetrahedron subtended by facet f (m^3)
f | facet index |
Definition at line 529 of file SBGATMassProperties.cpp.
|
inline |
Compute and return the weighting factors for the maximum unit normal component (MUNC).
Definition at line 80 of file SBGATMassProperties.hpp.
|
inline |
Definition at line 81 of file SBGATMassProperties.hpp.
|
inline |
Definition at line 82 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the max cell area in m^2
Definition at line 99 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the min cell area in m^2
Definition at line 93 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the normalized inertia moments assuming uniform density distribution across the shape, sorted from the smallest inertia to the largest. The normalization applied to the inertia tensor is I_norm = I / (mass * r_avg ^ 2) where r_avg = cbrt(3/4*Volume/pi)
Definition at line 184 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the dimensionless inertia tensor evaluated in the frame of origin assuming a constant density distribution across the shape. The normalization applied to the inertia tensor is I_norm = I / (mass * r_avg ^ 2) where r_avg = cbrt(3/4*Volume/pi)
Definition at line 145 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the normalized shape index. This characterizes the deviation of the shape of an object from a sphere. A sphere's NSI is one. This number is always >= 1.0.
Definition at line 114 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the dcm orienting the principal axes of the small body relative to the body coordinates frame. That is, denoting P the principal frame and B the frame in which the coordinates of the body are currently expressed, this method returns [PB]
Definition at line 165 of file SBGATMassProperties.hpp.
|
inline |
Computes and returns the principal dimensions (m) of the ellipsoid associated with the inertia tensor tensor, sorted from the longest (smallest inertia) to shortest (largest inertia)
Definition at line 174 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the area in m^2
Definition at line 87 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the inertia moments assuming uniform unit density distribution across the shape, sorted from the smallest inertia to the largest.
Definition at line 192 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the dimensionless inertia tensor evaluated in the frame of origin assuming a constant density distribution across the shape. The normalization applied to the inertia tensor is I_norm = I / (rho) where rho is the density
Definition at line 154 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the volume (m^3)
Definition at line 57 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the projected volume. Typically you should compare this volume to the value returned by GetVolume if you get an error (GetVolume()-GetVolumeProjected())*10000 that is greater than GetVolume() this should identify a problem:
Definition at line 67 of file SBGATMassProperties.hpp.
|
inline |
Compute and return the volume projected on to each axis aligned plane.
Definition at line 72 of file SBGATMassProperties.hpp.
|
inline |
Definition at line 73 of file SBGATMassProperties.hpp.
|
inline |
Definition at line 74 of file SBGATMassProperties.hpp.
|
static |
Constructs with initial values of zero.
|
privatedelete |
|
override |
Definition at line 367 of file SBGATMassProperties.cpp.
|
override |
Definition at line 376 of file SBGATMassProperties.cpp.
|
override |
Definition at line 370 of file SBGATMassProperties.cpp.
|
overrideprotected |
Definition at line 85 of file SBGATMassProperties.cpp.
void SBGATMassProperties::SaveMassProperties | ( | std::string | path | ) | const |
Save the computed mass properties to a JSON file
path | savepath (ex: "mass_properties.json") |
Definition at line 413 of file SBGATMassProperties.cpp.
SBGATMassProperties::vtkTypeMacro | ( | SBGATMassProperties | , |
vtkPolyDataAlgorithm | |||
) |
|
protected |
Definition at line 261 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 262 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 284 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 279 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 280 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 281 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 273 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 272 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 265 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 282 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 263 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 268 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 283 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 271 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 267 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 266 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 274 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 275 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 276 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 277 of file SBGATMassProperties.hpp.
|
protected |
Definition at line 278 of file SBGATMassProperties.hpp.