QwAnalysis
QwTrackingWorker Class Reference

Controls all the routines involved in finding tracks in an event. More...

#include <QwTrackingWorker.h>

+ Collaboration diagram for QwTrackingWorker:

Public Member Functions

 QwTrackingWorker (QwOptions &options, const QwGeometry &geometry)
 Default constructor with name. More...
 
virtual ~QwTrackingWorker ()
 Destructor. More...
 
void ProcessOptions (QwOptions &options)
 Process command line and config file options. More...
 
int GetDebugLevel () const
 Get the debug level. More...
 
void SetDebugLevel (const int debug)
 Set the debug level. More...
 
const QwGeometry GetGeometry () const
 Get the geometry. More...
 
void SetGeometry (const QwGeometry &geometry)
 Set the geometry. More...
 
double GetMagneticFieldCurrent () const
 Get the magnetic field current. More...
 
void SetMagneticFieldCurrent (const double current)
 Set the magnetic field current. More...
 
void ProcessEvent (const QwSubsystemArrayTracking *detectors, QwEvent *event)
 Process the hit list and construct the event. More...
 

Static Public Member Functions

static void DefineOptions (QwOptions &options)
 Define command line and config file options. More...
 

Data Fields

int ngood
 number of good events More...
 
int nbad
 number of bad events More...
 
int R2Good
 
int R2Bad
 
int R3Good
 
int R3Bad
 
int nbridged
 number of beidged tracks More...
 

Private Member Functions

void InitTree (const QwGeometry &geometry)
 Initialize the pattern search tree. More...
 
 QwTrackingWorker ()
 

Private Attributes

QwGeometry fGeometry
 Detector geometry. More...
 
int fDebug
 Debug level. More...
 
bool fRegenerate
 Regenerate the search tree. More...
 
int fLevelsR2
 Region 2 levels. More...
 
int fLevelsR3
 Region 3 levels. More...
 
std::vector
< QwTrackingTreeRegion * > 
fSearchTrees
 Internal list of search trees created by QwTrackingWorker. More...
 
Momentum determination bridging methods
QwBridgingTrackFilterfBridgingTrackFilter
 Track filter. More...
 
QwRayTracerfRayTracer
 Ray tracing bridging method. More...
 
Parsed configuration options
bool fDisableTracking
 Disable all tracking. More...
 
bool fPrintPatternDatabase
 Print the pattern database. More...
 
bool fShowEventPattern
 Show all event patterns. More...
 
bool fShowMatchingPattern
 Show matching event patterns. More...
 
bool fDisableMomentum
 Disable momentum reconstruction. More...
 
bool fDisableMatrixLookup
 Disable matrix lookup momentum reconstruction. More...
 
bool fDisableRayTracer
 Disable ray tracer momentum reconstruction. More...
 
bool fMismatchPkg
 Indicates if the pkg for R2 and R3 is differnt at the same octant. More...
 
std::string fFilenameFieldmap
 Filename of the fieldmap in QW_FIELDMAP. More...
 
std::string fFilenameLookupTable
 Filename of the lookup table in QW_LOOKUP. More...
 
Helper modules
QwTrackingTreeSearchfTreeSearch
 Module that handles the tree search. More...
 
QwTrackingTreeCombinefTreeCombine
 Module that combines treelines and partial tracks. More...
 
QwTrackingTreeSortfTreeSort
 Module that sorts lists of treelines and partial tracks. More...
 
QwTrackingTreeMatchfTreeMatch
 Module that matches up VDC front and back treelines. More...
 

Detailed Description

Controls all the routines involved in finding tracks in an event.

The tracking worker coordinates the track finding based on the hit list provided to the ProcessHits() method.

Definition at line 46 of file QwTrackingWorker.h.

Constructor & Destructor Documentation

QwTrackingWorker::QwTrackingWorker ( QwOptions options,
const QwGeometry geometry 
)

Default constructor with name.

Definition at line 119 of file QwTrackingWorker.cc.

References QwLog::endl(), fBridgingTrackFilter, fDebug, fDisableMomentum, fDisableRayTracer, fDisableTracking, fRayTracer, fShowMatchingPattern, fTreeCombine, fTreeMatch, fTreeSearch, fTreeSort, QwOptions::GetValue(), InitTree(), nbad, nbridged, ngood, ProcessOptions(), QwDebug, R2Bad, R2Good, R3Bad, R3Good, QwTrackingTreeSort::SetDebugLevel(), QwTrackingTreeCombine::SetDebugLevel(), QwTrackingTreeSearch::SetDebugLevel(), QwTrackingTreeMatch::SetDebugLevel(), QwTrackingTreeCombine::SetGeometry(), SetGeometry(), QwTrackingTreeCombine::SetMaxRoad(), QwTrackingTreeCombine::SetMaxXRoad(), and QwTrackingTreeSearch::SetShowMatchingPatterns().

120 {
121  QwDebug << "###### Calling QwTrackingWorker::QwTrackingWorker ()" << QwLog::endl;
122 
123  // Debug level
124  fDebug = 0;
125 
126  // Process options
127  ProcessOptions(options);
128 
129  // Store geometry
130  SetGeometry(geometry);
131 
132  // If tracking is disabled, stop here
133  if (fDisableTracking) return;
134 
135  // Initialize the pattern database
136  InitTree(geometry);
137  // Initialize a bridging filter
139 
140  // Initialize a ray tracer bridging method
141  if (! fDisableMomentum && ! fDisableRayTracer) {
142  fRayTracer = new QwRayTracer(options);
143  // or set to null if disabled
144  } else fRayTracer = 0;
145 
146  /* Reset counters of number of good and bad events */
147  ngood = 0;
148  nbad = 0;
149 
150  R2Good = 0;
151  R2Bad = 0;
152  R3Good = 0;
153  R3Bad = 0;
154 
155  nbridged = 0;
156 
157  // Initialize the tracking modules
158  // Tracking functionality is provided by these four sub-blocks.
163 
164  // Debug level of the various tracking modules
169 
170  // Pass the geometry
171  fTreeCombine->SetGeometry(geometry);
172 
173  // Process the options and set the respective flags in the modules
175  fTreeCombine->SetMaxRoad(options.GetValue<float>("QwTracking.R2.maxroad"));
176  fTreeCombine->SetMaxXRoad(options.GetValue<float>("QwTracking.R2.maxxroad"));
177 
178  QwDebug << "###### Leaving QwTrackingWorker::QwTrackingWorker ()" << QwLog::endl;
179 }
void SetDebugLevel(int debug)
Set the debug level.
bool fDisableTracking
Disable all tracking.
void SetShowMatchingPatterns(bool show=true)
Set the flag to show matching patterns when they are found.
void SetMaxRoad(const double maxroad)
Set the maximum road width.
void SetMaxXRoad(const double maxxroad)
Set the maximum X road width (?)
int nbad
number of bad events
QwTrackingTreeSort * fTreeSort
Module that sorts lists of treelines and partial tracks.
void SetGeometry(const QwGeometry &geometry)
Set the geometry.
void SetGeometry(const QwGeometry &geometry)
Set the geometry.
QwRayTracer * fRayTracer
Ray tracing bridging method.
int fDebug
Debug level.
QwTrackingTreeMatch * fTreeMatch
Module that matches up VDC front and back treelines.
Track filter for the bridging methods.
T GetValue(const std::string &key)
Get a templated value.
Definition: QwOptions.h:240
#define QwDebug
Predefined log drain for debugging output.
Definition: QwLog.h:60
void ProcessOptions(QwOptions &options)
Process command line and config file options.
QwTrackingTreeSearch * fTreeSearch
Module that handles the tree search.
QwTrackingTreeCombine * fTreeCombine
Module that combines treelines and partial tracks.
void SetDebugLevel(const int debuglevel)
Set the debug level.
Module that matches track segments for pairs of wire planes.
bool fShowMatchingPattern
Show matching event patterns.
QwBridgingTrackFilter * fBridgingTrackFilter
Track filter.
void SetDebugLevel(int debug)
Set the debug level.
This module is used to identify good track segments versus ghost tracks/hits.
void SetDebugLevel(const int debuglevel)
int ngood
number of good events
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
bool fDisableRayTracer
Disable ray tracer momentum reconstruction.
void InitTree(const QwGeometry &geometry)
Initialize the pattern search tree.
Combines track segments and performs line fitting.
int nbridged
number of beidged tracks
bool fDisableMomentum
Disable momentum reconstruction.

