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
SBGATObsRadar Class Reference

Computes range/range-rate Doppler images over the surface of provided small body. More...

#include <SBGATObsRadar.hpp>

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

Public Member Functions

 vtkTypeMacro (SBGATObsRadar, 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 CollectMeasurements (SBGATRadarObsSequence &measurements_sequence, const double &time, const int &N, const arma::vec &radar_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 BinObservations (const SBGATRadarObsSequence &measurements_sequence, const double &r_bin, const double &rr_bin)
 
void SaveImages (std::string savepath)
 
std::vector< vtkSmartPointer< vtkImageData > > GetImages () const
 
void ClearImages ()
 
- 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 SBGATObsRadarNew ()
 
- Static Public Member Functions inherited from SBGATObs
static SBGATObsNew ()
 

Protected Member Functions

 SBGATObsRadar ()
 
 ~SBGATObsRadar () override
 
void reverse_ray_trace (SBGATRadarObsSequence &measurements_sequence, const std::vector< std::vector< int > > &facets_in_view, const arma::vec &radar_dir, const int N, const bool penalize_indicence, const std::vector< arma::mat > &BN_dcms_vec, const std::vector< arma::vec > &positions_vec, const std::vector< arma::vec > &velocities_vec, const std::vector< arma::vec > &omega_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
 

Protected Attributes

std::vector< vtkSmartPointer< vtkImageData > > images
 
double max_value
 
- 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
 

Private Member Functions

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

Detailed Description

Computes range/range-rate Doppler images over the surface of provided small body.

Author
Benjamin Bercovici
Jay McMahon
Date
October 2018

Computes range/range-rate Doppler images over the surface of provided small body

Definition at line 56 of file SBGATObsRadar.hpp.

Constructor & Destructor Documentation

◆ SBGATObsRadar() [1/2]

SBGATObsRadar::SBGATObsRadar ( )
protected

Definition at line 60 of file SBGATObsRadar.cpp.

◆ ~SBGATObsRadar()

SBGATObsRadar::~SBGATObsRadar ( )
overrideprotected

Definition at line 64 of file SBGATObsRadar.cpp.

◆ SBGATObsRadar() [2/2]

SBGATObsRadar::SBGATObsRadar ( const SBGATObsRadar )
privatedelete

Member Function Documentation

◆ BinObservations()

void SBGATObsRadar::BinObservations ( const SBGATRadarObsSequence measurements_sequence,
const double &  r_bin,
const double &  rr_bin 
)

Bins the provided measurements sequence into a series of 2d-histogram Will throw an std::runtime_error exception if

  • either of the provided bin sizes are invalid (i.e <= 0)
  • the provided bin sizes yield empty histogram dimensions
    Parameters
    measurements_sequencereference to MeasurementsSequence, holding collected range/range-rate measurements at each observation time
    r_binrange bin size (m)
    rr_binrange-rate bin size (m/s)

Definition at line 225 of file SBGATObsRadar.cpp.

◆ ClearImages()

void SBGATObsRadar::ClearImages ( )
inline

Clears images, if any

Definition at line 128 of file SBGATObsRadar.hpp.

◆ CollectMeasurements()

void SBGATObsRadar::CollectMeasurements ( SBGATRadarObsSequence measurements_sequence,
const double &  time,
const int &  N,
const arma::vec &  radar_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 
)

Collects range/range-rate samples over the surface of the specified system of small bodies at a given time after the epoch. The radar source is positionned at 1e6 * l meters from the target's center of mass, where l is a measure of the first object's diagonal.

Parameters
measurements_sequencereference to MeasurementsSequence, holding collected range/range-rate measurements at each observation time
timeobservation timestamp
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
radar_dirunit direction towards radar 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 69 of file SBGATObsRadar.cpp.

◆ GetImages()

std::vector< vtkSmartPointer< vtkImageData > > SBGATObsRadar::GetImages ( ) const

Return vector holding the computed images

Returns
vector of radar images

Definition at line 352 of file SBGATObsRadar.cpp.

◆ New()

static SBGATObsRadar* SBGATObsRadar::New ( )
static

Constructs with initial values of zero.

◆ operator=()

void SBGATObsRadar::operator= ( const SBGATObsRadar )
privatedelete

◆ PrintHeader()

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

Definition at line 357 of file SBGATObsRadar.cpp.

◆ PrintSelf()

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

Definition at line 366 of file SBGATObsRadar.cpp.

◆ PrintTrailer()

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

Definition at line 360 of file SBGATObsRadar.cpp.

◆ reverse_ray_trace()

void SBGATObsRadar::reverse_ray_trace ( SBGATRadarObsSequence measurements_sequence,
const std::vector< std::vector< int > > &  facets_in_view,
const arma::vec &  radar_dir,
const int  N,
const bool  penalize_indicence,
const std::vector< arma::mat > &  BN_dcms_vec,
const std::vector< arma::vec > &  positions_vec,
const std::vector< arma::vec > &  velocities_vec,
const std::vector< arma::vec > &  omega_vec 
)
protected

Ray traces the facets in view to the radar and measurements sequence with range/range-rate data if sampled point was in view of the radar

Parameters
measurements_sequenceSequence of measurements to add new image data to
facets_in_viewreference to a vector of vector holding indices of (maybe) illuminated facets for all considered bodies
radar_dirradar 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
velocities_vecvector holding the velocities vector of the CM of each body
omega_vecvector holding the angular velocities vector of each body

Definition at line 117 of file SBGATObsRadar.cpp.

◆ SaveImages()

void SBGATObsRadar::SaveImages ( std::string  savepath)

Save the binned radar images to PNGs in the prescribed folder. The images will be normalized by the largest value in the observation sequence

Parameters
savepathpath to folder where images will be saved (ex: "output/")

Definition at line 323 of file SBGATObsRadar.cpp.

◆ vtkTypeMacro()

SBGATObsRadar::vtkTypeMacro ( SBGATObsRadar  ,
vtkPolyDataAlgorithm   
)

Member Data Documentation

◆ images

std::vector<vtkSmartPointer<vtkImageData> > SBGATObsRadar::images
protected

Definition at line 159 of file SBGATObsRadar.hpp.

◆ max_value

double SBGATObsRadar::max_value
protected

Definition at line 161 of file SBGATObsRadar.hpp.


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