MyMarlinTPC
170316
|
Track finder based on Fast Hough Transformation. More...
#include <RowBasedFastHoughTransformationProcessor.h>
Public Member Functions | |
virtual Processor * | newProcessor () |
RowBasedFastHoughTransformationProcessor () | |
Construct processor. More... | |
virtual void | init () |
Initialize processor. More... | |
virtual void | processRunHeader (EVENT::LCRunHeader *run) |
virtual void | processEvent (EVENT::LCEvent *evt) |
Process event. More... | |
virtual void | check (EVENT::LCEvent *evt) |
virtual void | end () |
Protected Attributes | |
std::string | _inputColName |
Name of the input collection. More... | |
std::string | _outputColName |
Name of the output collection. More... | |
double | _bfieldScaleFactor |
scale factor for magnetic field (default: 1.0) More... | |
bool | _useXY |
Use XY (anode plane) measurement. More... | |
bool | _useXZ |
Use XZ (drift time) measurement. More... | |
double | _centerXY |
Center of XY (anode plane) measurement. More... | |
double | _rangeXY |
Range of XY (anode plane) measurement. More... | |
double | _centerXZ |
Center of XZ (drift time) measurement. More... | |
double | _rangeXZ |
Range of XZ (drift time) measurement. More... | |
int | _minRow |
Minimal number of (correlated) rows. More... | |
float | _fracRow |
Relative minimum cube content. More... | |
int | _maxRowDiff |
Row correation parameter. More... | |
int | _maxCube |
Maximum number of (hyper) cubes to check. More... | |
int | _minLevel |
Minimum number of (hyper) cubes splittings. More... | |
int | _maxLevel |
Maximum number of (hyper) cubes splittings. More... | |
float | _effCut |
Minimum hit density (hits/tracklength) More... | |
float | _purCut |
Maximum hit density (hits/tracklength) More... | |
int | _maxHitsPerRow |
Maximum number hits per row. More... | |
double | _unusedCut |
Chi2 cut for matching unused hits in XY and Z. More... | |
int | _maxGap |
Cut for (row) distance to segment for matching unused hits in XY and Z. More... | |
bool | _encodedModuleID |
Module ID is encoded in CellID0. More... | |
bool | _refAtPCA |
Use Pca as reference point (else 1. hit) More... | |
bool | _timePixFlag |
Is timepix data. More... | |
Private Attributes | |
int | _scaleBits |
number of scale bits (2^(_scaleBits) = 1.0) More... | |
double | _Bzc |
magnetic field strength (Bz*c) More... | |
double | _Xcenter |
TPC center, X coordinate. More... | |
double | _Ycenter |
TPC center, Y coordinate. More... | |
Track finder based on Fast Hough Transformation.
Fast Hough transformation to construct tracks from space points.
The average track direction is in (local) X direction. The measured coordinates have to be shifted and rescaled to: X->u in [-1, +1], Y,Z -> v,w in [-0.5, 0.5]. The resolution in v and w should be comparable. The track projections v(u), w(u) are parametrized by series of Legendre polynomials L_i, in the (XY) bending plane up to order 2 and in the (XZ) plane with order 1.:
The track parameters a_i, b_i are in [-0.5, +0.5] (the initial hypercube). Each space point defines 2 hyperplanes (v(u)-v_meas=0, w(u)-w_meas=0) in the 5 (or 4) dimensional parameter space (a_i, b_i) with normal vectors (L_i(u), 0.) and (0., L_i(u)). The unit normal vector is used to calculate the distance of a hyperplane to the center of a hypercube. They intersect if the absolute value of the distance times the largest component of the unit normal vector is smaller than 0.5. The root cube is recursivly split for each dimension in two child cubes. Only childs containing a minimum number of rows are considered further until the hit density in a cube is compatible with a single track. The search stops with the first found track candidate.
Extended for timepix data (allow multiple hits per row).
"InputHits":string | The name of the input collection of TPC hits (default: "TPCHits") |
"OutputTracks":string | The name of the output collection with the found tracks (default: "FHTTracks") |
"CenterXYMeasurement":double | Center of XY (anode plane) measurement (default 60.) |
"RangeXYMeasurement":double | Range of XY (anode plane) measurement (default 200.) |
"CenterXZMeasurement":double | Center of XY (drift time) measurement (default 300.) |
"RangeXZMeasurement":double | Range of XY (drift time) measurement (default 600.) |
"BFieldScaleFactor":double | Optional parameter, scales magnetic field (map), use 1.0 (default) for field ON or 0.0 for field OFF |
"UseXYMeasurement":bool | Optional parameter, flag for using XY (anode plane) measurement (default true) |
"UseXZMeasurement":bool | Optional parameter, flag for using XZ (drift time) measurement (default true) |
"MinimumRows":int | Optional parameter, minimal number of (correlated) rows (default 8) |
"MaximumRowDifference":int | Optional parameter, row correation parameter (max. distance of correlated rows) (default 1) |
"FractionOfRows":float | Optional parameter, relative minimum cube content (rows >= fraction * (correlated) rows) (default 0.8) |
"MaximumCubes":int | Optional parameter, maximum number of (hyper) cubes to check (default 250) |
"MinimumLevel":int | Optional parameter, minimum number of (hyper) cubes splittings (default 5) |
"MaximumLevel":int | Optional parameter, maximum number of (hyper) cubes splittings (default 8) |
"EfficienyCut":float | Optional parameter, minimum hit density (hits/tracklength) (default 0.8) |
"PurityCut":float | Optional parameter, maximum hit density (hits/tracklength) (default 1.1) |
"MaxHitsPerRow":int | Optional parameter, mximum number of hits per row (default 1) |
"UnusedHitMatchingChi2Cut":double | Optional parameter, Chi2 cut for matching unused hits in XY and Z (default 20.) |
"UnusedHitMaxGapCut:int | Optional parameter, maximal (row) gap for matching unused hits in XY and Z (default 4) |
"EncodedModuleID":bool | Optional parameter, flag for encoding of module ID in CellID0 (default true) |
"ReferencePointAtPca":bool | Optional parameter, use PCA as reference point, else 1. hit (default false) |
"TimePixFlag":bool | Optional parameter, is timepix data(default false) |
A detailed description is available as lcnote LC-TOOL-2014-006.
Definition at line 92 of file RowBasedFastHoughTransformationProcessor.h.
marlintpc::RowBasedFastHoughTransformationProcessor::RowBasedFastHoughTransformationProcessor | ( | ) |
Construct processor.
Definition at line 42 of file RowBasedFastHoughTransformationProcessor.cc.
References _bfieldScaleFactor, _centerXY, _centerXZ, _effCut, _encodedModuleID, _fracRow, _inputColName, _maxCube, _maxGap, _maxHitsPerRow, _maxLevel, _maxRowDiff, _minLevel, _minRow, _outputColName, _purCut, _rangeXY, _rangeXZ, _refAtPCA, _timePixFlag, _unusedCut, _useXY, and _useXZ.
Referenced by newProcessor().
|
virtual |
Definition at line 408 of file RowBasedFastHoughTransformationProcessor.cc.
Referenced by newProcessor().
|
virtual |
Definition at line 412 of file RowBasedFastHoughTransformationProcessor.cc.
Referenced by newProcessor().
|
virtual |
Initialize processor.
Definition at line 97 of file RowBasedFastHoughTransformationProcessor.cc.
Referenced by newProcessor().
|
inlinevirtual |
Definition at line 95 of file RowBasedFastHoughTransformationProcessor.h.
References check(), end(), init(), processEvent(), processRunHeader(), and RowBasedFastHoughTransformationProcessor().
|
virtual |
Process event.
Definition at line 107 of file RowBasedFastHoughTransformationProcessor.cc.
References _bfieldScaleFactor, _Bzc, _centerXY, _centerXZ, _effCut, _encodedModuleID, _fracRow, _inputColName, _maxCube, _maxGap, _maxHitsPerRow, _maxLevel, _maxRowDiff, _minLevel, _minRow, _outputColName, _purCut, _rangeXY, _rangeXZ, _refAtPCA, _scaleBits, _timePixFlag, _unusedCut, _useXY, _useXZ, _Xcenter, _Ycenter, marlintpc::rb_Segment::addHit(), marlintpc::rb_HyperCube::divide(), marlintpc::rb_Segment::fillRowMap(), marlintpc::rb_Segment::getFirstRow(), marlintpc::rb_Segment::getHitList(), marlintpc::rb_Segment::getLastRow(), marlintpc::rb_Segment::getNdf(), and marlintpc::rb_Segment::match().
Referenced by newProcessor().
|
virtual |
Definition at line 102 of file RowBasedFastHoughTransformationProcessor.cc.
Referenced by newProcessor().
|
protected |
scale factor for magnetic field (default: 1.0)
Definition at line 115 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
private |
magnetic field strength (Bz*c)
Definition at line 139 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent().
|
protected |
Center of XY (anode plane) measurement.
Definition at line 118 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Center of XZ (drift time) measurement.
Definition at line 120 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Minimum hit density (hits/tracklength)
Definition at line 128 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Module ID is encoded in CellID0.
Definition at line 133 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Relative minimum cube content.
Definition at line 123 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Name of the input collection.
Definition at line 113 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Maximum number of (hyper) cubes to check.
Definition at line 125 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Cut for (row) distance to segment for matching unused hits in XY and Z.
Definition at line 132 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Maximum number hits per row.
Definition at line 130 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Maximum number of (hyper) cubes splittings.
Definition at line 127 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Row correation parameter.
Definition at line 124 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Minimum number of (hyper) cubes splittings.
Definition at line 126 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Minimal number of (correlated) rows.
Definition at line 122 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Name of the output collection.
Definition at line 114 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Maximum hit density (hits/tracklength)
Definition at line 129 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Range of XY (anode plane) measurement.
Definition at line 119 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Range of XZ (drift time) measurement.
Definition at line 121 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Use Pca as reference point (else 1. hit)
Definition at line 134 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
private |
number of scale bits (2^(_scaleBits) = 1.0)
Definition at line 138 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent().
|
protected |
Is timepix data.
Definition at line 135 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Chi2 cut for matching unused hits in XY and Z.
Definition at line 131 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Use XY (anode plane) measurement.
Definition at line 116 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
protected |
Use XZ (drift time) measurement.
Definition at line 117 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent(), and RowBasedFastHoughTransformationProcessor().
|
private |
TPC center, X coordinate.
Definition at line 140 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent().
|
private |
TPC center, Y coordinate.
Definition at line 141 of file RowBasedFastHoughTransformationProcessor.h.
Referenced by processEvent().