JG ROOT Tools libraries  version 5.0 - august 2014
RHoughTransform3D Class Reference

#include <RHoughTransform3D.hh>

Public Member Functions

Constructors, affectation, destructor
 RHoughTransform3D (TH3 *hptr=NULL, u_short div=2, u_int nxy=0, double xymax=0., double thr=0., const string &id="")
 
virtual ~RHoughTransform3D ()
 
Transform analysis options
u_int GetOptions () const
 
bool GetOptionThreshold () const
 
bool GetOptionAmplWeight () const
 
bool GetOptionThresholdAmplWeight () const
 
bool GetOptionSmooth () const
 
virtual void SetOptionAmplWeight (bool thr=false)
 
virtual void SetOptionUnitWeight ()
 
virtual void SetOptionSmooth ()
 
virtual void SetOptionNoSmooth ()
 
virtual void SetOptionDegree ()
 
virtual void SetOptionRadian ()
 
virtual bool IsOptionDegree () const
 
virtual bool IsOptionRadian () const
 
Transform analysis parameters
TH3 * GetSourceHisto () const
 
TH3 * GetWeightHisto (const string &hname="") const
 
virtual double GetWeight (double pval) const
 
const TVector3 & GetSourceCenter () const
 
TVector3 & GetSourceCenter ()
 
virtual void SetSourceCenter (const TVector3 &ctr)
 
virtual int SetSourceHisto (TH3 *hptr, double thr=0.)
 
virtual void SetThreshold (double thr)
 
virtual void ResetSourceData ()
 
virtual void InitSourceData ()
 
virtual void UpdateSourceData ()
 
virtual void ResetTesselation ()
 
virtual const RSphereTesselationSetTesselation (u_short div)
 
int GetNpoints () const
 
int GetNsegments () const
 
const RSphereTesselationGetSphereTesselation () const
 
const RTesselationPointGetPoint (int ip) const
 
const RTesselationSegmentGetSegment (int is) const
 
virtual void SetTransformId (const string &id)
 
virtual void ResetTransform ()
 
virtual void ResetTransformMain ()
 
virtual void ResetTransformLocal ()
 
virtual TH2 * SetTransformMain (u_int nxy, double xymax)
 
virtual TH2 * SetTransformLocal (u_int nxy, double dxy, u_int na, double da)
 
Transform analysis functions
virtual int SetAnalysisMain ()
 
virtual int SetAnalysisLocal ()
 
virtual int ComputeTransformXY (u_int ip)
 
virtual int ComputeTransformXY (double theta, double phi, bool deg=false)
 
virtual int ComputeTransformXY (TVector3 ux, TVector3 uy)
 
TH2 * GetLocalHistoAngle ()
 
TH2 * GetLocalHistoUxUy ()
 
virtual void ResetAnalysis ()
 
virtual int AnalyseTransformMain ()
 
virtual int AnalyseTransformLocal ()
 
double GetMaxScore () const
 
double GetMaxScore (double &theta, double &phi, double &ux, double &uy) const
 
virtual void SuppressSourceLine (double theta, double phi, double ux, double uy, double width, u_short mode=modeCut, u_short sdiv=1)
 
virtual void SuppressSourceLine (TVector3 Pt, TVector3 dir, double width, u_short mode=modeCut, u_short sdiv=1)
 
virtual u_int SearchLineSuppress (double hthr, u_int nlmax=1, double *thtab=NULL, double *rtab=NULL, double *htab=NULL, u_short sdiv=1)
 
ROOT related functions
virtual void SetLineColor (Color_t c)
 
virtual void SetLineWidth (Width_t w)
 
virtual void SetLineStyle (Style_t s)
 
virtual void SetLineAtt (TAttLine &att)
 
const TAttLine & GetLineAtt () const
 
TAttLine & GetLineAtt ()
 
TPolyMarker3D * CreatePoints3D () const
 
TPolyLine3D * CreateSegment3D (int is) const
 
TList * CreateLines3D (const TAttLine latt=TAttLine(kRed+1, 2, 1)) const
 
TList * CreateScorePoints3D (Color_t col=kAzure+1, Style_t sty=20, double siz=1.0) const
 
TList * CreateLocalPoints3D (Color_t col=kRed+1, Style_t sty=20, double siz=1.0) const
 
TPolyLine3D * CreateLine3D () const
 
TPolyLine3D * CreateLine3D (double theta, double phi, double ux, double uy) const
 
 ClassDef (RHoughTransform3D, 0)
 

Protected Attributes

RSphereTesselation sphere
 Discretized representation of directions on a sphere.
 
int * sph_pts
 Array of selected points for unique representation of 3D directions.
 
int * sph_seg
 Array of segments connecting the selected points.
 
TVector3 * sph_ux
 Unit X vector for discretization directions.
 
TVector3 * sph_uy
 Unit X vector for discretization directions.
 
double * sph_score
 Score for discretization directions.
 
u_int npts
 Number of selected points.
 
u_int nseg
 Number of selected segments.
 
TH3 * src_histo
 Pointer to the histogram to analyse.
 
TVector3 src_center
 Image histogram center.
 
double threshold
 Threshold value.
 
u_int src_data_dim
 Reduced source data dimension.
 
u_int src_data_num
 Number of source data points.
 
TVector3 * src_data_pos
 Array of source data point positions.
 
double * src_data_weight
 Array of source data point weight.
 
TH2 * h_anal_xy
 Pointer to the Ux-Uy weight histogram for the analysed direction.
 
string hxy_id
 Identifier for created histogram name.
 
double hxy_max
 Maximum value for current Ux-Uy weight histogram.
 
double hxy_max_x
 X value for analysis maximum.
 
double hxy_max_y
 Y value for analysis maximum.
 
double hxy_max_theta
 Theta value for analysis maximum.
 
double hxy_max_phi
 Phi value for analysis maximum.
 
bool analysis_done
 Whether the analysis has been performed.
 
double hxy_last_score
 Maximum score for current direction.
 
double hxy_last_x
 X value for current direction maximum.
 
double hxy_last_y
 Y value for current direction maximum.
 
double hxy_last_theta
 Theta value for current direction maximum.
 
double hxy_last_phi
 Phi value for current direction maximum.
 
TH2 * h_main_xy
 Pointer to the Ux-Uy score histogram for the main analysis (all directions)
 
TH2 * h_local_xy
 Pointer to the Ux-Uy score histogram for the improved local analysis.
 
TH2 * h_local_ang
 Pointer to the Theta-Phi score histogram for the improved local analysis.
 
TVector3 * pt_local
 Array of vectors from local analysis.
 
bool local_analysis
 Whether the ongoing analysis is local.
 
TAttLine line_att
 Lines draw attributes (ROOT)
 
u_int options
 Transform options.
 

Private Member Functions

 GObject (RHoughTransform3D)
 

Static Private Attributes

static double epsilon = 0.0001
 
static const u_int optThreshold = 0x00000002
 
static const u_int optAmplitude = 0x00000004
 
static const u_int optSmooth = 0x00000008
 
static const u_int optdegree = 0x00000010
 
static const u_short modeCut = 0x0000
 
static const u_short modeGaus = 0x0001
 
static const double rad2deg = 180./TMath::Pi()
 
static const double deg2rad = TMath::Pi()/180.
 

Detailed Description

This class is a 3D adaptation of the Hough transform analysis, in order to search for lines in a 3D histogram.

The discretization of 3D directions is based on the tesselation of a sphere starting from an icosahedron (see class RSphereTesselation). The analysis is limited to the directions on the unit sphere defining unique directions, then only half of the tesselation points are considered.

Line detection main procedure