+ Here is the call graph for this function:

QwTrackingWorker::~QwTrackingWorker ( )
virtual

Destructor.

Definition at line 183 of file QwTrackingWorker.cc.

References QwLog::endl(), fBridgingTrackFilter, fRayTracer, fSearchTrees, fTreeCombine, fTreeMatch, fTreeSearch, fTreeSort, QwObjectCounter< QwHitPattern >::GetObjectsAlive(), QwObjectCounter< QwTreeLine >::GetObjectsAlive(), QwObjectCounter< QwEvent >::GetObjectsAlive(), QwObjectCounter< QwPartialTrack >::GetObjectsAlive(), QwObjectCounter< QwTrack >::GetObjectsAlive(), QwObjectCounter< QwEventHeader >::GetObjectsAlive(), QwObjectCounter< QwHit >::GetObjectsAlive(), QwObjectCounter< QwTrack >::GetObjectsCreated(), QwObjectCounter< QwTreeLine >::GetObjectsCreated(), QwObjectCounter< QwPartialTrack >::GetObjectsCreated(), QwObjectCounter< QwHitPattern >::GetObjectsCreated(), QwObjectCounter< QwEventHeader >::GetObjectsCreated(), QwObjectCounter< QwEvent >::GetObjectsCreated(), QwObjectCounter< QwHit >::GetObjectsCreated(), and QwMessage.

