FairRoot/PandaRoot
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
PndPatternMatcher Class Reference

#include <PndPatternMatcher.h>

Inheritance diagram for PndPatternMatcher:

Public Member Functions

 PndPatternMatcher ()
 
virtual ~PndPatternMatcher ()
 
void SetOutputBranchname (TString name)
 
void SetPersistence (Bool_t val)
 
void SetMatchRatio (float ratio)
 

Protected Member Functions

virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void FinishTask ()
 

Private Member Functions

void LoadPatternDB (TString filename)
 
void FindMatch ()
 
void CreateTrackCandFromMatch (std::set< int > partialMatch)
 

Private Attributes

PndGeoSttParfSttParameters
 
TClonesArray * fEventHeader
 
TClonesArray * fSTTHitArray
 
TFile * dbFile
 
TTree * patternTree
 
std::vector< std::set< int > > dbTubeIDs
 
std::vector< PndTrackCandfPartialCand
 
TClonesArray * fPartialMatchCand
 
TString fTrackCandName
 
Bool_t fPersistence
 
bool verboseMatches
 
float fMatchRatio
 

Detailed Description

Definition at line 18 of file PndPatternMatcher.h.

Constructor & Destructor Documentation

PndPatternMatcher::PndPatternMatcher ( )

Definition at line 21 of file PndPatternMatcher.cxx.

21  : fPartialMatchCand(0) {
22  dbFile = NULL;
23  patternTree = NULL;
24 // treeReader = NULL;
25  fSttParameters = NULL;
26  fEventHeader = NULL;
27  fSTTHitArray = NULL;
28 
29 
30  fTrackCandName = "PartialMatchCand";
31  fPersistence = kTRUE;
32 
33  verboseMatches = false;
34  fMatchRatio = 0.5;
35 }
TClonesArray * fEventHeader
PndGeoSttPar * fSttParameters
TClonesArray * fSTTHitArray
TClonesArray * fPartialMatchCand
PndPatternMatcher::~PndPatternMatcher ( )
virtual

Definition at line 37 of file PndPatternMatcher.cxx.

37  {
38 }

Member Function Documentation

void PndPatternMatcher::CreateTrackCandFromMatch ( std::set< int >  partialMatch)
private

Definition at line 125 of file PndPatternMatcher.cxx.

References PndTrackCand::AddHit(), fPartialCand, fSTTHitArray, and PndSttHit::GetTubeID().

Referenced by FindMatch().

125  {
126 // Create a PndTrackCand from the partial match found in FindMatch()
127  PndTrackCand trackCand;
128  int nSttHits = fSTTHitArray->GetEntriesFast();
129 // Loop over all STT hits
130  for (int iHit = 0; iHit < nSttHits; ++iHit) {
131  PndSttHit *sttHit = (PndSttHit*) fSTTHitArray->At(iHit);
132  int tubeID = sttHit->GetTubeID();
133 // verify if this hit is part of the partial match
134  auto search = partialMatch.find(tubeID);
135 // if hit is found in partial match, add it to track cand
136  if (search != partialMatch.end()) {
137  trackCand.AddHit(sttHit->GetEntryNr(),iHit);
138  }
139  }
140 // Add new trackCand to vector of all trackCands for this event
141  fPartialCand.push_back(trackCand);
142 }
std::vector< PndTrackCand > fPartialCand
void AddHit(UInt_t detId, UInt_t hitId, Double_t rho)
TClonesArray * fSTTHitArray
Int_t GetTubeID() const
Definition: PndSttHit.h:75
void PndPatternMatcher::Exec ( Option_t *  opt)
protectedvirtual

Definition at line 67 of file PndPatternMatcher.cxx.

References FindMatch(), fPartialCand, and fPartialMatchCand.

67  {
68  fPartialCand.clear();
69  fPartialMatchCand->Delete();
70  FindMatch();
71 
72 // Write trackCands into the data stream
73  for (auto const& trackCand : fPartialCand) {
74  new ((*fPartialMatchCand)[fPartialMatchCand->GetEntriesFast()]) PndTrackCand(trackCand);
75  }
76 }
std::vector< PndTrackCand > fPartialCand
TClonesArray * fPartialMatchCand
void PndPatternMatcher::FindMatch ( )
private