The line searching principle is the following:

  • each point of the sphere defines a (theta,phi) direction, with direction vector Uz;
  • for each direction, a perpendicular plane (Ux,Uy) is defined, and a corresponding 2D histogram is created (the histogram is created only once, and reused for each direction);
    • for each bin of the source 3D histogram, the (ux,uy,uz) coordinates are computed to score the (Ux,Uy) histogram;
    • the maximum value (score) for the (Ux,Uy) histogram defines the score for the direction (theta,phi);
  • the direction with the maximum score is considered for a line detection:
    • if the score is higher than the detection threshold, the detection is valid

Main procedure parameters:

  • directions discretization level (sphere tesselation)
  • 3D source (signal) histogram center, as origin point for the projections
  • 2D score histogram: size and binning along Ux and Uy directions (the histogram dimensions must be large enough to cover the projections of the source histogram)

Line detection improvement

The main detection procedure precision is limited to the sphere tessalation. Once a line is detected, a local analysis allows to refine the line parameters around the result of the main procedure.

Multiple lines search: detected line suppression

Configuration functions

Configuration at construction time:

TH3 * hsrc = ... // define the source histogram
RHoughTransform3D h3d ( hsrc, u_short div = 2,
u_int nxy = 100, double xymax = 100.,
double thr = 0. );

Configuration functions:

TH3 * hsrc = ... // define the source histogram
RHoughTransform3D h3d; // empty analysis
double thr = 0.2; // threshold for source histogram
h3d.SetSourceHisto ( hsrc, thr );
u_short div = 2; // number of division for sphere tesselation (3D directions discretization)
h3d.SetTesselation ( 2 );
u_int nxy = 100; // number of bins for (Ux,Uy) score histogram
double xymax = 100.; // (Ux,Uy) score histogram axis limits: [-xymax;+xymax]
h3d.SetTransformMain ( nxy, xymax, "UxUy_Score" );
Line search functions

The main function to perform the line analysis search is RHoughTransform3D::AnalyseTransformMain().

This function computes the score histogram (Ux,Uy) for each direction of the 3D directions discretization. The direction with maximum score is memorized as the result of the analysis.

The score histogram for 1 direction is filled with the RHoughTransform3D::ComputeTransformXY() function.

Constructor & Destructor Documentation

RHoughTransform3D::RHoughTransform3D ( TH3 *  hptr = NULL,
u_short  div = 2,
u_int  nxy = 0,
double  xymax = 0.,
double  thr = 0.,
const string &  id = "" 
)

Constructor.

Parameters
hptrsource (image) histogram
divsphere tesselation division
nxynumber of bins of Ux-Uy weight histogram
xymaxmaximum distance to consider
thrthreshold for points selection
idtransform identifier (name)
RHoughTransform3D::~RHoughTransform3D ( )
virtual

Destructor.

References ResetSourceData(), and ResetTransform().

Member Function Documentation

RHoughTransform3D::GObject ( RHoughTransform3D  )
private

Macro from GCpp library that defines the following functions:

  • ClassName(): return the real class name of the object
  • StaticClassName(): return the used class name of the object (that may be a base class of the real object).
  • IsInstanceOf<T>(): return true if the current object is an instance of the template class name argument
  • Clone(): return an allocated copy of the object.
u_int RHoughTransform3D::GetOptions ( ) const
inline

Return the option flags for the transform.

References options.

bool RHoughTransform3D::GetOptionThreshold ( ) const
inline

Return true if the option for threshold validation is set.

References options, and optThreshold.

Referenced by GetWeight().

bool RHoughTransform3D::GetOptionAmplWeight ( ) const
inline

Return true if the option for weight from amplitude is set.

References optAmplitude, and options.

Referenced by GetWeight().

bool RHoughTransform3D::GetOptionThresholdAmplWeight ( ) const
inline

Return true if the options for weight from amplitude and threshold validation are both set.

References optAmplitude, options, and optThreshold.

bool RHoughTransform3D::GetOptionSmooth ( ) const
inline

Return true if the option for smoothing the Hough transform histogram is set.

