Small Body Geophysical Analysis Tool (SBGAT)
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
SBGATPolyhedronGravityModel Class Reference

Evaluation of potential, acceleration caused by a constant-density polyhedron. More...

#include <SBGATPolyhedronGravityModel.hpp>

Inheritance diagram for SBGATPolyhedronGravityModel:
Inheritance graph
[legend]
Collaboration diagram for SBGATPolyhedronGravityModel:
Collaboration graph
[legend]

Public Member Functions

 vtkTypeMacro (SBGATPolyhedronGravityModel, 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 GetPotential (double const *point) const
 
double GetPotential (const arma::vec::fixed< 3 > &point) const
 
void GetPotentialAcceleration (double const *point, double &potential, arma::vec::fixed< 3 > &acc) const
 
void GetPotentialAcceleration (const arma::vec::fixed< 3 > &point, double &potential, arma::vec::fixed< 3 > &acc) const
 
void GetPotentialAccelerationGravityGradient (double const *point, double &potential, arma::vec::fixed< 3 > &acc, arma::mat::fixed< 3, 3 > &gravity_gradient_mat) const
 
arma::vec::fixed< 3 > GetBodyFixedAccelerationf (const int &f) const
 
void GetPotentialAccelerationGravityGradient (const arma::vec::fixed< 3 > &point, double &potential, arma::vec::fixed< 3 > &acc, arma::mat::fixed< 3, 3 > &gravity_gradient_mat) const
 
arma::mat::fixed< 3, 3 > GetGravityGradient (const arma::vec::fixed< 3 > &point) const
 
arma::vec::fixed< 3 > GetAcceleration (const arma::vec::fixed< 3 > &point) const
 
arma::vec::fixed< 3 > GetAcceleration (double const *point) const
 
bool Contains (double const *point, double tol=1e-8) const
 
bool Contains (const arma::vec::fixed< 3 > &point, double tol=1e-8) const
 
double GetSlope (const int &f) const
 
double GetMass () const
 
double GetPerformanceFactor (const arma::vec::fixed< 3 > &pos, const int &f) const
 
double GetPerformanceFactor (const double *pos, const int &f) const
 
double GetLe (const arma::vec::fixed< 3 > &pos, const int &e) const
 
double GetLe (const double *pos, const int &e) const
 
arma::vec::fixed< 10 > GetXe (const arma::vec::fixed< 3 > &pos, const int &e) const
 
arma::vec::fixed< 10 > GetXf (const arma::vec::fixed< 3 > &pos, const int &f) const
 
arma::vec::fixed< 6 > GetEeParam (const int &e) const
 
arma::vec::fixed< 6 > GetFfParam (const int &f) const
 
void SetOmega (arma::vec::fixed< 3 > Omega)
 
arma::vec::fixed< 3 > GetOmega () const
 
- Public Member Functions inherited from SBGATMassProperties
 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
 
- Public Member Functions inherited from SBGATFilter
 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 SBGATPolyhedronGravityModelNew ()
 
static void ComputeSurfacePGM (vtkSmartPointer< vtkPolyData > selected_shape, const std::vector< unsigned int > &queried_elements, bool is_in_meters, double density, const arma::vec::fixed< 3 > &omega, std::vector< double > &slopes, std::vector< double > &inertial_potentials, std::vector< double > &body_fixed_potentials, std::vector< double > &inertial_acc_magnitudes, std::vector< double > &body_fixed_acc_magnitudes)
 
static void SaveSurfacePGM (vtkSmartPointer< vtkPolyData > selected_shape, const std::vector< unsigned int > &queried_elements, bool is_in_meters, const double &mass, const arma::vec::fixed< 3 > &omega, const std::vector< double > &slopes, const std::vector< double > &inertial_potentials, const std::vector< double > &body_fixed_potentials, const std::vector< double > &inertial_acc_magnitudes, const std::vector< double > &body_fixed_acc_magnitudes, std::string path)
 
static void SaveSurfacePGM (vtkSmartPointer< vtkPolyData > selected_shape, const std::vector< unsigned int > &queried_elements, bool is_in_meters, const double &mass, const arma::vec::fixed< 3 > &omega, const std::vector< double > &slopes, const std::vector< double > &inertial_potentials, const std::vector< double > &body_fixed_potentials, const std::vector< double > &inertial_acc_magnitudes, const std::vector< double > &body_fixed_acc_magnitudes, const std::vector< double > &slope_sds, std::string path)
 
static void LoadSurfacePGM (double &mass, arma::vec::fixed< 3 > &omega, std::vector< double > &slopes, std::vector< double > &inertial_potentials, std::vector< double > &body_fixed_potentials, std::vector< double > &inertial_acc_magnitudes, std::vector< double > &body_fixed_acc_magnitudes, std::vector< double > &slope_sds, std::string path)
 
static double GetUe (const arma::vec::fixed< 10 > &Xe)
 
static double GetUf (const arma::vec::fixed< 10 > &Xf)
 
static arma::vec::fixed< 3 > GetAe (const arma::vec::fixed< 10 > &Xe)
 
static arma::vec::fixed< 3 > GetAf (const arma::vec::fixed< 10 > &Xf)
 
- Static Public Member Functions inherited from SBGATMassProperties
static SBGATMassPropertiesNew ()
 
static void ComputeAndSaveMassProperties (vtkSmartPointer< vtkPolyData > shape, std::string path)
 
- Static Public Member Functions inherited from SBGATFilter
static SBGATFilterNew ()
 

Protected Member Functions

 SBGATPolyhedronGravityModel ()
 
 ~SBGATPolyhedronGravityModel () override
 
void Clear ()
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
- Protected Member Functions inherited from SBGATMassProperties
 SBGATMassProperties ()
 
 ~SBGATMassProperties () override
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
- Protected Member Functions inherited from SBGATFilter
 SBGATFilter ()
 
 ~SBGATFilter () override
 
void Clear ()
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 

Protected Attributes

double ** facet_dyads
 
double ** edge_dyads
 
arma::vec::fixed< 3 > Omega
 
- Protected Attributes inherited from SBGATMassProperties
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
 
- Protected Attributes inherited from SBGATFilter
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

 SBGATPolyhedronGravityModel (const SBGATPolyhedronGravityModel &)=delete
 
void operator= (const SBGATPolyhedronGravityModel &)=delete
 

Detailed Description

Evaluation of potential, acceleration caused by a constant-density polyhedron.

Author
Benjamin Bercovici
Date
October 2018

Computes the potential, acceleration caused by a polyhedron of constant density by evaluating the so called Polyhedron Gravity Model as derived by Werner and Scheeres. The input must be a topologically-closed polyhedron. This class will always use results expressed in meters as their distance unit (e.g accelerations in m/s^2, potentials in m^2/s^2,...) . Unit consistency is enforced through the use of the SetScaleMeters() and SetScaleKiloMeters() method.

See Werner, R. A., & Scheeres, D. J. (1997). Exterior gravitation of a polyhedron derived and compared with harmonic and mascon gravitation representations of asteroid 4769 Castalia. Celestial Mechanics and Dynamical Astronomy, 65(3), 313–344. https://doi.org/10.1007/BF00053511 for further details. Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen

Definition at line 44 of file SBGATPolyhedronGravityModel.hpp.

Constructor & Destructor Documentation

◆ SBGATPolyhedronGravityModel() [1/2]

SBGATPolyhedronGravityModel::SBGATPolyhedronGravityModel ( )
protected

Definition at line 73 of file SBGATPolyhedronGravityModel.cpp.

◆ ~SBGATPolyhedronGravityModel()

SBGATPolyhedronGravityModel::~SBGATPolyhedronGravityModel ( )
overrideprotected

Definition at line 81 of file SBGATPolyhedronGravityModel.cpp.

◆ SBGATPolyhedronGravityModel() [2/2]

SBGATPolyhedronGravityModel::SBGATPolyhedronGravityModel ( const SBGATPolyhedronGravityModel )
privatedelete

Member Function Documentation

◆ Clear()

void SBGATPolyhedronGravityModel::Clear ( )
protected

Definition at line 499 of file SBGATPolyhedronGravityModel.cpp.

◆ ComputeSurfacePGM()

void SBGATPolyhedronGravityModel::ComputeSurfacePGM ( vtkSmartPointer< vtkPolyData >  selected_shape,
const std::vector< unsigned int > &  queried_elements,
bool  is_in_meters,
double  density,
const arma::vec::fixed< 3 > &  omega,
std::vector< double > &  slopes,
std::vector< double > &  inertial_potentials,
std::vector< double > &  body_fixed_potentials,
std::vector< double > &  inertial_acc_magnitudes,
std::vector< double > &  body_fixed_acc_magnitudes 
)
static

Evaluates the Polyhedron Gravity Model at the surface of the specified surface elements in the provided shape

Parameters
[in]selected_shapeshape for which the surface polyhedron gravity model must be computed
[in]queried_elementsvector of elements indices where the polyhedron gravity model should be evaluated
[in]is_in_meterstrue if the shape coordinates were expressed in meters, false if they were expressed in kilometers
[in]densityshape bulk density in kg/m^3
[in]omegafixed angular velocity of shape expressed in rad/s
[out]slopesvector storing the gravitational slopes (degrees) evaluated at the center of each queried element
[out]inertial_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[out]body_fixed_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[out]inertial_acc_magnitudesvector storing the inertial gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[out]body_fixed_acc_magnitudesvector storing the body-fixed gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element

Definition at line 540 of file SBGATPolyhedronGravityModel.cpp.

◆ Contains() [1/2]

bool SBGATPolyhedronGravityModel::Contains ( double const *  point,
double  tol = 1e-8 
) const

Determines whether the provided point lies inside or outside the shape

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata
tolerance
Returns
true if the polydata contains the point, false otherwise

Definition at line 224 of file SBGATPolyhedronGravityModel.cpp.

◆ Contains() [2/2]

bool SBGATPolyhedronGravityModel::Contains ( const arma::vec::fixed< 3 > &  point,
double  tol = 1e-8 
) const

Determines whether the provided point lies inside or outside the shape

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata
tolerance
Returns
true if the polydata contains the point, false otherwise

Definition at line 220 of file SBGATPolyhedronGravityModel.cpp.

◆ GetAcceleration() [1/2]

arma::vec::fixed< 3 > SBGATPolyhedronGravityModel::GetAcceleration ( const arma::vec::fixed< 3 > &  point) const

Evaluates the Polyhedron Gravity Model acceleration at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
Returns
PGM acceleration evaluated at the queried point (m / s ^2)

Definition at line 244 of file SBGATPolyhedronGravityModel.cpp.

◆ GetAcceleration() [2/2]

arma::vec::fixed< 3 > SBGATPolyhedronGravityModel::GetAcceleration ( double const *  point) const

Evaluates the Polyhedron Gravity Model acceleration at the specified point assuming a constant density

Parameters
pointpointer to coordinates of queried point, expressed in the same frame as the input polydata
Returns
PGM acceleration evaluated at the queried point (m / s ^2)

Definition at line 248 of file SBGATPolyhedronGravityModel.cpp.

◆ GetAe()

arma::vec::fixed< 3 > SBGATPolyhedronGravityModel::GetAe ( const arma::vec::fixed< 10 > &  Xe)
static

Return the contribution of a specific edge to the acceleration at a specified field point

Parameters
Xevector of dyadic coefficients for edge e at the prescribed fieldpoint
Returns
contribution to the potential of this specific edge at the prescribed fieldpoint

Definition at line 1127 of file SBGATPolyhedronGravityModel.cpp.

◆ GetAf()

arma::vec::fixed< 3 > SBGATPolyhedronGravityModel::GetAf ( const arma::vec::fixed< 10 > &  Xf)
static

Return the contribution of a specific facet to the acceleration at a specified field point

Parameters
Xfvector of dyadic coefficients for facet f at the prescribed fieldpoint
Returns
contribution to the potential of this specific facet at the prescribed fieldpoint

Definition at line 1203 of file SBGATPolyhedronGravityModel.cpp.

◆ GetBodyFixedAccelerationf()

arma::vec::fixed< 3 > SBGATPolyhedronGravityModel::GetBodyFixedAccelerationf ( const int &  f) const

Get the body-fixed acceleration at the center of the specified facet

Parameters
ffacet index
Returns
body-fixed acceleration (m/s^2)

Definition at line 1224 of file SBGATPolyhedronGravityModel.cpp.

◆ GetEeParam()

arma::vec::fixed< 6 > SBGATPolyhedronGravityModel::GetEeParam ( const int &  e) const

Return the parametrization of the designated edge dyad Ee. This dyad is stored in a flattened double container holding nine values and ordered like so

E = { {E[0], E[1], E[2]}, {E[3], E[4], E[5]}, {E[6], E[7], E[8]} }; Due to the symmetric nature of E, its parametrization is {E[0],E[4],E[8],E[1],E[2],E[5]};

Parameters
eedge index
Returns
Ee dyad parametrization

Definition at line 1098 of file SBGATPolyhedronGravityModel.cpp.

◆ GetFfParam()

arma::vec::fixed< 6 > SBGATPolyhedronGravityModel::GetFfParam ( const int &  f) const

Return the parametrization of the designated facet dyad Ff. This dyad is stored in a flattened double container holding nine values and ordered like so

F = { {F[0], F[1], F[2]}, {F[3], F[4], F[5]}, {F[6], F[7], F[8]} }; Due to the symmetric nature of E, its parametrization is {F[0],F[4],F[8],F[1],F[2],F[5]};

Parameters
ffacet index
Returns
Ff dyad parametrization

Definition at line 1163 of file SBGATPolyhedronGravityModel.cpp.

◆ GetGravityGradient()

arma::mat::fixed< 3, 3 > SBGATPolyhedronGravityModel::GetGravityGradient ( const arma::vec::fixed< 3 > &  point) const

Evaluates the Polyhedron Gravity Model gravity gradient matrix at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
[out]gravity_gradient_matPGM gravity gradient matrix evaluated at the queried point (1 / s ^2)

Definition at line 1237 of file SBGATPolyhedronGravityModel.cpp.

◆ GetLe() [1/2]

double SBGATPolyhedronGravityModel::GetLe ( const arma::vec::fixed< 3 > &  pos,
const int &  e 
) const

Return the wire potential of the e-th edge at the specified position

Parameters
posposition of field point
eedge index
Returns
L_e

Definition at line 1052 of file SBGATPolyhedronGravityModel.cpp.

◆ GetLe() [2/2]

double SBGATPolyhedronGravityModel::GetLe ( const double *  pos,
const int &  e 
) const

Return the wire potential of the e-th edge at the specified position

Parameters
posposition of field point
eedge index
Returns
L_e

Definition at line 1057 of file SBGATPolyhedronGravityModel.cpp.

◆ GetMass()

double SBGATPolyhedronGravityModel::GetMass ( ) const
inline

Get polyhedron mass. Must have set the density and call

Returns
mass (kg)

Definition at line 197 of file SBGATPolyhedronGravityModel.hpp.

◆ GetOmega()

arma::vec::fixed<3> SBGATPolyhedronGravityModel::GetOmega ( ) const
inline

Return the angular velocity vector

Returns
Omega angular velocity expressed in body-frame (rad/s)

Definition at line 435 of file SBGATPolyhedronGravityModel.hpp.

◆ GetPerformanceFactor() [1/2]

double SBGATPolyhedronGravityModel::GetPerformanceFactor ( const arma::vec::fixed< 3 > &  pos,
const int &  f 
) const

Return the performance factor of the f-th facet at the specified position

Parameters
posposition of field point
ffacet index
Returns
omega_f

Definition at line 1017 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPerformanceFactor() [2/2]

double SBGATPolyhedronGravityModel::GetPerformanceFactor ( const double *  pos,
const int &  f 
) const

Return the performance factor of the f-th facet at the specified position

Parameters
posposition of field point
ffacet index
Returns
omega_f

Definition at line 1021 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotential() [1/2]

double SBGATPolyhedronGravityModel::GetPotential ( double const *  point) const

Evaluates the Polyhedron Gravity Model potential at the specified point assuming a constant density

Parameters
pointpointer to coordinates of queried point, expressed in the same frame as the polydata
Returns
PGM potential evaluated at the queried point (m ^ 2/ s ^2)

Definition at line 196 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotential() [2/2]

double SBGATPolyhedronGravityModel::GetPotential ( const arma::vec::fixed< 3 > &  point) const

Evaluates the Polyhedron Gravity Model potential at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata
Returns
PGM potential evaluated at the queried point (m ^ 2 / s ^2)

Definition at line 191 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotentialAcceleration() [1/2]

void SBGATPolyhedronGravityModel::GetPotentialAcceleration ( double const *  point,
double &  potential,
arma::vec::fixed< 3 > &  acc 
) const

Evaluates the Polyhedron Gravity Model potential and acceleration at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
[out]potentialPGM potential evaluated at the queried point (m ^ 2 / s ^2)
[out]accPGM acceleration evaluated at the queried point (m / s ^2)

Definition at line 308 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotentialAcceleration() [2/2]

void SBGATPolyhedronGravityModel::GetPotentialAcceleration ( const arma::vec::fixed< 3 > &  point,
double &  potential,
arma::vec::fixed< 3 > &  acc 
) const

Evaluates the Polyhedron Gravity Model potential and acceleration at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
[out]potentialPGM potential evaluated at the queried point (m ^ 2 / s ^2)
[out]accPGM acceleration evaluated at the queried point (m / s ^2)

Definition at line 299 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotentialAccelerationGravityGradient() [1/2]

void SBGATPolyhedronGravityModel::GetPotentialAccelerationGravityGradient ( double const *  point,
double &  potential,
arma::vec::fixed< 3 > &  acc,
arma::mat::fixed< 3, 3 > &  gravity_gradient_mat 
) const

Evaluates the Polyhedron Gravity Model potential, acceleration and gravity gradient matrix at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
[out]potentialPGM potential evaluated at the queried point (m ^ 2 / s ^2)
[out]accPGM acceleration evaluated at the queried point (m / s ^2)
[out]gravity_gradient_matPGM gravity gradient matrix evaluated at the queried point (1 / s ^2)

Definition at line 397 of file SBGATPolyhedronGravityModel.cpp.

◆ GetPotentialAccelerationGravityGradient() [2/2]

void SBGATPolyhedronGravityModel::GetPotentialAccelerationGravityGradient ( const arma::vec::fixed< 3 > &  point,
double &  potential,
arma::vec::fixed< 3 > &  acc,
arma::mat::fixed< 3, 3 > &  gravity_gradient_mat 
) const

Evaluates the Polyhedron Gravity Model potential, acceleration and gravity gradient matrix at the specified point assuming a constant density

Parameters
pointcoordinates of queried point, expressed in the same frame as the polydata used to construct the PGM
[out]potentialPGM potential evaluated at the queried point (m ^ 2 / s ^2)
[out]accPGM acceleration evaluated at the queried point (m / s ^2)
[out]gravity_gradient_matPGM gravity gradient matrix evaluated at the queried point (1 / s ^2)

Definition at line 390 of file SBGATPolyhedronGravityModel.cpp.

◆ GetSlope()

double SBGATPolyhedronGravityModel::GetSlope ( const int &  f) const

Computes the slope as the center of the designated facet

Parameters
ffacet index
Returns
slope (rad)

Definition at line 1232 of file SBGATPolyhedronGravityModel.cpp.

◆ GetUe()

double SBGATPolyhedronGravityModel::GetUe ( const arma::vec::fixed< 10 > &  Xe)
static

Return the contribution of a specific edge to the potential at a specified field point

Parameters
Xevector of dyadic coefficients for edge e at the prescribed fieldpoint
Returns
contribution to the potential of this specific edge at the prescribed fieldpoint

Definition at line 1109 of file SBGATPolyhedronGravityModel.cpp.

◆ GetUf()

double SBGATPolyhedronGravityModel::GetUf ( const arma::vec::fixed< 10 > &  Xf)
static

Return the contribution of a specific facet to the potential at a specified field point

Parameters
Xfvector of dyadic coefficients for facet f at the prescribed fieldpoint
Returns
contribution to the potential of this specific facet at the prescribed fieldpoint

Definition at line 1184 of file SBGATPolyhedronGravityModel.cpp.

◆ GetXe()

arma::vec::fixed< 10 > SBGATPolyhedronGravityModel::GetXe ( const arma::vec::fixed< 3 > &  pos,
const int &  e 
) const

Return the Xe^E vector holding the e-th edge dyadic factors (Le,r_ie_0^T,Ee^T)^T

Parameters
posfield point
eedge index
Returns
Xe^E vector holding the e-th edge dyadic factors

Definition at line 1083 of file SBGATPolyhedronGravityModel.cpp.

◆ GetXf()

arma::vec::fixed< 10 > SBGATPolyhedronGravityModel::GetXf ( const arma::vec::fixed< 3 > &  pos,
const int &  f 
) const

Return the Xf^F vector holding the f-th facet dyadic factors (omega_f,r_if_0^T,Ff^T)^T

Parameters
posfield point
ffacet index
Returns
Xf^F vector holding the f-th facet dyadic factors

Definition at line 1148 of file SBGATPolyhedronGravityModel.cpp.

◆ LoadSurfacePGM()

void SBGATPolyhedronGravityModel::LoadSurfacePGM ( double &  mass,
arma::vec::fixed< 3 > &  omega,
std::vector< double > &  slopes,
std::vector< double > &  inertial_potentials,
std::vector< double > &  body_fixed_potentials,
std::vector< double > &  inertial_acc_magnitudes,
std::vector< double > &  body_fixed_acc_magnitudes,
std::vector< double > &  slope_sds,
std::string  path 
)
static

Loads a previously computed surface Polyhedron Gravity Model from a file

Parameters
[out]massmass of shape model (kg)
[out]omegafixed angular velocity of shape (rad/s)
[out]slopesvector storing the gravitational slopes (degrees) evaluated at the center of each queried element
[out]inertial_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[out]body_fixed_potentialsvector storing the body-fixed gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[out]inertial_acc_magnitudesvector storing the inertial gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[out]body_fixed_acc_magnitudesvector storing the body-fixed gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[out]slope_sdsstandard deviation in the surface slopes (deg)
[in]pathload path (ex: "pgm_surface.json")

Definition at line 889 of file SBGATPolyhedronGravityModel.cpp.

◆ New()

static SBGATPolyhedronGravityModel* SBGATPolyhedronGravityModel::New ( )
static

Constructs with initial values of zero.

◆ operator=()

void SBGATPolyhedronGravityModel::operator= ( const SBGATPolyhedronGravityModel )
privatedelete

◆ PrintHeader()

void SBGATPolyhedronGravityModel::PrintHeader ( std::ostream &  os,
vtkIndent  indent 
)
override

Definition at line 491 of file SBGATPolyhedronGravityModel.cpp.

◆ PrintSelf()

void SBGATPolyhedronGravityModel::PrintSelf ( std::ostream &  os,
vtkIndent  indent 
)
override

Definition at line 523 of file SBGATPolyhedronGravityModel.cpp.

◆ PrintTrailer()

void SBGATPolyhedronGravityModel::PrintTrailer ( std::ostream &  os,
vtkIndent  indent 
)
override

Definition at line 494 of file SBGATPolyhedronGravityModel.cpp.

◆ RequestData()

int SBGATPolyhedronGravityModel::RequestData ( vtkInformation *  request,
vtkInformationVector **  inputVector,
vtkInformationVector *  outputVector 
)
overrideprotected

Definition at line 90 of file SBGATPolyhedronGravityModel.cpp.

◆ SaveSurfacePGM() [1/2]

void SBGATPolyhedronGravityModel::SaveSurfacePGM ( vtkSmartPointer< vtkPolyData >  selected_shape,
const std::vector< unsigned int > &  queried_elements,
bool  is_in_meters,
const double &  mass,
const arma::vec::fixed< 3 > &  omega,
const std::vector< double > &  slopes,
const std::vector< double > &  inertial_potentials,
const std::vector< double > &  body_fixed_potentials,
const std::vector< double > &  inertial_acc_magnitudes,
const std::vector< double > &  body_fixed_acc_magnitudes,
std::string  path 
)
static

Saves the provided Polyhedron Gravity Model to a file

Parameters
[in]selected_shapeshape for which the surface polyhedron gravity model must be computed
[in]queried_elementsshape indices of elements where the polyhedron gravity model should be evaluated
[in]is_in_meterstrue if the shape coordinates were expressed in meters, false if they were expressed in kilometers
[in]massmass of shape model (kg)
[in]omegafixed angular velocity of shape (rad/s)
[in]slopesvector storing the gravitational slopes (degrees) evaluated at the center of each queried element
[in]inertial_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[in]body_fixed_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[in]inertial_acc_magnitudesvector storing the inertial gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[in]body_fixed_acc_magnitudesvector storing the body-fixed gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[in]pathsave path (ex: "pgm_surface.json")

Definition at line 630 of file SBGATPolyhedronGravityModel.cpp.

◆ SaveSurfacePGM() [2/2]

void SBGATPolyhedronGravityModel::SaveSurfacePGM ( vtkSmartPointer< vtkPolyData >  selected_shape,
const std::vector< unsigned int > &  queried_elements,
bool  is_in_meters,
const double &  mass,
const arma::vec::fixed< 3 > &  omega,
const std::vector< double > &  slopes,
const std::vector< double > &  inertial_potentials,
const std::vector< double > &  body_fixed_potentials,
const std::vector< double > &  inertial_acc_magnitudes,
const std::vector< double > &  body_fixed_acc_magnitudes,
const std::vector< double > &  slope_sds,
std::string  path 
)
static

Saves the provided Polyhedron Gravity Model to a file

Parameters
[in]selected_shapeshape for which the surface polyhedron gravity model must be computed
[in]queried_elementsshape indices of elements where the polyhedron gravity model should be evaluated
[in]is_in_meterstrue if the shape coordinates were expressed in meters, false if they were expressed in kilometers
[in]massmass of shape model (kg)
[in]omegafixed angular velocity of shape (rad/s)
[in]slopesvector storing the gravitational slopes (degrees) evaluated at the center of each queried element
[in]inertial_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[in]body_fixed_potentialsvector storing the inertial gravitational potentials (m^2/s^2) evaluated at the center of each queried element
[in]inertial_acc_magnitudesvector storing the inertial gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[in]body_fixed_acc_magnitudesvector storing the body-fixed gravitational acceleration magnitudes (m/s^2) evaluated at the center of each queried element
[in]slope_sdsvector storing the standard deviation in the gravitational slopes (degrees) evaluated at the center of each queried element
[in]pathsave path (ex: "pgm_surface.json")

Definition at line 748 of file SBGATPolyhedronGravityModel.cpp.

◆ SetOmega()

void SBGATPolyhedronGravityModel::SetOmega ( arma::vec::fixed< 3 >  Omega)
inline

Set the angular velocity vector

Parameters
Omegaangular velocity expressed in body-frame (rad/s)

Definition at line 429 of file SBGATPolyhedronGravityModel.hpp.

◆ vtkTypeMacro()

SBGATPolyhedronGravityModel::vtkTypeMacro ( SBGATPolyhedronGravityModel  ,
vtkPolyDataAlgorithm   
)

Member Data Documentation

◆ edge_dyads

double** SBGATPolyhedronGravityModel::edge_dyads
protected

Definition at line 448 of file SBGATPolyhedronGravityModel.hpp.

◆ facet_dyads

double** SBGATPolyhedronGravityModel::facet_dyads
protected

Definition at line 447 of file SBGATPolyhedronGravityModel.hpp.

◆ Omega

arma::vec::fixed<3> SBGATPolyhedronGravityModel::Omega
protected

Definition at line 450 of file SBGATPolyhedronGravityModel.hpp.


The documentation for this class was generated from the following files: