MyMarlinTPC  170316
Public Member Functions | Private Attributes | List of all members
marlintpc::rb_HyperCube Class Reference

Row based hypercubes. More...

#include <RowBasedFastHoughTransformationProcessor.h>

Public Member Functions

 rb_HyperCube (unsigned int *, hyperPlaneListType &, intListType &, unsigned int=0)
 Construct row based hypercube. More...
 
void print () const
 Print hypercube. More...
 
hitListType divide (int, int, int, int, float, float)
 Divide hypercube (recursively). More...
 

Private Attributes

unsigned int _setup [2]
 setup (number of parameters for pad plane (3), drift time measurement (2)) More...
 
std::vector< rb_HyperPlane * > _planes
 list of hyperplanes More...
 
intListType _distances
 list distances More...
 
int _level
 (splitting) level More...
 
int _numChild
 number of child cubes More...
 
unsigned int _dimension
 (measurement) dimension (1 or 2) More...
 

Detailed Description

Row based hypercubes.

Has entries with distances and hyperplanes. Subdivided recursively until first track candidate found or maximal (splitting) level reached.

Definition at line 183 of file RowBasedFastHoughTransformationProcessor.h.

Constructor & Destructor Documentation

◆ rb_HyperCube()

marlintpc::rb_HyperCube::rb_HyperCube ( unsigned int *  setup,
hyperPlaneListType planes,
intListType distances,
unsigned int  level = 0 
)

Construct row based hypercube.

Parameters
[in]setupnumber of parameters per measurement (pad plane, drift time)
[in]planeslist of hyperplanes (pointer)
[in]distanceslist of distances
[in]levellevel (number of cube splittings)

Definition at line 492 of file RowBasedFastHoughTransformationProcessor.cc.

References _dimension, _numChild, _setup, and marlintpc::nCube.

Member Function Documentation

◆ divide()

hitListType marlintpc::rb_HyperCube::divide ( int  maxCube,
int  minRows,
int  minLevel,
int  maxLevel,
float  effCut,
float  purCut 
)

Divide hypercube (recursively).

Parameters
[in]maxCubemaximal number of cubes to check.
[in]minRowsminimum number of rows in child.
[in]minLevelminimal level (number of splittings)
[in]maxLevelmaximal level (number of splittings)
[in]effCutefficiency cut (minimal hit density)
[in]purCutpurity cut (maximal hit density)
Returns
list of hits for track candidate (from final cube)
  • For all child cubes count rows from intersecting hyperplanes.
  • Sort child cubes by (decreasing) number of rows and (increasing) distance spread (<d**2>).
  • Check (orderd) child cubes
  1. Accept cube with proper hit density (and minimum level) as track candidate.
  2. Subdivide cube if not reached maximum level.

Definition at line 518 of file RowBasedFastHoughTransformationProcessor.cc.

References _dimension, _distances, _level, _numChild, _planes, _setup, divide(), marlintpc::rb_HyperPlane::getCut(), marlintpc::rb_HyperPlane::getRow(), marlintpc::rb_HyperPlane::getStep(), and marlintpc::nCube.

Referenced by divide(), and marlintpc::RowBasedFastHoughTransformationProcessor::processEvent().

◆ print()

void marlintpc::rb_HyperCube::print ( ) const

Print hypercube.

Definition at line 503 of file RowBasedFastHoughTransformationProcessor.cc.

References _level, _numChild, and marlintpc::nCube.

Member Data Documentation

◆ _dimension

unsigned int marlintpc::rb_HyperCube::_dimension
private

(measurement) dimension (1 or 2)

Definition at line 195 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide(), and rb_HyperCube().

◆ _distances

intListType marlintpc::rb_HyperCube::_distances
private

list distances

Definition at line 192 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide().

◆ _level

int marlintpc::rb_HyperCube::_level
private

(splitting) level

Definition at line 193 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide(), and print().

◆ _numChild

int marlintpc::rb_HyperCube::_numChild
private

number of child cubes

Definition at line 194 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide(), print(), and rb_HyperCube().

◆ _planes

std::vector<rb_HyperPlane*> marlintpc::rb_HyperCube::_planes
private

list of hyperplanes

Definition at line 191 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide().

◆ _setup

unsigned int marlintpc::rb_HyperCube::_setup[2]
private

setup (number of parameters for pad plane (3), drift time measurement (2))

Definition at line 190 of file RowBasedFastHoughTransformationProcessor.h.

Referenced by divide(), and rb_HyperCube().


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