Computes range/range-rate Doppler images over the surface of provided small body.
More...
#include <SBGATObsRadar.hpp>
|
| 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 () |
|
| 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 |
|
|
| 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) |
|
| 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 |
|
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
- Copyright
- MIT License, Benjamin Bercovici and Jay McMahon
Definition at line 56 of file SBGATObsRadar.hpp.
◆ SBGATObsRadar() [1/2]
SBGATObsRadar::SBGATObsRadar |
( |
| ) |
|
|
protected |
◆ ~SBGATObsRadar()
SBGATObsRadar::~SBGATObsRadar |
( |
| ) |
|
|
overrideprotected |
◆ SBGATObsRadar() [2/2]
◆ 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_sequence | reference to MeasurementsSequence, holding collected range/range-rate measurements at each observation time |
r_bin | range bin size (m) |
rr_bin | range-rate bin size (m/s) |
Definition at line 225 of file SBGATObsRadar.cpp.
◆ ClearImages()
void SBGATObsRadar::ClearImages |
( |
| ) |
|
|
inline |
◆ 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_sequence | reference to MeasurementsSequence, holding collected range/range-rate measurements at each observation time |
time | observation timestamp |
N | minimum 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_dir | unit direction towards radar from target in inertial frame |
positions_vec | vector of positions of each target's center-of-mass expressed in the primary's body frame |
velocities_vec | vector of inertial velocities of each target's center-of-mass expressed in the primary's body frame |
mrps_vec | vector of MRPs defining the inertial-to-body DCM [BN] |
omegas_vec | vector of angular velocities of each body expressed in the inertial frame |
penalize_incidence | if 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()
Constructs with initial values of zero.
◆ operator=()
◆ PrintHeader()
void SBGATObsRadar::PrintHeader |
( |
std::ostream & |
os, |
|
|
vtkIndent |
indent |
|
) |
| |
|
override |
◆ PrintSelf()
void SBGATObsRadar::PrintSelf |
( |
std::ostream & |
os, |
|
|
vtkIndent |
indent |
|
) |
| |
|
override |
◆ PrintTrailer()
void SBGATObsRadar::PrintTrailer |
( |
std::ostream & |
os, |
|
|
vtkIndent |
indent |
|
) |
| |
|
override |
◆ 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_sequence | Sequence of measurements to add new image data to |
facets_in_view | reference to a vector of vector holding indices of (maybe) illuminated facets for all considered bodies |
radar_dir | radar direction expressed in inertial frame |
BN_dcms_vec | vector holding the DCMs orienting the body frame of each body w/r to inertial |
positions_vec | vector holding the position vector of the CM of each body w/r to the primary |
velocities_vec | vector holding the velocities vector of the CM of each body |
omega_vec | vector 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
-
savepath | path to folder where images will be saved (ex: "output/") |
Definition at line 323 of file SBGATObsRadar.cpp.
◆ vtkTypeMacro()
SBGATObsRadar::vtkTypeMacro |
( |
SBGATObsRadar |
, |
|
|
vtkPolyDataAlgorithm |
|
|
) |
| |
◆ images
std::vector<vtkSmartPointer<vtkImageData> > SBGATObsRadar::images |
|
protected |
◆ max_value
double SBGATObsRadar::max_value |
|
protected |
The documentation for this class was generated from the following files:
- /Users/bbercovici/GDrive/CUBoulder/Research/code/SBGAT/SbgatCore/include/SbgatCore/SBGATObsRadar.hpp
- /Users/bbercovici/GDrive/CUBoulder/Research/code/SBGAT/SbgatCore/source/SBGATObsRadar.cpp