Definition at line 81 of file PndPatternMatcher.cxx.

References CreateTrackCandFromMatch(), dbTubeIDs, fMatchRatio, fSTTHitArray, PndSttHit::GetTubeID(), and verboseMatches.

Referenced by Exec().

81  {
82 // Get a list of all tubeIDs in the event
83  int nSttHits = fSTTHitArray->GetEntriesFast();
84  std::set<int> tubeIDs;
85  for (int iHit = 0; iHit < nSttHits; ++iHit) {
86  PndSttHit *sttHit = (PndSttHit*) fSTTHitArray->At(iHit);
87  int tubeID = sttHit->GetTubeID();
88  tubeIDs.insert(tubeID);
89  }
90 
91 // std::cout << dbTubeIDs.at(1).size() << std::endl;
92 
93 
94  int matches = 0;
95 // Compare each pattern in the data base with the current set of tubeIDs
96  for (auto const& dbPattern : dbTubeIDs) {
97  std::set<int> intersection;
98  intersection.clear();
99 // Determine the intersection between the current database pattern and the event's set of tubeIDs
100  std::set_intersection(tubeIDs.begin(), tubeIDs.end(), dbPattern.begin(), dbPattern.end(), std::inserter(intersection,intersection.begin()));
101 // Determine the ratio to which there is an overlap with the database pattern
102  float matchRatio = (float)intersection.size() / (float)dbPattern.size();
103 // If the matching ratio is large enough, proceed to create trackCand from intersection
104  if (matchRatio > fMatchRatio) {
105 // std::cout<< "match ratio: " << matchRatio << std::endl;
106  matches++;
107 
108 // Detailed output in case it is needed (set verboseMatches to true)
109  if (verboseMatches) {
110  if (intersection.size() > 0) {
111  std::cout << "matching tubes: ";
112  for (auto const& tube : intersection) {
113  std::cout << tube << " ";
114  }
115  std::cout << std::endl;
116  }
117  }
118 
119  CreateTrackCandFromMatch(intersection);
120 
121  }
122  }
123  if (verboseMatches) std::cout << "partial matches found: " << matches << std::endl;
124 }
TClonesArray * fSTTHitArray
Int_t GetTubeID() const
Definition: PndSttHit.h:75
std::vector< std::set< int > > dbTubeIDs
void CreateTrackCandFromMatch(std::set< int > partialMatch)
void PndPatternMatcher::FinishTask ( )
protectedvirtual

Definition at line 77 of file PndPatternMatcher.cxx.

77  {
78 
79 }
InitStatus PndPatternMatcher::Init ( )
protectedvirtual

Definition at line 44 of file PndPatternMatcher.cxx.

References fEventHeader, fPartialMatchCand, fPersistence, fSTTHitArray, fTrackCandName, and LoadPatternDB().

44  {
45  LoadPatternDB("patternDB.root");
46 
47 // Get the running instance of the FairRootManager to access tree branches
48  FairRootManager *ioman = FairRootManager::Instance();
49  if(!ioman) {
50  std::cout << "-E- PatternDBGenerator::Init: FairRootManager not instantiated!" << std::endl;
51  return kFATAL;
52  }
53 
54  fEventHeader = (TClonesArray*) ioman->GetObject("EventHeader.");
55  if (!fEventHeader) {
56  std::cout << "-E- PatternDBGenerator:Init: EventHeader not instantiated!" << std::endl;
57  }
58 
59 // Access the STTHit branch
60  fSTTHitArray = (TClonesArray*) ioman->GetObject("STTHit");
61 
62 // Set the output branch
63  fPartialMatchCand = ioman->Register(fTrackCandName,"PndTrackCand","STT",fPersistence);
64 
65  return kSUCCESS;
66 }
TClonesArray * fEventHeader
void LoadPatternDB(TString filename)
TClonesArray * fSTTHitArray
TClonesArray * fPartialMatchCand
void PndPatternMatcher::LoadPatternDB ( TString  filename)
private

Definition at line 143 of file PndPatternMatcher.cxx.