References options, and optSmooth.

Referenced by ComputeTransformXY().

void RHoughTransform3D::SetOptionAmplWeight ( bool  thr = false)
virtual

Set the weight of image histogram bins to be the amplitude. All bins above threshold are taken into account if the threshold is defined.

Parameters
thrindicates whether the threshold is considered

References optAmplitude, options, and optThreshold.

void RHoughTransform3D::SetOptionUnitWeight ( )
virtual

Set the weight of image histogram bins to be 1. In this case, the threshold is automatically taken into account for bins to be taken into account.

References optAmplitude, options, and optThreshold.

void RHoughTransform3D::SetOptionSmooth ( )
inlinevirtual

Set the option for smoothing the Hough transform histogram.

References options, and optSmooth.

void RHoughTransform3D::SetOptionNoSmooth ( )
inlinevirtual

Unset the option for smoothing the Hough transform histogram.

References options, and optSmooth.

void RHoughTransform3D::SetOptionDegree ( )
inlinevirtual

Set the option for angles given in degrees.

References optdegree, and options.

void RHoughTransform3D::SetOptionRadian ( )
inlinevirtual

Set the option for angles given in radians.

References optdegree, and options.

bool RHoughTransform3D::IsOptionDegree ( ) const
inlinevirtual

Return true is the angles are given in degrees.

References optdegree, and options.

Referenced by AnalyseTransformLocal(), and IsOptionRadian().

bool RHoughTransform3D::IsOptionRadian ( ) const
inlinevirtual

Return true is the angles are given in radians.

References IsOptionDegree().

TH3 * RHoughTransform3D::GetSourceHisto ( ) const
inline

Return the pointer to the source image histogram.

References src_histo.

TH3 * RHoughTransform3D::GetWeightHisto ( const string &  hname = "") const

Create an histogram from the source image histo with weight of each bin according to current options.

Parameters
hnamecreated histogram name

References GetWeight(), and src_histo.

double RHoughTransform3D::GetWeight ( double  pval) const
virtual

Return the weight of a point according to the current options of the transform class.

Parameters
pvalamplitude of the point signal

References GetOptionAmplWeight(), GetOptionThreshold(), and threshold.

Referenced by GetWeightHisto(), and UpdateSourceData().

const TVector3 & RHoughTransform3D::GetSourceCenter ( ) const
inline

Return the center point of the source histogram.

References src_center.

TVector3 & RHoughTransform3D::GetSourceCenter ( )
inline

Return the center point of the source histogram.

References src_center.

void RHoughTransform3D::SetSourceCenter ( const TVector3 &  ctr)
virtual

Define the histogram center for transformations.

Parameters
ctrcenter point

References src_center, src_histo, and UpdateSourceData().

int RHoughTransform3D::SetSourceHisto ( TH3 *  hptr,
double  thr = 0. 
)
virtual

Set the source (image) histogram. The function returns 0 if no error occured.

Parameters
hptrsource (image) histogram
thrthreshold for points selection

References SetThreshold(), src_histo, and UpdateSourceData().

void RHoughTransform3D::SetThreshold ( double  thr)
inlinevirtual

Set the threshold value to take a point in the analysis.

Parameters
thrthreshold value

References threshold.

Referenced by SetSourceHisto().

void RHoughTransform3D::ResetSourceData ( )
virtual

Function that deletes the source data arrays built from the source histogram.

References src_data_dim, src_data_num, src_data_pos, and src_data_weight.

Referenced by InitSourceData(), and ~RHoughTransform3D().

void RHoughTransform3D::InitSourceData ( )
virtual

Function that creates the source data arrays built from the source histogram, if needed.

References ResetSourceData(), src_data_dim, src_data_num, src_data_pos, src_data_weight, and src_histo.

Referenced by UpdateSourceData().

void RHoughTransform3D::UpdateSourceData ( )
virtual

Function that fills the source data arrays built from the source histogram, if needed.

