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

Evaluation of the formal uncertainty in the volume, center of mass, inertia tensor parametrization from a topologically-closed, constant density polyhedron. More...

#include <SBGATMassPropertiesUQ.hpp>

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

Public Member Functions

virtual void SetModel (vtkSmartPointer< SBGATFilter > model)
 
const arma::mat & GetPartialComPartialC () const
 
const arma::mat & GetPartialIPartialC () const
 
const arma::mat & GetPartialSigmaPartialC () const
 
const arma::rowvec & GetPartialVolumePartialC () const
 
virtual void ApplyDeviation (const arma::vec &delta_C)
 
void PrecomputeMassPropertiesPartials ()
 
const arma::mat & GetPartialUnitDensityMomentsPartialI () const
 
const arma::mat::fixed< 3, 6 > & GetPartialSigmaPartialI () const
 
- Public Member Functions inherited from SBGATFilterUQ
vtkSmartPointer< SBGATFilterGetModel ()
 
arma::mat GetCovarianceSquareRoot (std::string method="chol") const
 
void ComputeVerticesCovarianceGlobal (const double &standard_dev, const double &correl_distance)
 
void AddRadialUncertaintyRegionToCovariance (int region_center_index, const double &standard_dev, const double &correl_distance)
 
int RegularizeCovariance ()
 
void AddNormalUncertaintyRegionToCovariance (int region_center_index, const double &standard_dev, const double &correl_distance)
 
void SetCovarianceComponent (const arma::mat::fixed< 3, 3 > &P, const int &v0, const int &v1)
 
void SetCovariance (const arma::mat &P_CC)
 
void SaveVerticesCovariance (std::string path) const
 
arma::mat GetVerticesCovariance () const
 
void SaveNonZeroVerticesCovariance (std::string path) const
 
int LoadVerticesCovarianceFromJson (std::string path)
 
arma::sp_mat PartialTfPartialC (const int &f) const
 
arma::mat::fixed< 3, 9 > PartialNfPartialTf (const int &f) const
 
void TakeAndSaveSlice (int axis, std::string path, const double &c) const
 

Static Public Member Functions

static void TestPartials (std::string input, double tol, bool shape_in_meters)
 
static void RunMCUQVolumeCOMInertia (std::string path_to_shape, const double &density, const bool &shape_in_meters, const arma::mat &C_CC, const unsigned int &N_samples, std::string output_dir, int N_saved_shapes, arma::mat &deviations, arma::vec &all_volumes, arma::mat &all_com, arma::mat &all_inertia)
 
- Static Public Member Functions inherited from SBGATFilterUQ
static arma::mat GetCovarianceSquareRoot (arma::mat P_CC, std::string method)
 
static double KLDivergence (const arma::vec &m0, const arma::mat &m1, const arma::mat &P0, const arma::mat &P1)
 

Protected Member Functions

arma::mat::fixed< 3, 6 > PartialSigmaPartialI () const
 
arma::mat::fixed< 3, 6 > PartialUnitDensityMomentsPartialI () const
 
arma::rowvec::fixed< 9 > PartialEqDeltaIfErPartialTf (const int &f, const int &q, const int &r, const arma::vec::fixed< 9 > &Tf) const
 
arma::mat::fixed< 3, 9 > PartialDeltaComPartialTf (const int &f) const
 
arma::rowvec::fixed< 9 > PartialEqDeltaIOverDeltaVfErPartialTf (const arma::vec::fixed< 3 > &e_q, const arma::vec::fixed< 3 > &e_r, const arma::vec::fixed< 9 > &Tf) const
 
arma::rowvec::fixed< 9 > PartialDeltaVfPartialTf (const int &f) const
 
arma::mat::fixed< 6, 9 > PartialDeltaIfPartialTf (const int &f) const
 
arma::mat::fixed< 6, 9 > PartialDeltaIOverDeltaVPartialTf (const int &f) const
 
virtual void ApplyTfDeviation (arma::vec::fixed< 9 > delta_Tf, const int &f)
 
- Protected Member Functions inherited from SBGATFilterUQ
void SaveSlice (int axis, std::string path, const std::vector< std::vector< arma::vec > > &lines) const
 