184 {
185  // Delete helper objects
187  if (fRayTracer) delete fRayTracer;
188 
189  // Delete search trees
190  for (size_t i = 0; i < fSearchTrees.size(); i++)
191  delete fSearchTrees[i];
192  fSearchTrees.clear();
193 
194  // Delete tracking modules
195  if (fTreeSearch) delete fTreeSearch;
196  if (fTreeCombine) delete fTreeCombine;
197  if (fTreeSort) delete fTreeSort;
198  if (fTreeMatch) delete fTreeMatch;
199 
200  // Summary of tracking objects that are still alive
202  QwMessage << "Number of tracking objects still alive:" << QwLog::endl;
203  QwMessage << " QwEvent: " << QwEvent::GetObjectsAlive()
204  << " (of " << QwEvent::GetObjectsCreated() << ")" << QwLog::endl;
205  QwMessage << " QwEventHeader: " << QwEventHeader::GetObjectsAlive()
206  << " (of " << QwEventHeader::GetObjectsCreated() << ")" << QwLog::endl;
207  QwMessage << " QwHit: " << QwHit::GetObjectsAlive()
208  << " (of " << QwHit::GetObjectsCreated() << ")" << QwLog::endl;
209  QwMessage << " QwHitPattern: " << QwHitPattern::GetObjectsAlive()
210  << " (of " << QwHitPattern::GetObjectsCreated() << ")" << QwLog::endl;
211  QwMessage << " QwTreeLine: " << QwTreeLine::GetObjectsAlive()
212  << " (of " << QwTreeLine::GetObjectsCreated() << ")" << QwLog::endl;
213  QwMessage << " QwPartialTrack: " << QwPartialTrack::GetObjectsAlive()
214  << " (of " << QwPartialTrack::GetObjectsCreated() << ")" << QwLog::endl;
215  QwMessage << " QwTrack: " << QwTrack::GetObjectsAlive()
216  << " (of " << QwTrack::GetObjectsCreated() << ")" << QwLog::endl;
218 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
QwTrackingTreeSort * fTreeSort
Module that sorts lists of treelines and partial tracks.
std::vector< QwTrackingTreeRegion * > fSearchTrees
Internal list of search trees created by QwTrackingWorker.
QwRayTracer * fRayTracer
Ray tracing bridging method.
QwTrackingTreeMatch * fTreeMatch
Module that matches up VDC front and back treelines.
QwTrackingTreeSearch * fTreeSearch
Module that handles the tree search.
QwTrackingTreeCombine * fTreeCombine
Module that combines treelines and partial tracks.
QwBridgingTrackFilter * fBridgingTrackFilter
Track filter.
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
static size_t GetObjectsAlive()
Get number of objects still alive.
static size_t GetObjectsCreated()
Get number of objects ever created.

+ Here is the call graph for this function:

QwTrackingWorker::QwTrackingWorker ( )
private

Member Function Documentation

void QwTrackingWorker::DefineOptions ( QwOptions options)
static

Define command line and config file options.

Definition at line 222 of file QwTrackingWorker.cc.

References QwOptions::AddOptions(), and default_bool_value.

Referenced by DefineOptionsTracking().

223 {
224  // General options
225  options.AddOptions("Tracking options")("QwTracking.debug",
226  po::value<int>()->default_value(0),
227  "track reconstruction debug level");
228  options.AddOptions("Tracking options")("QwTracking.regenerate",
229  po::value<bool>()->default_bool_value(false),
230  "regenerate search trees");
231  options.AddOptions("Tracking options")("QwTracking.disable-tracking",
232  po::value<bool>()->default_bool_value(false),
233  "disable all tracking analysis");
234  options.AddOptions("Tracking options")("QwTracking.print-pattern-db",
235  po::value<bool>()->default_bool_value(false),
236  "print pattern database");
237  options.AddOptions("Tracking options")("QwTracking.showeventpattern",
238  po::value<bool>()->default_bool_value(false),
239  "show bit pattern for all events");
240  options.AddOptions("Tracking options")("QwTracking.showmatchingpattern",
241  po::value<bool>()->default_bool_value(false),
242  "show bit pattern for matching tracks");
243 
244  options.AddOptions("Tracking options")("QwTracking.package-mismatch",
245  po::value<bool>()->default_bool_value(false),
246  "if the package number is different for R2 and R3 at the same octant");
247 
248  // Region 2
249  options.AddOptions("Tracking options")("QwTracking.R2.levels",
250  po::value<int>()->default_value(8),
251  "number of search tree levels in region 2");
252  options.AddOptions("Tracking options")("QwTracking.R2.maxslope",
253  po::value<float>()->default_value(0.862),
254  "maximum allowed slope for region 2 tracks");
255  options.AddOptions("Tracking options")("QwTracking.R2.maxroad",
256  po::value<float>()->default_value(1.4),
257  "maximum allowed road width for region 2 tracks");
258  options.AddOptions("Tracking options")("QwTracking.R2.maxxroad",
259  po::value<float>()->default_value(25.0),
260  "maximum allowed X road width for region 2 tracks");
261  options.AddOptions("Tracking options")("QwTracking.R2.MaxMissedPlanes",
262  po::value<int>()->default_value(1),
263  "maximum number of missed planes");
264  options.AddOptions("Tracking options")("QwTracking.R2.DropWorstHit",
265  po::value<bool>()->default_bool_value(false),
266  "attempt partial track fit without worst hit");
267 
268  // Region 3
269  options.AddOptions("Tracking options")("QwTracking.R3.levels",
270  po::value<int>()->default_value(4),
271  "number of search tree levels in region 3");
272  options.AddOptions("Tracking options")("QwTracking.R3.MaxMissedWires",
273  po::value<int>()->default_value(4),
274  "maximum number of missed wires");
275  options.AddOptions("Tracking options")("QwTracking.R3.RotatorTilt",
276  po::value<bool>()->default_value(true),
277  "apply region 3 rotator tilt correction");
278 
279  // Momentum reconstruction
280  options.AddOptions("Tracking options")("QwTracking.disable-momentum",
281  po::value<bool>()->default_bool_value(false),
282  "disable the momentum reconstruction");
283  options.AddOptions("Tracking options")("QwTracking.disable-matrixlookup",
284  po::value<bool>()->default_bool_value(false),
285  "disable the use of the momentum lookup table");
286  options.AddOptions("Tracking options")("QwTracking.disable-raytracer",
287  po::value<bool>()->default_bool_value(false),
288  "disable the magnetic field map tracking");
289  options.AddOptions("Tracking options")("QwTracking.lookuptable",
290  po::value<std::string>()->default_value("QwTrajMatrix.root"),
291  "filename of the lookup table in QW_LOOKUP");
292 }
#define default_bool_value(b)
Definition: QwOptions.h:51
po::options_description_easy_init AddOptions(const std::string &blockname="Specialized options")
Add an option to a named block or create new block.
Definition: QwOptions.h:164

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

int QwTrackingWorker::GetDebugLevel ( ) const
inline

Get the debug level.

Definition at line 71 of file QwTrackingWorker.h.

References fDebug.

71 { return fDebug; };
int fDebug
Debug level.
const QwGeometry QwTrackingWorker::GetGeometry ( ) const
inline

Get the geometry.

Definition at line 76 of file QwTrackingWorker.h.

References fGeometry.

76 { return fGeometry; };
QwGeometry fGeometry
Detector geometry.
double QwTrackingWorker::GetMagneticFieldCurrent ( ) const
inline

Get the magnetic field current.

Definition at line 81 of file QwTrackingWorker.h.

References fRayTracer, and QwRayTracer::GetMagneticFieldCurrent().

Referenced by main().

81  {
83  else return 0.0;
84  }
QwRayTracer * fRayTracer
Ray tracing bridging method.
double GetMagneticFieldCurrent() const
Get magnetic field current.
Definition: QwRayTracer.h:61

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwTrackingWorker::InitTree ( const QwGeometry geometry)
private

Initialize the pattern search tree.

Region 2 contains 4 layers

Region 3 contains 8 layers

Create a new search tree

Set the maximum slope

Set the geometry

Set up the filename with the following format tree[numlayers]-[levels]-[u|l]-[1|2|3]-[n|u|v|x|y].tre

Each element of search tree will point to a pattern database

Definition at line 329 of file QwTrackingWorker.cc.

References QwLog::endl(), fLevelsR2, fLevelsR3, fPrintPatternDatabase, fRegenerate, fSearchTrees, QwDetectorInfo::GetActiveWidthZ(), QwDetectorInfo::GetElementSpacing(), getenv_safe_string(), QwDetectorInfo::GetNumberOfElements(), QwDetectorInfo::GetRegion(), QwOptions::GetValue(), gQwOptions, QwGeometry::in(), QwTrackingTree::inittree(), kDirectionV, kDirectionX, kPackage1, kRegionID2, kRegionID3, QwTrackingTreeRegion::PrintNodes(), QwDebug, QwTrackingTree::SetGeometry(), and QwTrackingTree::SetMaxSlope().

Referenced by QwTrackingWorker().

330 {
331  // Loop over packages
332  for (EQwDetectorPackage package = kPackage1;
333  package <= kPackage2; package++) {
334  // Loop over regions
335  for (EQwRegionID region = kRegionID2;
336  region <= kRegionID3; region++) {
337  // Loop over directions
338  for (EQwDirectionID direction = kDirectionX;
339  direction <= kDirectionV; direction++) {
340 
341  // Get the detectors, and skip if none are active
342  QwGeometry detectors(geometry.in(region).in(package).in(direction));
343  if (detectors.size() == 0) continue;
344 
345  // Get the first plane
346  QwDetectorInfo* det = detectors.front();
347 
348  int levels = 0;
349  int numlayers = 0;
350  double width = 0.0;
351 
352  /// Region 2 contains 4 layers
353  if (det->GetRegion() == kRegionID2) {
354  numlayers = 4; // TODO replace this with info from the geometry file
355  width = det->GetElementSpacing() * det->GetNumberOfElements();
356  levels = fLevelsR2;
357  }
358 
359  /// Region 3 contains 8 layers
360  if (det->GetRegion() == kRegionID3) {
361  numlayers = 8; // TODO replace this with info from the geometry file
362  width = det->GetActiveWidthZ();
363  levels = fLevelsR3;
364  }
365 
366  /// Create a new search tree
367  QwTrackingTree thetree(numlayers);
368 
369  /// Set the maximum slope
370  thetree.SetMaxSlope(gQwOptions.GetValue<float>("QwTracking.R2.maxslope"));
371 
372  /// Set the geometry
373  thetree.SetGeometry(detectors);
374 
375  /// Set up the filename with the following format
376  /// tree[numlayers]-[levels]-[u|l]-[1|2|3]-[n|u|v|x|y].tre
377  std::stringstream filename;
378  filename << getenv_safe_string("QW_SEARCHTREE")
379  << "/tree" << numlayers << "-" << levels
380  << "-" << package << "-" << region << "-" << direction << ".tre";
381  QwDebug << "Tree filename: " << filename.str() << QwLog::endl;
382 
383  /// Each element of search tree will point to a pattern database
384  QwTrackingTreeRegion* searchtree = thetree.inittree(filename.str(),
385  levels,
386  numlayers,
387  width,
388  det,
389  fRegenerate);
390  if (fPrintPatternDatabase && region == kRegionID2) searchtree->PrintNodes();
391  if (fPrintPatternDatabase && region == kRegionID3) searchtree->PrintTrees();
392  // Add tree to internal list for deletion
393  fSearchTrees.push_back(searchtree);
394 
395  // Set the detector identification
396  searchtree->SetRegion(region);
397  searchtree->SetPackage(package);
398  searchtree->SetDirection(direction);
399 
400  // Store pointer to tree in detector objects
401  for (size_t i = 0; i < detectors.size(); i++) {
402  detectors.at(i)->SetTrackingSearchTree(searchtree);
403  }
404 
405  } // end of loop over directions
406  } // end of loop over regions
407  } // end of loop over packages
408 }
const QwGeometry in(const EQwRegionID &r) const
Get detectors in given region.
Definition: QwGeometry.h:92
bool fPrintPatternDatabase
Print the pattern database.
A container for the pattern databases for each detector region.
void PrintNodes() const
Print the list of nodes.
EQwRegionID GetRegion() const
int fLevelsR2
Region 2 levels.
QwOptions gQwOptions
Definition: QwOptions.cc:27
EQwDetectorPackage
Definition: QwTypes.h:70
std::vector< QwTrackingTreeRegion * > fSearchTrees
Internal list of search trees created by QwTrackingWorker.
T GetValue(const std::string &key)
Get a templated value.
Definition: QwOptions.h:240
EQwRegionID
Definition: QwTypes.h:16
#define QwDebug
Predefined log drain for debugging output.
Definition: QwLog.h:60
Creates and manages the treesearch pattern database.
double GetActiveWidthZ() const
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
double GetElementSpacing() const
EQwDirectionID
Definition: QwTypes.h:41
const std::string getenv_safe_string(const char *name)
Definition: QwOptions.h:37
Collection of QwDetectorInfo pointers that specifies an experimental geometry.
Definition: QwGeometry.h:27
int fLevelsR3
Region 3 levels.
int GetNumberOfElements() const
bool fRegenerate
Regenerate the search tree.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwTrackingWorker::ProcessEvent ( const QwSubsystemArrayTracking detectors,
QwEvent event 
)

Process the hit list and construct the event.

ProcessHits() - this function is the main tracking function. It performs tracking in projections (u/v/x) to form treelines, it combines projection tracks to tracks in space and bridges tracks in space before and after the magnet to form recon- structed particle tracks. The function calculates the momentum of the tracks and afterwards reiterates the track parameters using the now known track position and momentum for iterating the hit positions in the tracking chambers.

inputs: (1) QwHitContainer &hitlist - pointer to the QwHitContainer object with the hit list.

outputs: (1) Event* ProcessHits() - pointer to the structure with all of the reconstruction information for this event.

If no hits, return

If tracking is disabled, return

Loop through all detector packages

Loop through the detector regions

Loop through the detector directions

-— 1st: check that the direction is tree-searchable -—

-— 2nd: do some variable initialization for the found tree line linked list -—

-— 3rd: create the bit patterns for the hits -—

Get the subhitlist of hits in this detector

-— TASK 2: Combine the treelines into partial tracks -—

-— TASK 3: Sort out the Partial Tracks -—

-— TASK 4: Hook up the partial track info to the event info -—

If momentum reconstruction is disabled

Definition at line 431 of file QwTrackingWorker.cc.

References QwHitContainer::Append(), QwRayTracer::Bridge(), QwLog::endl(), fBridgingTrackFilter, fDebug, fDisableMomentum, fDisableRayTracer, fDisableTracking, fGeometry, QwBridgingTrackFilter::Filter(), fLevelsR2, fLevelsR3, fMismatchPkg, fRayTracer, fShowEventPattern, fTreeCombine, fTreeMatch, fTreeSearch, fTreeSort, QwDetectorInfo::GetNumberOfElements(), QwEvent::GetNumberOfTreeLines(), QwDetectorInfo::GetPlane(), QwHitContainer::GetSubList_Plane(), QwEvent::GetTrack(), QwTrackingTreeRegion::GetWidth(), QwGeometry::in(), QwDetectorInfo::IsInactive(), QwTrackingTreeRegion::IsSearchable(), kDirectionU, kDirectionV, kDirectionX, kPackage1, kPackage2, kRegionID2, kRegionID3, QwTrackingTreeMatch::MatchRegion3(), MAX_LAYERS, nbad, nbridged, QwTreeLine::next, ngood, QwHitContainer::Print(), QwTreeLine::Print(), QwDebug, QwMessage, QwOut, QwVerbose, QwWarning, R2Bad, R2Good, R3Bad, R3Good, QwTrackingTreeSort::rcPartConnSort(), QwTrackingTreeSort::rcTreeConnSort(), QwTrackingTreeSearch::SearchTreeLines(), VQwTrackingElement::SetDirection(), QwHitPattern::SetHDCHitList(), VQwTrackingElement::SetPackage(), VQwTrackingElement::SetPlane(), VQwTrackingElement::SetRegion(), QwHitPattern::SetVDCHit(), QwTrackingTreeCombine::TlTreeCombine(), and QwTrackingTreeCombine::TlTreeLineSort().

Referenced by main().

433 {
434  // Get hit list from event
435  QwHitContainer* hitlist = event->GetHitContainer();
436 
437  // Print hitlist
438  if (fDebug) hitlist->Print();
439 
440  /// If no hits, return
441  if (hitlist->size() == 0) {
442  delete hitlist;
443  return;
444  }
445 
446  /// If tracking is disabled, return
447  if (fDisableTracking) {
448  delete hitlist;
449  return;
450  }
451 
452 
453  /// Loop through all detector packages
454  // Normally only two of these will generate actual tracks,
455  // since there are only two tracking packages for the eight
456  // possible positions
457  for (EQwDetectorPackage package = kPackage1;
458  package <= kPackage2;
459  package++) {
460  QwDebug << "[QwTrackingWorker::ProcessHits] Package: " << package << QwLog::endl;
461 
462  /// Loop through the detector regions
463  for (EQwRegionID region = kRegionID2;
464  region <= kRegionID3;
465  region++) {
466  QwDebug << "[QwTrackingWorker::ProcessHits] Region: " << region << QwLog::endl;
467 
468 
469  int dlayer = 0; // number of detector planes in the search
470  int tlayers = 0; // number of tracking layers
471 
472  /// Loop through the detector directions
473  for (EQwDirectionID dir = kDirectionX;
474  dir <= kDirectionV;
475  dir++) {
476  QwDebug << "[QwTrackingWorker::ProcessHits] Direction: " << dir << QwLog::endl;
477 
478  // Find these detectors
479  QwGeometry detectors(fGeometry.in(region).in(package).in(dir));
480  if (detectors.size() == 0) continue;
481 
482  /*! ---- 1st: check that the direction is tree-searchable ---- */
483 
484  // The search tree for this detector will be stored in searchtree
485  QwTrackingTreeRegion* searchtree = detectors.at(0)->GetTrackingSearchTree();
486 
487  // Check whether the search tree exists
488  if (! searchtree) {
489  QwDebug << "[QwTrackingWorker::ProcessHits] No such tree!" << QwLog::endl;
490  continue;
491  }
492 
493  // Check whether the search tree is searchable
494  if (! searchtree->IsSearchable()) {
495  event->fTreeLine[package][region][dir] = 0;
496  QwDebug << "[QwTrackingWorker::ProcessHits] Search tree not searchable!" << QwLog::endl;
497  continue;
498  }
499 
500 
501  /*! ---- 2nd: do some variable initialization for the found tree line
502  linked list ---- */
503 
504  // Start the search for this set of like-pitched planes
505  QwTreeLine* treelinelist = 0; // local list of found tree lines
506 
507  /*! ---- 3rd: create the bit patterns for the hits ---- */
508 
509  /* Region 3 VDC */
510  if (region == kRegionID3) {
511 
512  dlayer = 0; /* set "number of detectors" to zero */
513  QwTreeLine *treelinelist1 = 0, *treelinelist2 = 0;
514 
515  // Set hit pattern for this region
516  QwDebug << "Setting hit pattern (region 3)" << QwLog::endl;
517 
518  /* Loop over the like-pitched planes in a region */
519  for (std::vector<QwDetectorInfo*>::iterator iter = detectors.begin();
520  iter != detectors.end(); iter++) {
521  QwDetectorInfo* det = (*iter);
522 
523  // Get plane number
524  int plane = det->GetPlane();
525  // Get number of wires
526  int numwires = det->GetNumberOfElements();
527 
528  // Print detector info
529  if (fDebug) QwOut << *det << QwLog::endl;
530 
531  // If detector is inactive for tracking, skip it
532  if (det->IsInactive()) continue;
533 
534  /// Get the subhitlist of hits in this detector
535  QwHitContainer *subhitlist = hitlist->GetSubList_Plane(region, package, plane);
536  // If no hits in this detector, skip to the next detector.
537  if (subhitlist->size() == 0) {
538  delete subhitlist;
539  continue;
540  }
541  // Print the hit list for this detector
542  if (fDebug) {
543  std::cout << "region: " << region << " package: " << package << " plane: " << plane << std::endl;
544  subhitlist->Print();
545  }
546 
547  // Create a vector of hit patterns
548  std::vector<QwHitPattern> patterns(numwires+1,QwHitPattern(fLevelsR3)); // wires are counted from 1
549 
550  // Loop over the hits in the subhitlist
551  for (QwHitContainer::iterator hit = subhitlist->begin();
552  hit != subhitlist->end();
553  hit++) {
554 
555  // Construct the hit pattern for this set of hits
556  QwHitPattern hitpattern(fLevelsR3);
557  hitpattern.SetVDCHit(searchtree->GetWidth(), &(*hit));
558  // Add hit pattern to the vector
559  int wire = hit->GetElement();
560  patterns.at(wire) += hitpattern;
561 
562  } // end of loop over hits in this event
563 
564  // Print hit pattern, if requested
565  if (fShowEventPattern) {
566  QwOut << "event pattern: " << QwLog::endl;
567  for (size_t wire = 0; wire < patterns.size(); wire++)
568  if (patterns.at(wire).HasHits())
569  QwMessage << wire << ":" << patterns.at(wire) << QwLog::endl;
570  }
571 
572  // Copy the new hit patterns into the old array structure
573  // TODO This is temporary
574  char** channel = new char*[patterns.size()];
575  int** hashchannel = new int*[patterns.size()];
576  for (size_t wire = 0; wire < patterns.size(); wire++) {
577  channel[wire] = new char[patterns.at(wire).GetNumberOfBins()];
578  hashchannel[wire] = new int[patterns.at(wire).GetFinestBinWidth()];
579  patterns.at(wire).GetPattern(channel[wire]);
580  patterns.at(wire).GetPatternHash(hashchannel[wire]);
581  }
582 
583  // All hits in this detector (VDC) are added to the bit pattern.
584  // We can start the tree search now.
585  // NOTE Somewhere around here a memory leak lurks
586  QwDebug << "Searching for matching patterns (direction " << dir << ")" << QwLog::endl;
587  treelinelist = fTreeSearch->SearchTreeLines(searchtree,
588  channel, hashchannel, fLevelsR3,
589  numwires, MAX_LAYERS);
590 
591  // Delete the old array structures
592  for (size_t wire = 0; wire < patterns.size(); wire++) {
593  delete[] channel[wire];
594  delete[] hashchannel[wire];
595  }
596  delete[] channel;
597  delete[] hashchannel;
598 
599  // TODO These treelines should contain the region id etc
600  // We should set the QwDetectorInfo link here already,
601  // or manually set the QwDetectorID fields. Also, we
602  // should put QwDetectorInfo in the searchtree object.
603  for (QwTreeLine* treeline = treelinelist;
604  treeline; treeline = treeline->next) {
605  treeline->SetRegion(region);
606  treeline->SetPackage(package);
607  treeline->SetDirection(dir);
608  treeline->SetPlane(plane);
609  }
610 
611  // Print list of tree lines
612  if (fDebug) {
613  cout << "List of treelines:" << endl;
614  if (treelinelist) treelinelist->Print();
615  }
616 
617  QwDebug << "Calculate chi^2" << QwLog::endl;
618  double width = searchtree->GetWidth();
619  fTreeCombine->TlTreeLineSort(treelinelist, subhitlist,
620  package, region, dir,
621  1UL << (fLevelsR3 - 1), 0, dlayer, width);
622 
623  QwDebug << "Sort patterns" << QwLog::endl;
624  fTreeSort->rcTreeConnSort (treelinelist, region);
625 
626  if (plane < 3)
627  treelinelist1 = treelinelist;
628  else
629  treelinelist2 = treelinelist;
630 
631  dlayer++;
632 
633  // Delete subhitlist
634  delete subhitlist;
635 
636  } // end of loop over like-pitched planes in a region
637  event->AddTreeLineList(treelinelist1);
638  event->AddTreeLineList(treelinelist2);
639  treelinelist = 0;
640  // treelinelist 1 and treelinelist 2 are in the same dir
641  QwDebug << "Matching region 3 segments" << QwLog::endl;
642  if (treelinelist1 && treelinelist2) {
643  treelinelist = fTreeMatch->MatchRegion3 (treelinelist1, treelinelist2);
644  event->fTreeLine[package][region][dir] = treelinelist;
645  event->AddTreeLineList(treelinelist);
646 
647  if (fDebug) {
648  cout << "VDC1:" << endl;
649  if (treelinelist1) treelinelist1->Print();
650  cout << "VDC2:" << endl;
651  if (treelinelist2) treelinelist2->Print();
652  }
653  if (fDebug) {
654  cout << "VDC1+2:" << endl;
655  if (treelinelist) treelinelist->Print();
656  }
657 
658  }
659 
660  // Delete tree line lists after storing results in event structure
661  QwTreeLine* tl = 0;
662  QwTreeLine* tl_next = 0;
663  tl = treelinelist1;
664  while (tl) { tl_next = tl->next; delete tl; tl = tl_next; }
665  tl = treelinelist2;
666  while (tl) { tl_next = tl->next; delete tl; tl = tl_next; }
667 
668  tlayers = MAX_LAYERS; /* remember the number of tree-detector */
669 
670 
671  /* Region 2 HDC */
672  } else if (region == kRegionID2) {
673 
674  // Set hit pattern for this region
675  QwDebug << "Setting hit pattern (region 2)" << QwLog::endl;
676 
677  // Create a vector of hit patterns
678  std::vector<QwHitPattern> patterns;
679 
680  // Get the hit list for this package/region/direction
681  QwHitContainer *treelinehits = new QwHitContainer();
682 
683  /* Loop over the like-pitched planes in a region */
684  tlayers = 0;
685  for (std::vector<QwDetectorInfo*>::iterator iter = detectors.begin();
686  iter != detectors.end(); iter++, tlayers++) {
687  QwDetectorInfo* det = (*iter);
688 
689  // Get plane number
690  Int_t plane = det->GetPlane();
691 
692  // Print detector info
693  if (fDebug) QwOut << plane << ": " << *det << QwLog::endl;
694 
695  // Get the subhitlist of hits in this detector plane
696  QwHitContainer *planehits = hitlist->GetSubList_Plane(region, package, plane);
697  treelinehits->Append(planehits);
698  if (fDebug) planehits->Print();
699  // Construct the hit pattern for this set of hits
700  QwHitPattern hitpattern(fLevelsR2);
701  // Set the detector identification
702  hitpattern.SetRegion(region);
703  hitpattern.SetPackage(package);
704  hitpattern.SetDirection(dir);
705  hitpattern.SetPlane(plane);
706  // Set the hit pattern
707  hitpattern.SetHDCHitList(searchtree->GetWidth(), planehits);
708  // Add hit pattern to the vector
709  patterns.push_back(hitpattern);
710 
711  // Delete subhitlist
712  delete planehits;
713 
714  } // end of loop over like-pitched planes in a region
715 
716  // Print hit pattern, if requested
717  if (fShowEventPattern)
718  for (size_t layer = 0; layer < patterns.size(); layer++)
719  QwMessage << patterns.at(layer) << QwLog::endl;
720 
721  // Copy the new hit patterns into the old array structure
722  // TODO This is temporary
723  int levels = patterns.at(0).GetNumberOfLevels();
724  char** channel = new char*[patterns.size()];
725  int** hashchannel = new int*[patterns.size()];
726  for (size_t layer = 0; layer < patterns.size(); layer++) {
727  channel[layer] = new char[patterns.at(layer).GetNumberOfBins()];
728  hashchannel[layer] = new int[patterns.at(layer).GetFinestBinWidth()];
729  patterns.at(layer).GetPattern(channel[layer]);
730  patterns.at(layer).GetPatternHash(hashchannel[layer]);
731  }
732 
733  QwDebug << "Search for matching patterns (direction " << dir << ")" << QwLog::endl;
734  treelinelist = fTreeSearch->SearchTreeLines(searchtree,
735  channel, hashchannel, levels,
736  0, tlayers);
737 
738  // Delete the old array structures
739  for (size_t layer = 0; layer < patterns.size(); layer++) {
740  delete[] channel[layer];
741  delete[] hashchannel[layer];
742  }
743  delete[] channel;
744  delete[] hashchannel;
745 
746  // TODO These treelines should contain the region id etc
747  // We should set the QwDetectorInfo link here already,
748  // or manually set the QwDetectorID fields. Also, we
749  // should put QwDetectorInfo in the searchtree object.
750  for (QwTreeLine* treeline = treelinelist;
751  treeline; treeline = treeline->next) {
752  treeline->SetRegion(region);
753  treeline->SetPackage(package);
754  treeline->SetDirection(dir);
755  }
756 
757  QwDebug << "Calculate chi^2" << QwLog::endl;
758  if (searchtree) {
759 
760  double width = searchtree->GetWidth();
761  fTreeCombine->TlTreeLineSort (treelinelist, treelinehits,
762  package, region, dir,
763  1UL << (fLevelsR2 - 1),
764  tlayers, 0, width);
765  }
766 
767  //QwDebug << "Sort patterns" << QwLog::endl;
768  //fTreeSort->rcTreeConnSort (treelinelist, region);
769 
770  if (fDebug) {
771  cout << "List of treelines:" << endl;
772  if (treelinelist) treelinelist->Print();
773  }
774  event->fTreeLine[package][region][dir] = treelinelist;
775  event->AddTreeLineList(treelinelist);
776 
777  // Delete treelinehits
778  delete treelinehits;
779 
780  /* Any other region */
781  } else {
782  QwWarning << "[QwTrackingWorker::ProcessHits] Warning: no support for this detector." << QwLog::endl;
783  return;
784  }
785 
786 
787  } // end of loop over the three wire-pitch directions
788 
789 
790  /*! ---- TASK 2: Combine the treelines into partial tracks ---- */
791 
792  // List of partial tracks to return
793  std::vector<QwPartialTrack*> parttracklist;
794 
795  if (event->GetNumberOfTreeLines() > 0 && tlayers)
796  {
797  std::vector<QwTreeLine*> treelines_x =
798  event->fTreeLine[package][region][kDirectionX]->GetListAsVector();
799  std::vector<QwTreeLine*> treelines_u =
800  event->fTreeLine[package][region][kDirectionU]->GetListAsVector();
801  std::vector<QwTreeLine*> treelines_v =
802  event->fTreeLine[package][region][kDirectionV]->GetListAsVector();
803 
804  parttracklist = fTreeCombine->TlTreeCombine(
805  treelines_x, treelines_u, treelines_v,
806  package, region,
807  tlayers, dlayer);
808  } else continue;
809 
810 
811  // If we found partial tracks in this event
812  if (parttracklist.size() > 0) {
813 
814  /*! ---- TASK 3: Sort out the Partial Tracks ---- */
815 
816  fTreeSort->rcPartConnSort(parttracklist);
817 
818 
819  /*! ---- TASK 4: Hook up the partial track info to the event info ---- */
820 
821  // Add partial tracks to the event
822  event->AddPartialTrackList(parttracklist);
823 
824  // Debug output
825  if (fDebug) {
826  QwOut << "List of partial tracks:" << QwLog::endl;
827  for (size_t pt = 0; pt < parttracklist.size(); pt++)
828  QwOut << *parttracklist[pt] << QwLog::endl;
829  }
830 
831  // Delete partial tracks
832  for (size_t pt = 0; pt < parttracklist.size(); pt++)
833  delete parttracklist[pt];
834 
835 
836  // Count this as a good event
837  QwVerbose << "Found a good partial track in region " << region << QwLog::endl;
838  if (region == 2) ++R2Good;
839  if (region == 3) ++R3Good;
840  ngood++;
841 
842  } else {
843 
844  // Count this as a bad event
845  QwVerbose << "Couldn't find a good partial track in region " << region << QwLog::endl;
846  if (region == 2) ++R2Bad;
847  if (region == 3) ++R3Bad;
848  nbad++;
849  }
850 
851  } /* end of loop over the regions */
852 
853  } /* end of loop over the detector packages */
854 
855  // Delete local copy of the hit list
856  delete hitlist;
857 
858 
859  /* ==============================
860  * Correlate front and back
861  * tracks from x, y and y' infor-
862  * mation
863  * ============================== */
864 
865  /// If momentum reconstruction is disabled
866  if (fDisableMomentum) return;
867 
868  // Loop over packages in region 3
869  for (EQwDetectorPackage R3package = kPackage1;
870  R3package <= kPackage2; R3package++) {
871 
872  // Determine package in region 2 (considering possibility of reversed run)
873  EQwDetectorPackage R2package;
874  if (fMismatchPkg)
875  R2package = (R3package == kPackage1)? kPackage2: kPackage1;
876  else
877  R2package = (R3package == kPackage1)? kPackage1: kPackage2;
878 
879  // Get the lists of partial tracks in the front and back detectors
880  std::vector<QwPartialTrack*> frontlist = event->GetListOfPartialTracks(kRegionID2,R2package);
881  std::vector<QwPartialTrack*> backlist = event->GetListOfPartialTracks(kRegionID3,R3package);
882 
883  // Loop over all good front and back partial tracks
884  for (size_t ifront = 0; ifront < frontlist.size(); ifront++) {
885  QwPartialTrack* front = frontlist[ifront];
886 
887  for (size_t iback = 0; iback < backlist.size(); iback++) {
888  QwPartialTrack* back = backlist[iback];
889 
890  // Filter reasonable pairs
891  int status = fBridgingTrackFilter->Filter(front, back);
892  status = 0;
893  if (status != 0) {
894  QwMessage << "Tracks did not pass filter." << QwLog::endl;
895  continue;
896  }
897 
898  // Attempt to bridge tracks using ray-tracing
899  if (! fDisableRayTracer) {
900  const QwTrack* track = fRayTracer->Bridge(front, back);
901  if (track) {
902  event->AddTrack(track);
903  delete track;
904  }
905  }
906 
907  } // end of loop over back tracks
908 
909  } // end of loop over front tracks
910 
911  } // end of loop over packages
912 
913 
914  // Calculate kinematics
915  int num_of_bridged = event->GetNumberOfTracks();
916  if ( num_of_bridged > 0) {
917  nbridged += num_of_bridged;
918  event->CalculateKinematics(event->GetTrack(0));
919  }
920 
921  // Print the result
922  if (fDebug) event->Print();
923 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
bool IsSearchable() const
Is this region searchable?
#define QwOut
Predefined log drain for explicit output.
Definition: QwLog.h:35
bool fDisableTracking
Disable all tracking.
const QwGeometry in(const EQwRegionID &r) const
Get detectors in given region.
Definition: QwGeometry.h:92
Int_t GetNumberOfTreeLines() const
Get the number of tree lines.
Definition: QwEvent.h:254
int GetPlane() const
void SetRegion(EQwRegionID region)
Set the region.
A container for the pattern databases for each detector region.
Hit patterns used in the tracking tree search.
Definition: QwHitPattern.h:44
#define QwVerbose
Predefined log drain for verbose messages.
Definition: QwLog.h:55
int fLevelsR2
Region 2 levels.
int rcTreeConnSort(QwTreeLine *head, EQwRegionID region)
The best (by chi^2) treelines are select.
int nbad
number of bad events
void TlTreeLineSort(QwTreeLine *tl, QwHitContainer *hl, EQwDetectorPackage package, EQwRegionID region, EQwDirectionID dir, unsigned long bins, int tlayer, int dlayer, double width)
EQwDetectorPackage
Definition: QwTypes.h:70
QwTrackingTreeSort * fTreeSort
Module that sorts lists of treelines and partial tracks.
EStatus Filter(const QwPartialTrack *front, const QwPartialTrack *back) const
Filter front and back track combinations.
QwRayTracer * fRayTracer
Ray tracing bridging method.
int fDebug
Debug level.
QwTrackingTreeMatch * fTreeMatch
Module that matches up VDC front and back treelines.
std::vector< QwPartialTrack * > TlTreeCombine(const std::vector< QwTreeLine * > &treelines_x, const std::vector< QwTreeLine * > &treelines_u, const std::vector< QwTreeLine * > &treelines_v, EQwDetectorPackage package, EQwRegionID region, int tlayer, int dlayer)
EQwRegionID
Definition: QwTypes.h:16
#define QwDebug
Predefined log drain for debugging output.
Definition: QwLog.h:60
QwTrackingTreeSearch * fTreeSearch
Module that handles the tree search.
QwTrackingTreeCombine * fTreeCombine
Module that combines treelines and partial tracks.
void Print(const Option_t *options=0) const
Definition: QwTreeLine.cc:317
Contains the complete track as a concatenation of partial tracks.
Definition: QwTrack.h:30
int rcPartConnSort(std::vector< QwPartialTrack * > head)
const QwTrack * Bridge(const QwPartialTrack *front, const QwPartialTrack *back)
Bridge from the front to back partial track.
Definition: QwRayTracer.cc:178
QwBridgingTrackFilter * fBridgingTrackFilter
Track filter.
QwTreeLine * MatchRegion3(QwTreeLine *frontlist, QwTreeLine *backlist)
Match the tree lines in two planes in region 3.
QwGeometry fGeometry
Detector geometry.
bool fShowEventPattern
Show all event patterns.
One-dimensional (u, v, or x) track stubs and associated hits.
Definition: QwTreeLine.h:51
const QwTrack * GetTrack(const int t) const
Get the specified track.
Definition: QwEvent.cc:545
int ngood
number of good events
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
bool fDisableRayTracer
Disable ray tracer momentum reconstruction.
void Print(const Option_t *option=0) const
EQwDirectionID
Definition: QwTypes.h:41
Collection of QwDetectorInfo pointers that specifies an experimental geometry.
Definition: QwGeometry.h:27
double GetWidth() const
Get the width.
int fLevelsR3
Region 3 levels.
QwTreeLine * next
Definition: QwTreeLine.h:231
QwTreeLine * SearchTreeLines(QwTrackingTreeRegion *searchtree, char *pattern[16], int *hashpat[16], int maxlevel, int numwires, int numlayers)
Search for the tree lines consistent with the hit pattern.
int nbridged
number of beidged tracks
#define MAX_LAYERS
Definition: globals.h:8
#define QwWarning
Predefined log drain for warnings.
Definition: QwLog.h:45
void GetSubList_Plane(EQwRegionID region, EQwDetectorPackage package, Int_t plane, std::vector< QwHit > &sublist)
void Append(const QwHitContainer &mylist)
int GetNumberOfElements() const
Contains the straight part of a track in one region only.
bool IsInactive() const
bool fMismatchPkg
Indicates if the pkg for R2 and R3 is differnt at the same octant.
bool fDisableMomentum
Disable momentum reconstruction.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwTrackingWorker::ProcessOptions ( QwOptions options)

Process command line and config file options.

Definition at line 296 of file QwTrackingWorker.cc.

References fDebug, fDisableMatrixLookup, fDisableMomentum, fDisableRayTracer, fDisableTracking, fFilenameLookupTable, fLevelsR2, fLevelsR3, fMismatchPkg, fPrintPatternDatabase, fRegenerate, fShowEventPattern, fShowMatchingPattern, and QwOptions::GetValue().

Referenced by QwTrackingWorker().

297 {
298  // Enable tracking debug flag
299  fDebug = options.GetValue<int>("QwTracking.debug");
300  fRegenerate = options.GetValue<bool>("QwTracking.regenerate");
301 
302  // Disable tracking and/or momentu reconstruction
303  fDisableTracking = options.GetValue<bool>("QwTracking.disable-tracking");
304  fDisableMomentum = options.GetValue<bool>("QwTracking.disable-momentum");
305  fMismatchPkg=options.GetValue<bool>("QwTracking.package-mismatch");
306 
307  // Set the flags for printing the pattern database
308  fPrintPatternDatabase = options.GetValue<bool>("QwTracking.print-pattern-db");
309 
310  // Set the flags for showing the matching event patterns
311  fShowEventPattern = options.GetValue<bool>("QwTracking.showeventpattern");
312  fShowMatchingPattern = options.GetValue<bool>("QwTracking.showmatchingpattern");
313 
314  // Enable/disable the lookup table and raytracer momentum reconstruction methods
315  fDisableMatrixLookup = options.GetValue<bool>("QwTracking.disable-matrixlookup");
316  fDisableRayTracer = options.GetValue<bool>("QwTracking.disable-raytracer");
317 
318  // Lookup table filename
319  fFilenameLookupTable = options.GetValue<std::string>("QwTracking.lookuptable");
320 
321  // Number of levels (search tree depth) for region 2
322  fLevelsR2 = options.GetValue<int>("QwTracking.R2.levels");
323  // Number of levels (search tree depth) for region 3
324  fLevelsR3 = options.GetValue<int>("QwTracking.R3.levels");
325 }
bool fDisableTracking
Disable all tracking.
bool fPrintPatternDatabase
Print the pattern database.
int fLevelsR2
Region 2 levels.
int fDebug
Debug level.
bool fDisableMatrixLookup
Disable matrix lookup momentum reconstruction.
T GetValue(const std::string &key)
Get a templated value.
Definition: QwOptions.h:240
bool fShowMatchingPattern
Show matching event patterns.
bool fShowEventPattern
Show all event patterns.
bool fDisableRayTracer
Disable ray tracer momentum reconstruction.
std::string fFilenameLookupTable
Filename of the lookup table in QW_LOOKUP.
int fLevelsR3
Region 3 levels.
bool fMismatchPkg
Indicates if the pkg for R2 and R3 is differnt at the same octant.
bool fRegenerate
Regenerate the search tree.
bool fDisableMomentum
Disable momentum reconstruction.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwTrackingWorker::SetDebugLevel ( const int  debug)
inline

Set the debug level.

Definition at line 73 of file QwTrackingWorker.h.

References fDebug.

73 { fDebug = debug; };
int fDebug
Debug level.
void QwTrackingWorker::SetGeometry ( const QwGeometry geometry)
inline

Set the geometry.

Definition at line 78 of file QwTrackingWorker.h.

References fGeometry.

Referenced by QwTrackingWorker().

78 { fGeometry = geometry; };
QwGeometry fGeometry
Detector geometry.

+ Here is the caller graph for this function:

void QwTrackingWorker::SetMagneticFieldCurrent ( const double  current)
inline

Set the magnetic field current.

Definition at line 86 of file QwTrackingWorker.h.

References fRayTracer, and QwRayTracer::SetMagneticFieldCurrent().

Referenced by main().

86  {
88  }
QwRayTracer * fRayTracer
Ray tracing bridging method.
void SetMagneticFieldCurrent(const double current)
Set magnetic field current.
Definition: QwRayTracer.h:66

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Field Documentation

QwBridgingTrackFilter* QwTrackingWorker::fBridgingTrackFilter
private

Track filter.

Definition at line 115 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), QwTrackingWorker(), and ~QwTrackingWorker().

int QwTrackingWorker::fDebug
private

Debug level.

Definition at line 99 of file QwTrackingWorker.h.

Referenced by GetDebugLevel(), ProcessEvent(), ProcessOptions(), QwTrackingWorker(), and SetDebugLevel().

bool QwTrackingWorker::fDisableMatrixLookup
private

Disable matrix lookup momentum reconstruction.

Definition at line 128 of file QwTrackingWorker.h.

Referenced by ProcessOptions().

bool QwTrackingWorker::fDisableMomentum
private

Disable momentum reconstruction.

Definition at line 127 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), ProcessOptions(), and QwTrackingWorker().