References epsilon, GetWeight(), InitSourceData(), src_center, src_data_dim, src_data_num, src_data_pos, src_data_weight, and src_histo.

Referenced by SetSourceCenter(), SetSourceHisto(), and SuppressSourceLine().

void RHoughTransform3D::ResetTesselation ( )
virtual

Reset the tesselation data.

References npts, nseg, sph_pts, sph_score, sph_seg, sph_ux, and sph_uy.

Referenced by SetTesselation().

int RHoughTransform3D::GetNpoints ( ) const
inline

Return the number of effective points for unique discretization of 3D directions.

References npts.

int RHoughTransform3D::GetNsegments ( ) const
inline

Return the number of effective segments from the unique discretization of 3D directions.

References nseg.

const RSphereTesselation & RHoughTransform3D::GetSphereTesselation ( ) const
inline

Return the objet defining the sphere tesselation for the 3D space directions discretization.

References sphere.

const RTesselationPoint & RHoughTransform3D::GetPoint ( int  ip) const
inline

Return a point (by constant reference).

Parameters
ippoint index

References RSolidTesselation::GetPoint(), sph_pts, and sphere.

const RTesselationSegment & RHoughTransform3D::GetSegment ( int  is) const
inline

Return a segment (by constant reference).

Parameters
issegment index

References RSolidTesselation::GetSegment(), sph_seg, and sphere.

void RHoughTransform3D::SetTransformId ( const string &  id)
inlinevirtual

Define the base name for created histogram. The analysis transforms (main and local) must be resetted to recreate the histograms. The effective histogram names are:

  • "\<id\>_MainXY" for the (Ux,Uy) score histogram for main analysis
  • "\<id\>_ImprXY" for the (Ux,Uy) score histogram for local analysis
  • "\<id\>_ImprAng" for the (Theta,Phi) score histogram for local analysis
    Parameters
    idhistograms base name identifier

References hxy_id.

void RHoughTransform3D::ResetTransform ( )
virtual

Reset the transformation data for the main and improved analysis.

References ResetTransformLocal(), and ResetTransformMain().

Referenced by ~RHoughTransform3D().

void RHoughTransform3D::ResetTransformMain ( )
virtual

Reset the transformation data for the main analysis.

References h_main_xy.

Referenced by ResetTransform(), and SetTransformMain().

void RHoughTransform3D::ResetTransformLocal ( )
virtual

Reset the transformation data for the improved analysis.

References h_local_ang, h_local_xy, and pt_local.

Referenced by ResetTransform(), and SetTransformLocal().

TH2 * RHoughTransform3D::SetTransformMain ( u_int  nxy,
double  xymax 
)
virtual

Function defining the score histogram for the global transform analysis. The source histogram pointer must be defined.

If the maximum distance is 0, it is taken from the image histogram diagonal.

If the transform is created successfully, the Ux-Uy weight histogram pointer is returned, otherwise a NULL pointer.

Parameters
nxynumber of bins of Ux-Uy weight histogram
xymaxmaximum distance to consider

References h_main_xy, hxy_id, R_H3D_MAIN_XY, ResetTransformMain(), and src_histo.

TH2 * RHoughTransform3D::SetTransformLocal ( u_int  nxy,
double  dxy,
u_int  na,
double  da 
)
virtual

Function defining the score histogram for the local transform analysis (see local transform analysis in RHoughTransform3D class documentation).

The global transform data must have been previously defined.

The function creates the (Ux,Uy) score histogram, and also a (theta,phi) histogram that stores the maximum results for (Ux,Uy) for each angle. The axis of angles score histogram are in degrees/radians according to the corresponding option.

If the transform is created successfully, the (Ux,Uy) score histogram pointer is returned, otherwise a NULL pointer.

Parameters
nxynumber of bins of (Ux,Uy) weight histogram
dxy(Ux,Uy) score histogram bin size
nanumber of bins of (theta,phi) score histogram
da(theta,phi) score histogram bin size