void TakeSlice (int axis, std::vector< std::vector< arma::vec > > &lines, const double &c) const
 

Static Protected Member Functions

static arma::mat::fixed< 3, 9 > PartialDeltaCMfPartialTf ()
 
static void TestPartialDeltaVfPartialTf (std::string input, double tol, bool shape_in_meters)
 
static void TestPartialDeltaIOverDeltaVPartialTf (std::string input, double tol, bool shape_in_meters)
 
static void TestPartialDeltaIfPartialTf (std::string input, double tol, bool shape_in_meters)
 
static void TestPartialDeltaVPartialC (std::string input, double tol, bool shape_in_meters)
 
static void TestGetPartialVolumePartialC (std::string input, double tol, bool shape_in_meters)
 
static void TestGetPartialComPartialC (std::string input, double tol, bool shape_in_meters)
 
static void TestGetPartialIPartialC (std::string input, double tol, bool shape_in_meters)
 
static void TestGetPartialAllInertiaPartialC (std::string input, double tol, bool shape_in_meters)
 
static void TestPartialEqDeltaIfErPartialTf (std::string input, double tol, bool shape_in_meters)
 
static void TestGetPartialSigmaPartialC (std::string input, double tol, bool shape_in_meters)
 

Protected Attributes

arma::rowvec precomputed_partialVpartialC
 
arma::mat precomputed_partialGpartialC
 
arma::mat precomputed_partialSigmapartialC
 
arma::mat precomputed_partialIpartialC
 
arma::mat::fixed< 3, 6 > precomputed_partialUnitDensityMomentsPartialI
 
arma::mat::fixed< 3, 6 > precomputed_partialSigmapartialI
 
- Protected Attributes inherited from SBGATFilterUQ
vtkSmartPointer< SBGATFiltermodel
 
arma::mat P_CC
 
arma::sp_mat P_CC_sparse
 

Detailed Description

Evaluation of the formal uncertainty in the volume, center of mass, inertia tensor parametrization from a topologically-closed, constant density polyhedron.

Author
Benjamin Bercovici
Date
January 2019

Definition at line 22 of file SBGATMassPropertiesUQ.hpp.

Member Function Documentation

◆ ApplyDeviation()

void SBGATMassPropertiesUQ::ApplyDeviation ( const arma::vec &  delta_C)
virtual

Applies prescribed deviation to all the N_vertices control points and updates model

Parameters
delta_Cdeviation (3 * N_vertices x 1)

Implements SBGATFilterUQ.

Reimplemented in SBGATPolyhedronGravityModelUQ.

Definition at line 852 of file SBGATMassPropertiesUQ.cpp.

◆ ApplyTfDeviation()

void SBGATMassPropertiesUQ::ApplyTfDeviation ( arma::vec::fixed< 9 >  delta_Tf,
const int &  f 
)
protectedvirtual

Applies deviation to the coordinates of the vertices in the prescribed facet and updates the pgm

Parameters
delta_Tfdeviation
ffacet index

Implements SBGATFilterUQ.

Reimplemented in SBGATPolyhedronGravityModelUQ.

Definition at line 293 of file SBGATMassPropertiesUQ.cpp.

◆ GetPartialComPartialC()

const arma::mat& SBGATMassPropertiesUQ::GetPartialComPartialC ( ) const
inline

Return the partial derivative of the shape's center of mass with respect to the shape's vertices coordinates

Returns
partial derivative of the shape's center of mass with respect to the shape's vertices

Definition at line 47 of file SBGATMassPropertiesUQ.hpp.

◆ GetPartialIPartialC()

const arma::mat& SBGATMassPropertiesUQ::GetPartialIPartialC ( ) const
inline

Return the partial derivative of the 6 unique components of the inertia tensor {I(0,0),I(1,1),I(2,2),I(0,1),I(0,2),I(1,2)} with respect to the shape coordinates

Returns
partial derivative

Definition at line 55 of file SBGATMassPropertiesUQ.hpp.

◆ GetPartialSigmaPartialC()

const arma::mat& SBGATMassPropertiesUQ::GetPartialSigmaPartialC ( ) const
inline

