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

Defines the SBGATObsLightcurve class. More...

#include <SBGATObsLightcurve.hpp>

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

Public Member Functions

 vtkTypeMacro (SBGATObsLightcurve, SBGATObs)
 
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 CollectMeasurements (std::vector< std::array< double, 2 > > &measurements, const double &time, const int &N, const arma::vec &sun_dir, const arma::vec &observer_dir, const std::vector< arma::vec > &positions_vec, const std::vector< arma::vec > &velocities_vec, const std::vector< arma::vec > &mrps_vec, const std::vector< arma::vec > &omegas_vec, const bool &penalize_indicence)
 
void SaveLightCurveData (const std::vector< std::array< double, 2 > > &measurements, std::string savepath)
 
- Public Member Functions inherited from SBGATObs
 vtkTypeMacro (SBGATObs, 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 ()
 
int get_number_of_bodies () const
 

Static Public Member Functions

static SBGATObsLightcurveNew ()
 
- Static Public Member Functions inherited from SBGATObs
static SBGATObsNew ()
 

Protected Member Functions

 SBGATObsLightcurve ()
 
 ~SBGATObsLightcurve () override
 
void reverse_ray_trace (std::array< double, 2 > &measurements_temp, const std::vector< std::vector< int > > &facets_in_view, const arma::vec &sun_dir, const arma::vec &observer_dir, const int N, const bool penalize_indicence, const std::vector< arma::mat > &BN_dcms_vec, const std::vector< arma::vec > &positions_vec)
 
- Protected Member Functions inherited from SBGATObs
 SBGATObs ()
 
 ~SBGATObs () override
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
int FillInputPortInformation (int port, vtkInformation *info) VTK_OVERRIDE
 
void prefind_facets_inview (std::vector< int > &facets_in_view, const unsigned int &body_index, const std::vector< arma::vec > &dir_to_check_vec, const std::vector< arma::mat > &BN_dcms_vec, const std::vector< arma::vec > &positions_vec)
 
void find_min_facet_surface_area ()
 
bool check_line_for_intersect (const int &origin_body_index, const arma::vec &start_point_origin_body, const arma::vec &end_point_inertial, const std::vector< arma::mat > &BN_dcms_vec, const std::vector< arma::vec > &positions_vec, const double &tol) const
 

Private Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from SBGATObs
std::vector< vtkSmartPointer< vtkModifiedBSPTree > > bspTree_vec
 
std::vector< vtkPolyData * > polydata_vec
 
std::vector< arma::vec > center_of_mass_vec
 
double scaleFactor = 1
 
double min_area
 
int number_of_bodies
 

Detailed Description

Defines the SBGATObsLightcurve class.

Author
Benjamin Bercovici
Jay McMahon
Date
October 2018

Defines the SBGATObsLightcurve class, used to generate lightcurves over a single/binary asteroid system at a fixed phase angle

Definition at line 45 of file SBGATObsLightcurve.hpp.

Constructor & Destructor Documentation

◆ SBGATObsLightcurve() [1/2]

SBGATObsLightcurve::SBGATObsLightcurve ( )
protected

Definition at line 41 of file SBGATObsLightcurve.cpp.

◆ ~SBGATObsLightcurve()

SBGATObsLightcurve::~SBGATObsLightcurve ( )
overrideprotected

Definition at line 46 of file SBGATObsLightcurve.cpp.

◆ SBGATObsLightcurve() [2/2]

SBGATObsLightcurve::SBGATObsLightcurve ( const SBGATObsLightcurve )
privatedelete

Member Function Documentation

◆ CollectMeasurements()

void SBGATObsLightcurve::CollectMeasurements ( std::vector< std::array< double, 2 > > &  measurements,
const double &  time,
const int &  N,
const arma::vec &  sun_dir,
const arma::vec &  observer_dir,
const std::vector< arma::vec > &  positions_vec,
const std::vector< arma::vec > &  velocities_vec,
const std::vector< arma::vec > &  mrps_vec,
const std::vector< arma::vec > &  omegas_vec,
const bool &  penalize_indicence 
)

Computes collected luminosity over the surface of the considered small bodies at the specified time after epoch. Exposure is instantaneous. The luminosity is computed as the number of sample points in view of both the sun and the observer at each time (the "hit count"), normalized by the largest hit count in the observation sequence. The positions of each body is specified. The attitude of each body is supposed to derive from a constant-spin rotation such that each of the dcms BN is equal to identity when dt == 0.

Parameters
measurementsreference to a vector of std::arrays holding (times,luminosity)
timemeasurement time
Nminimum number of measurements to produce over the smallest facet in the shape. The number of samples for any other facet will be equal to N * facet_surface_area / smallest_facet_surface_area
sun_dirunit direction towards sun from target in inertial frame
observer_dirunit direction towards observer from target in inertial frame
positions_vecvector of positions of each target's center-of-mass expressed in the primary's body frame
velocities_vecvector of inertial velocities of each target's center-of-mass expressed in the primary's body frame
mrps_vecvector of MRPs defining the inertial-to-body DCM [BN]
omegas_vecvector of angular velocities of each body expressed in the inertial frame
penalize_incidenceif true, each measurement will be weighed by the cos(incidence) angle between the sampled point and the radar squared. If false, all measurements (in view of the radar and not blocked) are weighed equally have the same weight

Definition at line 51 of file SBGATObsLightcurve.cpp.

◆ New()

static SBGATObsLightcurve* SBGATObsLightcurve::New ( )
static

◆ operator=()

void SBGATObsLightcurve::operator= ( const SBGATObsLightcurve )
privatedelete

◆ PrintHeader()

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

Definition at line 140 of file SBGATObsLightcurve.cpp.

◆ PrintSelf()

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

Definition at line 151 of file SBGATObsLightcurve.cpp.

◆ PrintTrailer()

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

Definition at line 143 of file SBGATObsLightcurve.cpp.

◆ reverse_ray_trace()

void SBGATObsLightcurve::reverse_ray_trace ( std::array< double, 2 > &  measurements_temp,
const std::vector< std::vector< int > > &  facets_in_view,
const arma::vec &  sun_dir,
const arma::vec &  observer_dir,
const int  N,
const bool  penalize_indicence,
const std::vector< arma::mat > &  BN_dcms_vec,
const std::vector< arma::vec > &  positions_vec 
)
protected

Ray traces all of the facets in view to the sun/observer and increment measurement counter if in view

Parameters
measurements_tempreference to an std::array holding (times,luminosity)
facets_in_viewreference to a vector of vector holding indices of (maybe) illuminated facets for all considered bodies
sun_dirsun direction expressed in inertial frame
observer_dirobserver direction expressed in inertial frame
BN_dcms_vecvector holding the DCMs orienting the body frame of each body w/r to inertial
positions_vecvector holding the position vector of the CM of each body w/r to the primary

Definition at line 161 of file SBGATObsLightcurve.cpp.

◆ SaveLightCurveData()

void SBGATObsLightcurve::SaveLightCurveData ( const std::vector< std::array< double, 2 > > &  measurements,
std::string  savepath 
)

Save the raw data points from the unnormalized (times,luminosity) time series to file

Parameters
measurementsreference to a vector of std::arrays holding (times,luminosity)
savepathpath to file where raw lightcurve data will be saved (ex: "output/lightcurve.txt")

Definition at line 124 of file SBGATObsLightcurve.cpp.

◆ vtkTypeMacro()

SBGATObsLightcurve::vtkTypeMacro ( SBGATObsLightcurve  ,
SBGATObs   
)

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