References dbFile, dbTubeIDs, PndPattern::GetTubeIDs(), and patternTree.

Referenced by Init().

143  {
144  dbFile = new TFile(filename,"READ");
145  patternTree = (TTree*) dbFile->Get("trackPatterns");
146 
147  dbTubeIDs.clear();
148 
149  PndPattern *pattern;
150  patternTree->SetBranchAddress("pattern", &pattern);
151 
152  Long64_t nPatterns = patternTree->GetEntriesFast();
153  for (Long64_t iEntry = 0; iEntry < nPatterns; ++iEntry) {
154  patternTree->GetEntry(iEntry);
155  std::set<int> tubeIDs = pattern->GetTubeIDs();
156  dbTubeIDs.push_back(tubeIDs);
157  }
158  std::cout << "# DB patterns: " << dbTubeIDs.size() << std::endl;
159 
160 }
std::set< int > GetTubeIDs() const
Definition: PndPattern.h:32
std::vector< std::set< int > > dbTubeIDs
const string filename
void PndPatternMatcher::SetMatchRatio ( float  ratio)
inline

Definition at line 25 of file PndPatternMatcher.h.

References fMatchRatio.

25 {fMatchRatio = ratio;}
void PndPatternMatcher::SetOutputBranchname ( TString  name)
inline

Definition at line 23 of file PndPatternMatcher.h.

References fTrackCandName, and name.

TString name
void PndPatternMatcher::SetParContainers ( )
protectedvirtual

Definition at line 40 of file PndPatternMatcher.cxx.

References fSttParameters, and rtdb.

40  {
41  FairRuntimeDb *rtdb = FairRunAna::Instance()->GetRuntimeDb();
42  fSttParameters = (PndGeoSttPar*) rtdb->getContainer("PndGeoSttPar");
43 }
PndGeoSttPar * fSttParameters
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
void PndPatternMatcher::SetPersistence ( Bool_t  val)
inline

Definition at line 24 of file PndPatternMatcher.h.

References fPersistence, and val.

Referenced by patternMatcher().

24 {fPersistence = val;}
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11

Member Data Documentation

TFile* PndPatternMatcher::dbFile
private

Definition at line 43 of file PndPatternMatcher.h.

Referenced by LoadPatternDB().

std::vector<std::set<int> > PndPatternMatcher::dbTubeIDs
private

Definition at line 45 of file PndPatternMatcher.h.

Referenced by FindMatch(), and LoadPatternDB().

TClonesArray* PndPatternMatcher::fEventHeader
private

Definition at line 40 of file PndPatternMatcher.h.

Referenced by Init().

float PndPatternMatcher::fMatchRatio
private

Definition at line 55 of file PndPatternMatcher.h.

Referenced by FindMatch(), and SetMatchRatio().

std::vector<PndTrackCand> PndPatternMatcher::fPartialCand
private

Definition at line 48 of file PndPatternMatcher.h.

Referenced by CreateTrackCandFromMatch(), and Exec().

TClonesArray* PndPatternMatcher::fPartialMatchCand
private

Definition at line 49 of file PndPatternMatcher.h.

Referenced by Exec(), and Init().

Bool_t PndPatternMatcher::fPersistence
private

Definition at line 52 of file PndPatternMatcher.h.

Referenced by Init(), and SetPersistence().

TClonesArray* PndPatternMatcher::fSTTHitArray
private

Definition at line 41 of file PndPatternMatcher.h.

Referenced by CreateTrackCandFromMatch(), FindMatch(), and Init().

PndGeoSttPar* PndPatternMatcher::fSttParameters
private

Definition at line 39 of file PndPatternMatcher.h.

Referenced by SetParContainers().

TString PndPatternMatcher::fTrackCandName
private

Definition at line 50 of file PndPatternMatcher.h.

Referenced by Init(), and SetOutputBranchname().

TTree* PndPatternMatcher::patternTree
private

Definition at line 44 of file PndPatternMatcher.h.

Referenced by LoadPatternDB().

bool PndPatternMatcher::verboseMatches
private

Definition at line 54 of file PndPatternMatcher.h.

Referenced by FindMatch().


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