Return the partial derivative of the MRP orienting the body-frame (B) to principal-frame (P) dcm (PB) with respect to the shape vertices coordinates

Returns
partial derivative

Definition at line 63 of file SBGATMassPropertiesUQ.hpp.

◆ GetPartialSigmaPartialI()

const arma::mat::fixed<3,6>& SBGATMassPropertiesUQ::GetPartialSigmaPartialI ( ) const
inline

Return the partial derivative of the MRP orienting the body-frame (B) to principal-frame (P) dcm (PB) with respect to the inertia tensor parametrization

Returns
partial derivative

Definition at line 99 of file SBGATMassPropertiesUQ.hpp.

◆ GetPartialUnitDensityMomentsPartialI()

const arma::mat& SBGATMassPropertiesUQ::GetPartialUnitDensityMomentsPartialI ( ) const
inline

Return the partial derivative of the unit density moments relative to a change in the inertia tensor parametrization

Returns
partial derivative of the unit density moments relative to a change in the inertia tensor parametrization

Definition at line 91 of file SBGATMassPropertiesUQ.hpp.

◆ GetPartialVolumePartialC()

const arma::rowvec& SBGATMassPropertiesUQ::GetPartialVolumePartialC ( ) const
inline

Return the partial derivative of the volume with respect to the shape coordinates

Returns
partial derivative

Definition at line 71 of file SBGATMassPropertiesUQ.hpp.

◆ PartialDeltaCMfPartialTf()

arma::mat::fixed< 3, 9 > SBGATMassPropertiesUQ::PartialDeltaCMfPartialTf ( )
staticprotected

Return the partial derivative of the center of mass of the considered tetrahedron with respect to the facet coordinates

Returns
partial derivative of center of mass with respect to facet coordinates

Definition at line 37 of file SBGATMassPropertiesUQ.cpp.

◆ PartialDeltaComPartialTf()

arma::mat::fixed< 3, 9 > SBGATMassPropertiesUQ::PartialDeltaComPartialTf ( const int &  f) const
protected

Return the partial derivative of the shape's center-of-mass with respect to the f-th facet coordinates

Parameters
ffacet index
Returns
partial derivative of shape's center-of-mass with respect to facet coordinates

Definition at line 644 of file SBGATMassPropertiesUQ.cpp.

◆ PartialDeltaIfPartialTf()

arma::mat::fixed< 6, 9 > SBGATMassPropertiesUQ::PartialDeltaIfPartialTf ( const int &  f) const
protected

Return the partial derivative of the tetrahedron's inertia tensor parametrization relative to the facet coordinates

Parameters
ffacet index
Returns
partial derivative of the tetrahedron's inertia tensor parametrization relative to the facet coordinates

Definition at line 725 of file SBGATMassPropertiesUQ.cpp.

◆ PartialDeltaIOverDeltaVPartialTf()

arma::mat::fixed<6,9> SBGATMassPropertiesUQ::PartialDeltaIOverDeltaVPartialTf ( const int &  f) const
protected

Return the partial derivative of a tetrahedron's inertia-times-volume tensor parametrization with respect to the subtending facet's vertices coordinates

Parameters
ffacet index
Returns
partial derivative of the tetrahedron's inertia-times-volume tensor parametrization with respect to the subtending facet's vertices coordinates

◆ PartialDeltaVfPartialTf()

arma::rowvec::fixed< 9 > SBGATMassPropertiesUQ::PartialDeltaVfPartialTf ( const int &  f) const
protected

Return the partial derivative of the volume of the tetrahedron subtended by facet f with respect to the facet coordinates

Parameters
ffacet index
Returns
partial derivative of tetrahedron volume with respect to facet coordinates

Definition at line 16 of file SBGATMassPropertiesUQ.cpp.

◆ PartialEqDeltaIfErPartialTf()

arma::rowvec::fixed< 9 > SBGATMassPropertiesUQ::PartialEqDeltaIfErPartialTf ( const int &  f,
const int &  q,
const int &  r,
const arma::vec::fixed< 9 > &  Tf 
) const
protected