bool QwTrackingWorker::fDisableRayTracer
private

Disable ray tracer momentum reconstruction.

Definition at line 129 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), ProcessOptions(), and QwTrackingWorker().

bool QwTrackingWorker::fDisableTracking
private

Disable all tracking.

Definition at line 123 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), ProcessOptions(), and QwTrackingWorker().

std::string QwTrackingWorker::fFilenameFieldmap
private

Filename of the fieldmap in QW_FIELDMAP.

Definition at line 131 of file QwTrackingWorker.h.

std::string QwTrackingWorker::fFilenameLookupTable
private

Filename of the lookup table in QW_LOOKUP.

Definition at line 132 of file QwTrackingWorker.h.

Referenced by ProcessOptions().

QwGeometry QwTrackingWorker::fGeometry
private

Detector geometry.

Definition at line 96 of file QwTrackingWorker.h.

Referenced by GetGeometry(), ProcessEvent(), and SetGeometry().

int QwTrackingWorker::fLevelsR2
private

Region 2 levels.

Definition at line 105 of file QwTrackingWorker.h.

Referenced by InitTree(), ProcessEvent(), and ProcessOptions().

int QwTrackingWorker::fLevelsR3
private

Region 3 levels.

Definition at line 107 of file QwTrackingWorker.h.