References RSphereTesselation::GetAngleMin(), h_local_ang, h_local_xy, h_main_xy, hxy_id, optdegree, options, pt_local, R_H3D_LOCAL_ANG, R_H3D_LOCAL_XY, rad2deg, ResetTransformLocal(), and sphere.

int RHoughTransform3D::SetAnalysisMain ( )
virtual

Select score histogram for the main global analysis.

References h_anal_xy, h_main_xy, and local_analysis.

Referenced by AnalyseTransformMain().

int RHoughTransform3D::SetAnalysisLocal ( )
virtual

Select score histogram for the local (improved) analysis.

References h_anal_xy, h_local_ang, h_local_xy, h_main_xy, and local_analysis.

Referenced by AnalyseTransformLocal().

int RHoughTransform3D::ComputeTransformXY ( u_int  ip)
virtual

Function computing the Ux-Uy weight histogram with current parameters, for a point in the current discretization of space.

Parameters
ippoint index

References npts, sph_ux, and sph_uy.

Referenced by AnalyseTransformLocal(), AnalyseTransformMain(), and ComputeTransformXY().

int RHoughTransform3D::ComputeTransformXY ( double  theta,
double  phi,
bool  deg = false 
)
virtual

Function computing the Ux-Uy weight histogram with current parameters, for a (theta,phi) direction.

Parameters
thetaangle with Z axis (0:pi/2)
phiXY projection angle (2pi)
degwhether angle are in degree

References ComputeTransformXY(), RNormalUx(), and RNormalUy().

int RHoughTransform3D::ComputeTransformXY ( TVector3  ux,
TVector3  uy 
)
virtual

Function computing the Ux-Uy weight histogram with current parameters, for Ux-Uy plane vectors. The function returns -1 if a new maximum has been found, else 0, or a positive value in case of error.

Parameters
uxX unit vector
uyY unit vector

References GetOptionSmooth(), h_anal_xy, hxy_last_phi, hxy_last_score, hxy_last_theta, hxy_last_x, hxy_last_y, hxy_max, hxy_max_phi, hxy_max_theta, hxy_max_x, hxy_max_y, local_analysis, rad2deg, RGetString(), src_data_num, src_data_pos, src_data_weight, and src_histo.

TH2 * RHoughTransform3D::GetLocalHistoAngle ( )
inline

Return the (theta,phi) weight histogram from local analysis.

References h_local_ang.

TH2 * RHoughTransform3D::GetLocalHistoUxUy ( )
inline

Return the (Ux,Uy) weight histogram from local analysis.

References h_local_xy.

void RHoughTransform3D::ResetAnalysis ( )
virtual

Reset the maximum search analysis.

References analysis_done, hxy_max, hxy_max_phi, hxy_max_theta, hxy_max_x, hxy_max_y, and local_analysis.

Referenced by AnalyseTransformLocal(), and AnalyseTransformMain().

int RHoughTransform3D::AnalyseTransformMain ( )
virtual

For all direction in the space discretization, compute the Ux-Uy weight histograms to extract the one with maximum score. The function returns the point corresponding index.

References analysis_done, ComputeTransformXY(), hxy_last_score, hxy_max, npts, ResetAnalysis(), SetAnalysisMain(), and sph_score.

int RHoughTransform3D::AnalyseTransformLocal ( )
virtual

For all direction of the local analysis, compute the Ux-Uy weight histograms to extract the one with maximum score. The function returns 0 if no error occured.

The results of the analysis are stored in the data accessible with the RHoughTransform3D::GetMaxScore functions.

References analysis_done, ComputeTransformXY(), deg2rad, GetMaxScore(), h_local_ang, h_local_xy, hxy_last_score, hxy_max, IsOptionDegree(), local_analysis, pt_local, rad2deg, ResetAnalysis(), RNormalVectors(), and SetAnalysisLocal().

double RHoughTransform3D::GetMaxScore ( ) const