Return the partial derivative of the (q,r) component of the contribution of the f-facet to the shape's inertia tensor, relative to the f-facet vertices coordinates

Parameters
ffacet index
qfirst index
rsecond index
Tf9x1 vector holding coordinates of vertices in facet EXPRESSED IN METERS
Returns
partial derivative of the (q,r) component of the contribution of the f-facet to the shape's inertia tensor

Definition at line 661 of file SBGATMassPropertiesUQ.cpp.

◆ PartialEqDeltaIOverDeltaVfErPartialTf()

arma::rowvec::fixed< 9 > SBGATMassPropertiesUQ::PartialEqDeltaIOverDeltaVfErPartialTf ( const arma::vec::fixed< 3 > &  e_q,
const arma::vec::fixed< 3 > &  e_r,
const arma::vec::fixed< 9 > &  Tf 
) const
protected

Return the partial derivative of e_q.T * DeltaIOverDeltaVfEr * e_r with respect to the f-facet vertices coordinates

Parameters
e_qfirst 3x1 vector canonical unit vector
e_rfirst 3x1 vector canonical unit vector
Tf9x1 vector holding coordinates of vertices in facet EXPRESSED IN METERS
Returns
the partial derivative of e_q.T * DeltaIOverDeltaVfEr * e_r with respect to the f-facet vertices coordinates

Definition at line 688 of file SBGATMassPropertiesUQ.cpp.

◆ PartialSigmaPartialI()

arma::mat::fixed< 3, 6 > SBGATMassPropertiesUQ::PartialSigmaPartialI ( ) const
protected

Returns the partial derivative of the MRP orienting the principal axes with respect to the parametrization of the unit-density inertia tensor

Returns
partial derivative

Definition at line 880 of file SBGATMassPropertiesUQ.cpp.

◆ PartialUnitDensityMomentsPartialI()

arma::mat::fixed< 3, 6 > SBGATMassPropertiesUQ::PartialUnitDensityMomentsPartialI ( ) const
protected

Returns partial derivative of unit-density inertia moments relative to the parametrization of the unit-density inertia tensor

Definition at line 965 of file SBGATMassPropertiesUQ.cpp.

◆ PrecomputeMassPropertiesPartials()

void SBGATMassPropertiesUQ::PrecomputeMassPropertiesPartials ( )

Evaluates the partial of the volume, center of mass and mrp orienting the principal axes relative to the vertices coordinates and stores the computed partials in designated containers

Definition at line 1022 of file SBGATMassPropertiesUQ.cpp.

◆ RunMCUQVolumeCOMInertia()

void SBGATMassPropertiesUQ::RunMCUQVolumeCOMInertia ( std::string  path_to_shape,
const double &  density,
const bool &  shape_in_meters,
const arma::mat &  C_CC,
const unsigned int &  N_samples,
std::string  output_dir,
int  N_saved_shapes,
arma::mat &  deviations,
arma::vec &  all_volumes,
arma::mat &  all_com,
arma::mat &  all_inertia 
)
static

Runs a Monte Carlo on the shape and the volume, center-of-mass and inertia tensor parametrizaton

Parameters
[in]path_to_shapepath to reference shape
[in]densitysmall body density in kg/m^3
[in]shape_in_meterstrue if reference shape has its units expressed in meters, false otherwise
[in]C_CCsquare root of the covariance of the shape vertices coordinates. Must be of dimensions (3N_C * 3N_C) where N_C is the number of vertices in the reference shape
[in]N_samplesnumber of shape outcomes to draw
[in]output_dirpath ending in "/" where to save shape-related monte-carlo data. Only used if last argument is larger than 0
[in]N_saved_shapesnumber of shape outcomes to save. must be lesser or equal than N_samples
[out]deviationsholds N_samples 3*N_C column vectors storing the deviation applied to the coordinates of the reference shape at every sample
[out]all_volumesholds N_samples of the volume
[out]all_comholds N_samples of the center-of-mass
[out]all_inertiaholds N_samples of the inertia tensor parametrization

Definition at line 1058 of file SBGATMassPropertiesUQ.cpp.

◆ SetModel()