Referenced by InitTree(), ProcessEvent(), and ProcessOptions().

bool QwTrackingWorker::fMismatchPkg
private

Indicates if the pkg for R2 and R3 is differnt at the same octant.

Definition at line 130 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), and ProcessOptions().

bool QwTrackingWorker::fPrintPatternDatabase
private

Print the pattern database.

Definition at line 124 of file QwTrackingWorker.h.

Referenced by InitTree(), and ProcessOptions().

QwRayTracer* QwTrackingWorker::fRayTracer
private

Ray tracing bridging method.

Definition at line 117 of file QwTrackingWorker.h.

Referenced by GetMagneticFieldCurrent(), ProcessEvent(), QwTrackingWorker(), SetMagneticFieldCurrent(), and ~QwTrackingWorker().

bool QwTrackingWorker::fRegenerate
private

Regenerate the search tree.

Definition at line 102 of file QwTrackingWorker.h.

Referenced by InitTree(), and ProcessOptions().

std::vector<QwTrackingTreeRegion*> QwTrackingWorker::fSearchTrees
private

Internal list of search trees created by QwTrackingWorker.

Definition at line 146 of file QwTrackingWorker.h.

Referenced by InitTree(), and ~QwTrackingWorker().

bool QwTrackingWorker::fShowEventPattern
private