Return the maximum score found in current analysis.

References analysis_done, and hxy_max.

Referenced by AnalyseTransformLocal(), CreateLine3D(), and GetMaxScore().

double RHoughTransform3D::GetMaxScore ( double &  theta,
double &  phi,
double &  ux,
double &  uy 
) const

Return the maximum score found in current analysis. The arguments receive the corresponding coordinate values for this maximum score.

Parameters
thetatheta angle of maximum score
phiphi angle of maximum score
uxUx position in Ux-Uy weight histogram
uyUy position in Ux-Uy weight histogram

References GetMaxScore(), hxy_max_phi, hxy_max_theta, hxy_max_x, hxy_max_y, optdegree, options, and rad2deg.

void RHoughTransform3D::SuppressSourceLine ( double  theta,
double  phi,
double  ux,
double  uy,
double  width,
u_short  mode = modeCut,
u_short  sdiv = 1 
)
virtual

Create a copy of the source histogram, for which all the bins closer that r from the line defined by theta and dist arguments are set to 0.

Parameters
thetaline theta angle parameter (degree)
philine phi angle parameter (degree)
uxline position along Ux axis
uyline position along Uy axis
widthdistance to the line of points to remove (a vector for anisotropic width)
modesuppression mode:
  • 0 for sharp cut
  • 1 for gauss shape reduction cut
sdivhistogram bins subdivision

References deg2rad, RGetVector(), RNormalVectors(), src_center, and src_histo.

void RHoughTransform3D::SuppressSourceLine ( TVector3  Pt,
TVector3  dir,
double  width,
u_short  mode = modeCut,
u_short  sdiv = 1 
)
virtual

Create a copy of the source histogram, for which all the bins closer that r from the line defined by the pointPt and the direction dir are set to 0.

Parameters
Ptpoint of the line to suppress
dirdirection of the line to suppress
widthdistance to the line of points to remove (a vector for anisotropic width)
modesuppression mode:
  • 0 for sharp cut
  • 1 for gauss shape reduction cut
sdivhistogram bins subdivision

References src_histo, and UpdateSourceData().

u_int RHoughTransform3D::SearchLineSuppress ( double  hthr,
u_int  nlmax = 1,
double *  thtab = NULL,
double *  rtab = NULL,
double *  htab = NULL,
u_short  sdiv = 1 
)
virtual

Search for lines in the image histogram. When a line is found, it is suppressed from the image in order to search for the next line.

The function returns the number of lines found.

Be careful when setting no limit to the number of lines, since in case of problem, it may result in writing results out of the boundaries of thtab and rtab if they are defined.

Parameters
hthrthreshold in the Hough transform to accept a line detection
nlmaxmaximum number of lines to search (0 for no limit)
thtabarray of resulting angles (if not NULL)
rtabarray of resulting distances (if not NULL)
htabarray of resulting Hough transform maxima (if not NULL)
sdivimage histogram bins subdivision for line suppression

References h_anal_xy, and src_histo.

void RHoughTransform3D::SetLineColor ( Color_t  c)
inlinevirtual

Set the line drawing color.

Parameters
cline color

References line_att.

void RHoughTransform3D::SetLineWidth ( Width_t  w)
inlinevirtual

Set the line drawing width.

Parameters
wline width

References line_att.

void RHoughTransform3D::SetLineStyle ( Style_t  s)
inlinevirtual

Set the line drawing style.

Parameters
sline style

References line_att.

void RHoughTransform3D::SetLineAtt ( TAttLine &  att)
inlinevirtual

Set the line drawing attributes.

Parameters
attline attributes

References line_att.

const TAttLine & RHoughTransform3D::GetLineAtt ( ) const
inline

Set the line drawing attributes.

References line_att.

TAttLine & RHoughTransform3D::GetLineAtt ( )
inline

Set the line drawing attributes.

References line_att.