virtual void SBGATMassPropertiesUQ::SetModel ( vtkSmartPointer< SBGATFilter model)
inlinevirtual

Sets the model associated to this uncertainty quantification container and updates the partials of the mass properties relative to the shape

Parameters
[in]pgmpointer to valid SBGATFilter
[in]

Implements SBGATFilterUQ.

Definition at line 32 of file SBGATMassPropertiesUQ.hpp.

◆ TestGetPartialAllInertiaPartialC()

void SBGATMassPropertiesUQ::TestGetPartialAllInertiaPartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 473 of file SBGATMassPropertiesUQ.cpp.

◆ TestGetPartialComPartialC()

void SBGATMassPropertiesUQ::TestGetPartialComPartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 223 of file SBGATMassPropertiesUQ.cpp.

◆ TestGetPartialIPartialC()

void SBGATMassPropertiesUQ::TestGetPartialIPartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 334 of file SBGATMassPropertiesUQ.cpp.

◆ TestGetPartialSigmaPartialC()

void SBGATMassPropertiesUQ::TestGetPartialSigmaPartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 405 of file SBGATMassPropertiesUQ.cpp.

◆ TestGetPartialVolumePartialC()

void SBGATMassPropertiesUQ::TestGetPartialVolumePartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 153 of file SBGATMassPropertiesUQ.cpp.

◆ TestPartialDeltaIfPartialTf()

void SBGATMassPropertiesUQ::TestPartialDeltaIfPartialTf ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 570 of file SBGATMassPropertiesUQ.cpp.

◆ TestPartialDeltaIOverDeltaVPartialTf()

static void SBGATMassPropertiesUQ::TestPartialDeltaIOverDeltaVPartialTf ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

◆ TestPartialDeltaVfPartialTf()

void SBGATMassPropertiesUQ::TestPartialDeltaVfPartialTf ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 69 of file SBGATMassPropertiesUQ.cpp.

◆ TestPartialDeltaVPartialC()

static void SBGATMassPropertiesUQ::TestPartialDeltaVPartialC ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

◆ TestPartialEqDeltaIfErPartialTf()

void SBGATMassPropertiesUQ::TestPartialEqDeltaIfErPartialTf ( std::string  input,
double  tol,
bool  shape_in_meters 
)
staticprotected

Definition at line 748 of file SBGATMassPropertiesUQ.cpp.

◆ TestPartials()

void SBGATMassPropertiesUQ::TestPartials ( std::string  input,
double  tol,
bool  shape_in_meters 
)
static

Runs a finite-differencing based test of the implemented PGM partials

Parameters
inputpath to obj file used to test the partials
tolrelative tolerance
shape_in_meterstrue if tested shape has its coordinates expressed in meters

Definition at line 48 of file SBGATMassPropertiesUQ.cpp.

Member Data Documentation

◆ precomputed_partialGpartialC

arma::mat SBGATMassPropertiesUQ::precomputed_partialGpartialC
protected

Definition at line 249 of file SBGATMassPropertiesUQ.hpp.

◆ precomputed_partialIpartialC

arma::mat SBGATMassPropertiesUQ::precomputed_partialIpartialC
protected

Definition at line 251 of file SBGATMassPropertiesUQ.hpp.

◆ precomputed_partialSigmapartialC

arma::mat SBGATMassPropertiesUQ::precomputed_partialSigmapartialC
protected

Definition at line 250 of file SBGATMassPropertiesUQ.hpp.

◆ precomputed_partialSigmapartialI

arma::mat::fixed<3,6> SBGATMassPropertiesUQ::precomputed_partialSigmapartialI
protected

Definition at line 253 of file SBGATMassPropertiesUQ.hpp.

◆ precomputed_partialUnitDensityMomentsPartialI

arma::mat::fixed<3,6> SBGATMassPropertiesUQ::precomputed_partialUnitDensityMomentsPartialI
protected

Definition at line 252 of file SBGATMassPropertiesUQ.hpp.

◆ precomputed_partialVpartialC

arma::rowvec SBGATMassPropertiesUQ::precomputed_partialVpartialC
protected

Definition at line 248 of file SBGATMassPropertiesUQ.hpp.


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