Show all event patterns.

Definition at line 125 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), and ProcessOptions().

bool QwTrackingWorker::fShowMatchingPattern
private

Show matching event patterns.

Definition at line 126 of file QwTrackingWorker.h.

Referenced by ProcessOptions(), and QwTrackingWorker().

QwTrackingTreeCombine* QwTrackingWorker::fTreeCombine
private

Module that combines treelines and partial tracks.

Definition at line 140 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), QwTrackingWorker(), and ~QwTrackingWorker().

QwTrackingTreeMatch* QwTrackingWorker::fTreeMatch
private

Module that matches up VDC front and back treelines.

Definition at line 142 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), QwTrackingWorker(), and ~QwTrackingWorker().

QwTrackingTreeSearch* QwTrackingWorker::fTreeSearch
private

Module that handles the tree search.

Definition at line 139 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), QwTrackingWorker(), and ~QwTrackingWorker().

QwTrackingTreeSort* QwTrackingWorker::fTreeSort
private

Module that sorts lists of treelines and partial tracks.

Definition at line 141 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), QwTrackingWorker(), and ~QwTrackingWorker().

int QwTrackingWorker::nbad

number of bad events

Definition at line 51 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::nbridged

number of beidged tracks

Definition at line 58 of file QwTrackingWorker.h.

Referenced by main(), ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::ngood

number of good events

Definition at line 50 of file QwTrackingWorker.h.

Referenced by main(), ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::R2Bad

Definition at line 54 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::R2Good

Definition at line 53 of file QwTrackingWorker.h.

Referenced by main(), ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::R3Bad

Definition at line 56 of file QwTrackingWorker.h.

Referenced by ProcessEvent(), and QwTrackingWorker().

int QwTrackingWorker::R3Good

Definition at line 55 of file QwTrackingWorker.h.

Referenced by main(), ProcessEvent(), and QwTrackingWorker().


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