TPolyMarker3D * RHoughTransform3D::CreatePoints3D ( ) const
TPolyLine3D * RHoughTransform3D::CreateSegment3D ( int  is) const
TList * RHoughTransform3D::CreateLines3D ( const TAttLine  latt = TAttLine(kRed+1,2,1 )) const

Create a list of TPolyLine3D objects (ROOT) from segments.

Parameters
lattline attributes set to all segments

References CreateSegment3D(), and nseg.

TList * RHoughTransform3D::CreateScorePoints3D ( Color_t  col = kAzure+1,
Style_t  sty = 20,
double  sizmax = 1.0 
) const

Create a list of 1-point TPolyMarker3D objects (ROOT) from points from the sphere points analysis. The points size is proportional to the score.

Parameters
colpoints color
stypoints style
sizmaxpoint size for the maximum score

References analysis_done, RSolidTesselation::GetPoint(), RTesselationPoint::GetX(), RTesselationPoint::GetY(), RTesselationPoint::GetZ(), npts, sph_pts, sph_score, and sphere.

TList * RHoughTransform3D::CreateLocalPoints3D ( Color_t  col = kRed+1,
Style_t  sty = 20,
double  sizmax = 1.0 
) const

Create a list of 1-point TPolyMarker3D objects (ROOT) from points from the local analysis. The points size is proportional to the score.

Parameters
colpoints color
stypoints style
sizmaxpoint size for the maximum score

References analysis_done, h_local_ang, local_analysis, and pt_local.

TPolyLine3D * RHoughTransform3D::CreateLine3D ( ) const

Create a 3D line in the source histogram coordinates, with currently defined maximum score line parameters.

References GetMaxScore().

TPolyLine3D * RHoughTransform3D::CreateLine3D ( double  theta,
double  phi,
double  ux,
double  uy 
) const

Create a 3D line in the source histogram coordinates.

Parameters
thetatheta angle of maximum score
phiphi angle of maximum score
uxUx position in Ux-Uy weight histogram
uyUy position in Ux-Uy weight histogram

References deg2rad, optdegree, options, RLine3D(), RNormalVectors(), src_center, and src_histo.

RHoughTransform3D::ClassDef ( RHoughTransform3D  ,
 
)

For ROOT dictionary.

Member Data Documentation

double RHoughTransform3D::epsilon = 0.0001
staticprivate

Weight limit for ignoring point

Referenced by UpdateSourceData().

const u_int RHoughTransform3D::optThreshold = 0x00000002
staticprivate

Option for threshold on histogram amplitude for contributing to the analysis

Referenced by GetOptionThreshold(), GetOptionThresholdAmplWeight(), SetOptionAmplWeight(), and SetOptionUnitWeight().

const u_int RHoughTransform3D::optAmplitude = 0x00000004
staticprivate

Option for amplitude as weight of histogram point

Referenced by GetOptionAmplWeight(), GetOptionThresholdAmplWeight(), SetOptionAmplWeight(), and SetOptionUnitWeight().

const u_int RHoughTransform3D::optSmooth = 0x00000008
staticprivate

Option for smoothing the Hough transform histogram

Referenced by GetOptionSmooth(), SetOptionNoSmooth(), and SetOptionSmooth().

const u_int RHoughTransform3D::optdegree = 0x00000010
staticprivate
const u_short RHoughTransform3D::modeCut = 0x0000
staticprivate

Mode for line suppression: sharp cut

const u_short RHoughTransform3D::modeGaus = 0x0001
staticprivate

Mode for line suppression: gauss shape reduction

const double RHoughTransform3D::rad2deg = 180./TMath::Pi()
staticprivate

Conversion factor from radians to degrees

Referenced by AnalyseTransformLocal(), ComputeTransformXY(), GetMaxScore(), and SetTransformLocal().

const double RHoughTransform3D::deg2rad = TMath::Pi()/180.
staticprivate

Conversion factor from degrees to radians

Referenced by AnalyseTransformLocal(), CreateLine3D(), and SuppressSourceLine().


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