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

#include <PndSttMatchTracks.h>

Inheritance diagram for PndSttMatchTracks:
PndPersistencyTask

Public Member Functions

 PndSttMatchTracks ()
 
 PndSttMatchTracks (Int_t verbose)
 
 PndSttMatchTracks (const char *name, const char *title="Pnd Stt Match Tracks Task", Int_t verbose=1)
 
virtual ~PndSttMatchTracks ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void AddHitCollectionName (char *hitCollectionName, char *pointCollectionName)
 
PndSttHitGetHitFromCollections (Int_t hitCounter)
 
FairMCPoint * GetPointFromCollections (Int_t hitCounter)
 
void SetPersistence (Bool_t persistence)
 
void SetPersistency (Bool_t val=kTRUE)
 
Bool_t GetPersistency ()
 

Private Member Functions

void AddAllCollections ()
 
void AddHitCollection (char const *collectionName, char const *pointCollectionName)
 
 PndSttMatchTracks (const PndSttMatchTracks &L)
 
PndSttMatchTracksoperator= (const PndSttMatchTracks &)
 
 ClassDef (PndSttMatchTracks, 1)
 

Private Attributes

TClonesArray * fTrackCandidates
 
TClonesArray * fMatches
 
std::map< Int_t, Int_t > fMatchMap
 
Int_t fVerbose
 
Bool_t fCollectionsComplete
 
std::vector< std::string > fHitCollectionNames
 
std::vector< std::string > fPointCollectionNames
 
TList fHitCollectionList
 
TList fPointCollectionList
 

Detailed Description

Definition at line 30 of file PndSttMatchTracks.h.

Constructor & Destructor Documentation

PndSttMatchTracks::PndSttMatchTracks ( )

Default constructor

Definition at line 22 of file PndSttMatchTracks.cxx.

References fCollectionsComplete, fMatches, fVerbose, and PndPersistencyTask::SetPersistency().

23  : PndPersistencyTask("STT track match") {
24  fMatches = NULL;
25  fVerbose = 1;
26  fCollectionsComplete = kFALSE;
27  SetPersistency(kTRUE);
28 }
void SetPersistency(Bool_t val=kTRUE)
TClonesArray * fMatches
PndSttMatchTracks::PndSttMatchTracks ( Int_t  verbose)

Constructor with verbosity level

Definition at line 34 of file PndSttMatchTracks.cxx.

References fCollectionsComplete, fMatches, fVerbose, PndPersistencyTask::SetPersistency(), and verbose.

35  : PndPersistencyTask("STT track match") {
36  fMatches = NULL;
37  fVerbose = verbose;
38  fCollectionsComplete = kFALSE;
39  SetPersistency(kTRUE);
40 }
#define verbose
void SetPersistency(Bool_t val=kTRUE)
TClonesArray * fMatches
PndSttMatchTracks::PndSttMatchTracks ( const char *  name,
const char *  title = "Pnd Stt Match Tracks Task",
Int_t  verbose = 1 
)

Constructor with name, title and verbosity

Parameters
nameName of taks
titleTitle of task (default PndPersistencyTask)
verboseVerbosity level (default 1)

Definition at line 46 of file PndSttMatchTracks.cxx.

References fCollectionsComplete, fMatches, fVerbose, PndPersistencyTask::SetPersistency(), SetTitle(), and verbose.

49  fMatches = NULL;
50  fVerbose = verbose;
51  fCollectionsComplete = kFALSE;
52  SetPersistency(kTRUE);
53  SetTitle(title);
54 }
#define verbose
void SetPersistency(Bool_t val=kTRUE)
h_MC_angle SetTitle("MC truth: opening angle of #pi^{0}")
TString name
TClonesArray * fMatches
PndSttMatchTracks::~PndSttMatchTracks ( )
virtual

Destructor

Definition at line 60 of file PndSttMatchTracks.cxx.

References fHitCollectionNames, and fPointCollectionNames.

61 {
62  fHitCollectionNames.clear();
63  fPointCollectionNames.clear();
64 }
std::vector< std::string > fHitCollectionNames
std::vector< std::string > fPointCollectionNames
PndSttMatchTracks::PndSttMatchTracks ( const PndSttMatchTracks L)
private

Member Function Documentation

void PndSttMatchTracks::AddAllCollections ( )
private

Definition at line 336 of file PndSttMatchTracks.cxx.

References AddHitCollection(), counter, fCollectionsComplete, fHitCollectionNames, and fPointCollectionNames.

Referenced by Exec().

337 {
339  {
340  for (size_t counter = 0; counter < fHitCollectionNames.size(); counter++)
341  {
343  }
344  fCollectionsComplete = kTRUE;
345  }
346 }
std::vector< std::string > fHitCollectionNames
std::vector< std::string > fPointCollectionNames
int counter
Definition: ZeeAnalysis.C:59
void AddHitCollection(char const *collectionName, char const *pointCollectionName)
void PndSttMatchTracks::AddHitCollection ( char const *  collectionName,
char const *  pointCollectionName 
)
private

Definition at line 300 of file PndSttMatchTracks.cxx.

References fHitCollectionList, and fPointCollectionList.

Referenced by AddAllCollections().

301 {
302  // Get and check FairRootManager
303  FairRootManager
304  *ioman = FairRootManager::Instance();
305 
306  if (!ioman)
307  {
308  cout << "-E- PndSttFindTracks::AddHitCollection: "
309  << "RootManager not instantised!" << endl;
310  }
311 
312  // Get hit Array
313  TClonesArray
314  *fHitArray = (TClonesArray*) ioman->GetObject(hitCollectionName);
315 
316  if (!fHitArray)
317  {
318  cout << "-W- PndSttFindTracks::AddHitCollection: No " << hitCollectionName << " STT hit array!"
319  << endl;
320  }
321 
322  // Get point Array
323  TClonesArray
324  *fPointArray = (TClonesArray*) ioman->GetObject(pointCollectionName);
325 
326  if (!fPointArray)
327  {
328  cout << "-W- PndSttFindTracks::AddHitCollection: No " << pointCollectionName << " STT hit array!"
329  << endl;
330  }
331 
332  fHitCollectionList.Add(fHitArray);
333  fPointCollectionList.Add(fPointArray);
334 }
void PndSttMatchTracks::AddHitCollectionName ( char *  hitCollectionName,
char *  pointCollectionName 
)

Add an hit collection to perform trackfinding on

Definition at line 290 of file PndSttMatchTracks.cxx.

References fHitCollectionNames, and fPointCollectionNames.

Referenced by idealcomplete(), and RecoComplete().

291 {
292  string
293  newPointName(pointCollectionName),
294  newHitName(hitCollectionName);
295 
296  fHitCollectionNames.push_back(newHitName);
297  fPointCollectionNames.push_back(newPointName);
298 }
std::vector< std::string > fHitCollectionNames
std::vector< std::string > fPointCollectionNames
PndSttMatchTracks::ClassDef ( PndSttMatchTracks  ,
 
)
private
void PndSttMatchTracks::Exec ( Option_t *  opt)
virtual

Execution

Definition at line 100 of file PndSttMatchTracks.cxx.

References AddAllCollections(), Double_t, fHitCollectionList, fMatches, fMatchMap, fPointCollectionList, fTrackCandidates, fVerbose, GetHitFromCollections(), PndTrackCandHit::GetHitId(), PndTrackCand::GetNHits(), GetPointFromCollections(), PndTrackCand::GetSortedHit(), nHits, and point.

101 {
103 
104  if (fHitCollectionList.GetEntries() == 0)
105  {
106  cout << "-E- PndSttTrackFinderIdeal::DoFind: "
107  << "No hit arrays present, call AddHitCollection() first (at least once)! " << endl;
108  }
109 
110  if (fPointCollectionList.GetEntries() == 0)
111  {
112  cout << "-E- PndSttTrackFinderIdeal::DoFind: "
113  << "No point arrays present, call AddHitCollection() first (at least once)! " << endl;
114  }
115 
116  // Clear output array
117  fMatches->Delete();
118 
119  // Create some pointers and variables
120  PndTrackCand* trackCand = NULL;
121  PndSttHit* mHit = NULL;
122  FairMCPoint* point = NULL;
123 
124  Int_t nHits = 0;
125  Int_t nMCTracks = 0;
126  Int_t iPoint = 0;
127  //Int_t iFlag = 0; //[R.K. 01/2017] unused variable?
128  Int_t iMCTrack = 0;
129  Int_t nAll = 0;
130  Int_t nTrue = 0;
131  Int_t nWrong = 0;
132  Int_t nFake = 0;
133  Int_t nHitSum = 0;
134  Int_t nTrueSum = 0;
135  Int_t nWrongSum = 0;
136  Int_t nFakeSum = 0;
137  Int_t nMCTrackSum = 0;
138 
139 
140  // Loop over TracksCand
141  Int_t nTracks = fTrackCandidates->GetEntriesFast();
142 
143 
144  for (Int_t iTrack=0; iTrack<nTracks; iTrack++)
145  {
146  trackCand = (PndTrackCand*) fTrackCandidates->At(iTrack);
147 
148  if ( ! trackCand)
149  {
150  cout << "-W- PndSttMatchTracks::Exec: Empty STTTrackCand at "
151  << iTrack << endl;
152  continue;
153  }
154 
155  nHits = trackCand->GetNHits();
156 
157  nAll = nTrue = nWrong = nFake = nMCTracks = 0;
158 
159  fMatchMap.clear();
160  if (fVerbose > 2) cout << endl << "Track " << iTrack << ", Hits "
161  << nHits << endl;
162 
163  // Loop over Hits of track
164  for (Int_t iMHit=0; iMHit<nHits; iMHit++)
165  {
166 
167  PndTrackCandHit candhit = trackCand->GetSortedHit(iMHit);
168 
169  // alter here
170  mHit = GetHitFromCollections(candhit.GetHitId());
171 
172 
173 
174  if ( ! mHit )
175  {
176  cout << "-E- PndSttMatchTracks::Exec: "
177  << "No Hit " << iMHit << " for track " << iTrack << endl;
178  continue;
179  }
180 
181 
182  iPoint = mHit->GetRefIndex();
183 
184  if ( iPoint < 0 )
185  {
186  nFake++;
187  continue;
188  }
189 
190  // alter here
191  point = GetPointFromCollections(candhit.GetHitId());
192 
193 
194  if ( ! point )
195  {
196  cout << "-E- PndSttMatchTracks::Exec: "
197  << "Empty MCPoint " << iPoint << " from Hit " << iMHit
198  << " (track " << iTrack << ")" << endl;
199  continue;
200  }
201 
202  iMCTrack = point->GetTrackID();
203 
204  if ( fVerbose > 2 ) cout << "Track " << iTrack << ", hit "
205  << candhit.GetHitId()
206  << ", STTPoint " << iPoint << ", MCTrack "
207  << iMCTrack << endl;
208  fMatchMap[iMCTrack]++;
209  }
210 
211  // Search for best matching MCTrack
212  iMCTrack = -1;
213 
214  map<Int_t, Int_t>::iterator it;
215  for (it=fMatchMap.begin(); it!=fMatchMap.end(); it++)
216  {
217  if (fVerbose > 2) cout << it->second
218  << " common points wth MCtrack "
219  << it->first << endl;
220  nMCTracks++;
221  nAll += it->second;
222  if ( it->second > nTrue )
223  {
224  iMCTrack = it->first;
225  nTrue = it->second;
226  }
227  }
228 
229  nWrong = nAll - nTrue;
230  if (fVerbose>1) cout << "-I- PndSttMatchTracks: STTTrack " << iTrack
231  << ", MCTrack " << iMCTrack << ", true "
232  << nTrue << ", wrong " << nWrong << ", fake "
233  << nFake << ", #MCTracks " << nMCTracks << endl;
234 
235  // Create SttTrackMatch
236  new ((*fMatches)[iTrack]) PndSttTrackMatch(iMCTrack, nTrue,
237  nWrong, nFake,
238  nMCTracks);
239 
240  // Some statistics
241  nHitSum += nHits;
242  nTrueSum += nTrue;
243  nWrongSum += nWrong;
244  nFakeSum += nFake;
245  nMCTrackSum += nMCTracks;
246 
247  } // Track loop
248 
249  // Event statistics
250 
251  Double_t qTrue = 0.;
252  if ( nHitSum)
253  qTrue = Double_t(nTrueSum) / Double_t(nHitSum) * 100.;
254 
255  if (fVerbose == 3)
256  {
257  Double_t
258  qWrong = 0.,
259  qFake = 0.,
260  qMC = 0.;
261 
262  if (nHitSum)
263  {
264  qWrong = Double_t(nWrongSum) / Double_t(nHitSum) * 100.;
265  qFake = Double_t(nFakeSum) / Double_t(nHitSum) * 100.;
266  }
267  if (nTracks)
268  qMC = Double_t(nMCTrackSum) / Double_t(nTracks);
269 
270  cout << endl;
271  cout << "-------------------------------------------------------"
272  << endl;
273  cout << "-I- STT Track Matching -I-"
274  << endl;
275  cout << "Reconstructed STTTracks : " << nTracks << endl;;
276  cout << "True hit assignments : " << qTrue << " %" << endl;
277  cout << "Wrong hit assignments : " << qWrong << " %" << endl;
278  cout << "Fake hit assignments : " << qFake << " %" << endl;
279  cout << "MCTracks per STTTrack : " << qMC << endl;
280  cout << "--------------------------------------------------------"
281  << endl;
282  }
283  else if (fVerbose) cout << "-I- PndSttMatchTracks: rec. " << nTracks << ", quota "
284  << qTrue << " % " << endl;
285 
286 }
PndTrackCandHit GetSortedHit(UInt_t i)
Definition: PndTrackCand.h:54
FairMCPoint * GetPointFromCollections(Int_t hitCounter)
TClonesArray * fTrackCandidates
PndSttHit * GetHitFromCollections(Int_t hitCounter)
int nHits
Definition: RiemannTest.C:16
Double_t
TClonesArray * point
Definition: anaLmdDigi.C:29
TClonesArray * fMatches
UInt_t GetNHits() const
Definition: PndTrackCand.h:59
std::map< Int_t, Int_t > fMatchMap
Int_t GetHitId() const
void PndSttMatchTracks::Finish ( )
virtual

Finishing

Definition at line 405 of file PndSttMatchTracks.cxx.

405 { }
PndSttHit * PndSttMatchTracks::GetHitFromCollections ( Int_t  hitCounter)

Definition at line 348 of file PndSttMatchTracks.cxx.

References At, fHitCollectionList, and GetEntriesFast().

Referenced by Exec().

349 {
350  PndSttHit
351  *retval = NULL;
352 
353  Int_t
354  relativeCounter = hitCounter;
355 
356  for (Int_t collectionCounter = 0; collectionCounter < fHitCollectionList.GetEntries(); collectionCounter++)
357  {
358  Int_t
359  size = ((TClonesArray *)fHitCollectionList.At(collectionCounter))->GetEntriesFast();
360 
361  if (relativeCounter < size)
362  {
363  retval = (PndSttHit*) ((TClonesArray *)fHitCollectionList.At(collectionCounter))->At(relativeCounter);
364  break;
365  }
366  else
367  {
368  relativeCounter -= size;
369  }
370  }
371  return retval;
372 }
cout<< "POINTs for new FwEndCap == "<< tsim-> GetEntriesFast()
cout<< "blue = Monte Carlo "<< endl;cout<< "red = Helix Hit "<< endl;cout<< "green = Center Of Tubes "<< endl;for(Int_t k=0;k< track->GetEntriesFast();k++){PndSttTrack *stttrack=(PndSttTrack *) track-> At(k)
Definition: checkhelixhit.C:56
Bool_t PndPersistencyTask::GetPersistency ( )
inlineinherited

Definition at line 32 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency.

Referenced by PndLmdPixelHitProducerFast::GetPersistance(), PndMdtDigitization::Init(), PndMdtHitProducerIdeal::Init(), PndMdtClusterTask::Init(), PndFtsHitProducerRealFast::Init(), PndDiscTaskReconstruction::Init(), PndRichHitProducer::Init(), PndSttHitProducerRealFast::Init(), PndSttHelixHitProducer::Init(), PndDiscTaskPID::Init(), PndIdealTrackFinder::Init(), PndSttMvdGemTracking::Init(), PndMdtTrkProducer::Init(), PndFtsHitProducerRealFull::Init(), PndLmdPixelClusterTask::Init(), PndSttHitProducerRealFull::Init(), PndLmdStripClusterTask::Init(), PndEmcApdHitProducer::Init(), PndMissingPzCleanerTask::Init(), PndEmcMakeRecoHit::Init(), PndEmcMakeClusterOnline::Init(), PndTrackSmearTask::Init(), PndEmcFWEndcapTimebasedWaveforms::Init(), PndSttHitProducerIdeal::Init(), PndEmcFWEndcapDigi::Init(), PndFtsHitProducerIdeal::Init(), PndEmcMakeCluster::Init(), PndMdtPointsToWaveform::Init(), PndDiscTaskDigitization::Init(), PndEmcMakeDigi::Init(), PndSdsTimeWalkCorrTask::Init(), PndLmdPixelHitProducerFast::Init(), PndDrcHitFinder::Init(), PndRichHitFinder::Init(), PndEmcMakeCorr::Init(), PndFtofHitProducerIdeal::Init(), PndEmcHitsToWaveform::Init(), PndSciTDigiTask::Init(), PndDrcHitProducerIdeal::Init(), PndSdsHitProducerIdeal::Init(), PndSciTHitProducerIdeal::Init(), PndRecoMultiKalmanTask2::Init(), PndEmcHitProducer::Init(), PndDrcHitProducerReal::Init(), PndDskFLGHitProducerIdeal::Init(), PndEmcTmpWaveformToDigi::Init(), PndDrcDigiTask::Init(), PndEmcWaveformToDigi::Init(), Init(), PndEmcWaveformToCalibratedDigi::Init(), PndTrkTracking2::Init(), PndSttFindTracks::Init(), PndEmcMultiWaveformToCalibratedDigi::Init(), PndRecoKalmanTask2::Init(), PndDrcTimeDigiTask::Init(), PndEmcExpClusterSplitter::Init(), PndFtsHoughTrackerTask::Init(), PndSdsNoiseProducer::Init(), PndEmcPhiBumpSplitter::Init(), PndSdsIdealRecoTask::Init(), PndSdsHybridHitProducer::Init(), PndRecoMultiKalmanTask::Init(), PndSdsIdealClusterTask::Init(), PndRecoKalmanTask::Init(), PndSdsStripHitProducerDif::Init(), PndGemDigitize::Init(), PndSdsStripHitProducer::Init(), PndGemFindHits::Init(), PndSdsPixelClusterTask::Init(), PndSdsStripClusterTask::Init(), PndMvdGemTrackFinderOnHits::Init(), PndBarrelTrackFinder::Init(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcMakeBump::PndEmcMakeBump(), PndUnassignedHitsTask::RegisterBranches(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndEmcMakeBump::SetStorageOfData(), and PndEmcFullDigiTask::StoreDigi().

32 { return fPersistency; }
FairMCPoint * PndSttMatchTracks::GetPointFromCollections ( Int_t  hitCounter)

Definition at line 374 of file PndSttMatchTracks.cxx.

References At, fHitCollectionList, fPointCollectionList, and GetEntriesFast().

Referenced by Exec().

375 {
376  FairMCPoint
377  *retval = NULL;
378 
379  Int_t
380  relativeCounter = hitCounter;
381 
382  for (Int_t collectionCounter = 0; collectionCounter < fHitCollectionList.GetEntries(); collectionCounter++)
383  {
384  Int_t
385  size = ((TClonesArray *)fHitCollectionList.At(collectionCounter))->GetEntriesFast();
386 
387  if (relativeCounter < size)
388  {
389  Int_t
390  tmpHit = ((PndSttHit*) ((TClonesArray *)fHitCollectionList.At(collectionCounter))->At(relativeCounter))->GetRefIndex();
391 
392  retval = (FairMCPoint*) ((TClonesArray *)fPointCollectionList.At(collectionCounter))->At(tmpHit);
393 
394  break;
395  }
396  else
397  {
398  relativeCounter -= size;
399  }
400  }
401  return retval;
402 }
cout<< "POINTs for new FwEndCap == "<< tsim-> GetEntriesFast()
cout<< "blue = Monte Carlo "<< endl;cout<< "red = Helix Hit "<< endl;cout<< "green = Center Of Tubes "<< endl;for(Int_t k=0;k< track->GetEntriesFast();k++){PndSttTrack *stttrack=(PndSttTrack *) track-> At(k)
Definition: checkhelixhit.C:56
InitStatus PndSttMatchTracks::Init ( )
virtual

Intialisation at beginning of each event

Definition at line 71 of file PndSttMatchTracks.cxx.

References fMatches, fTrackCandidates, and PndPersistencyTask::GetPersistency().

71  {
72 
73  // Get FairRootManager
74  FairRootManager* ioman = FairRootManager::Instance();
75  if (! ioman) {
76  cout << "-E- PndSttMatchTracks::Init: "
77  << "RootManager not instantised!" << endl;
78  return kFATAL;
79  }
80 
81  // Get TrackCand Array
82  fTrackCandidates = (TClonesArray*) ioman->GetObject("STTTrackCand");
83  if ( ! fTrackCandidates ) {
84  cout << "-E- PndSttMatchTracks::Init: No STTTrackCand array!" << endl;
85  return kERROR;
86  }
87 
88  // Create and register SttTrackMatch array
89  fMatches = new TClonesArray("PndSttTrackMatch",100);
90  ioman->Register("STTTrackMatch", "STT", fMatches, GetPersistency());
91 
92  return kSUCCESS;
93 
94 }
TClonesArray * fTrackCandidates
TClonesArray * fMatches
PndSttMatchTracks& PndSttMatchTracks::operator= ( const PndSttMatchTracks )
inlineprivate

Definition at line 97 of file PndSttMatchTracks.h.

97 {return *this;};
void PndSttMatchTracks::SetPersistence ( Bool_t  persistence)
inline

set persistence flag

Definition at line 74 of file PndSttMatchTracks.h.

References PndPersistencyTask::SetPersistency().

74 { SetPersistency(persistence); }
void SetPersistency(Bool_t val=kTRUE)
void PndPersistencyTask::SetPersistency ( Bool_t  val = kTRUE)
inlineinherited

Definition at line 31 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency, and val.

Referenced by barrelTrackFinder(), digi_complete(), digi_complete_newSTT(), digiOnly_complete(), PndBarrelTrackFinder::PndBarrelTrackFinder(), PndCATracking::PndCATracking(), PndDrcHitFinder::PndDrcHitFinder(), PndEmc2DLocMaxFinder::PndEmc2DLocMaxFinder(), PndEmcExpClusterSplitter::PndEmcExpClusterSplitter(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcFWEndcapDigi::PndEmcFWEndcapDigi(), PndEmcFWEndcapTimebasedWaveforms::PndEmcFWEndcapTimebasedWaveforms(), PndEmcHitProducer::PndEmcHitProducer(), PndEmcHitsToWaveform::PndEmcHitsToWaveform(), PndEmcMakeBump::PndEmcMakeBump(), PndEmcMakeCluster::PndEmcMakeCluster(), PndEmcMakeClusterOnline::PndEmcMakeClusterOnline(), PndEmcMakeDigi::PndEmcMakeDigi(), PndEmcMakeRecoHit::PndEmcMakeRecoHit(), PndEmcMultiWaveformToCalibratedDigi::PndEmcMultiWaveformToCalibratedDigi(), PndEmcPhiBumpSplitter::PndEmcPhiBumpSplitter(), PndEmcTmpWaveformToDigi::PndEmcTmpWaveformToDigi(), PndEmcWaveformToCalibratedDigi::PndEmcWaveformToCalibratedDigi(), PndEmcWaveformToDigi::PndEmcWaveformToDigi(), PndFtofHitProducerIdeal::PndFtofHitProducerIdeal(), PndFtsCATracking::PndFtsCATracking(), PndFtsHitProducerIdeal::PndFtsHitProducerIdeal(), PndFtsHitProducerRealFast::PndFtsHitProducerRealFast(), PndFtsHitProducerRealFull::PndFtsHitProducerRealFull(), PndFtsHoughTrackerTask::PndFtsHoughTrackerTask(), PndGemDigitize::PndGemDigitize(), PndGemFindHits::PndGemFindHits(), PndIdealTrackFinder::PndIdealTrackFinder(), PndLmdPixelClusterTask::PndLmdPixelClusterTask(), PndLmdPixelHitProducerFast::PndLmdPixelHitProducerFast(), PndMdtClusterTask::PndMdtClusterTask(), PndMdtDigitization::PndMdtDigitization(), PndMdtHitProducerIdeal::PndMdtHitProducerIdeal(), PndMdtPointsToWaveform::PndMdtPointsToWaveform(), PndMdtTrkProducer::PndMdtTrkProducer(), PndMissingPzCleanerTask::PndMissingPzCleanerTask(), PndMvdGemTrackFinderOnHits::PndMvdGemTrackFinderOnHits(), PndMvdHitProducerIdeal::PndMvdHitProducerIdeal(), PndMvdPixelClusterTask::PndMvdPixelClusterTask(), PndMvdTimeWalkCorrTask::PndMvdTimeWalkCorrTask(), PndMvdToPix4ClusterTask::PndMvdToPix4ClusterTask(), PndRecoKalmanTask::PndRecoKalmanTask(), PndRecoKalmanTask2::PndRecoKalmanTask2(), PndRecoMultiKalmanTask::PndRecoMultiKalmanTask(), PndRecoMultiKalmanTask2::PndRecoMultiKalmanTask2(), PndRichHitFinder::PndRichHitFinder(), PndRichHitProducer::PndRichHitProducer(), PndSciTDigiTask::PndSciTDigiTask(), PndSciTHitProducerIdeal::PndSciTHitProducerIdeal(), PndSdsHitProducerIdeal::PndSdsHitProducerIdeal(), PndSdsHybridHitProducer::PndSdsHybridHitProducer(), PndSdsIdealClusterTask::PndSdsIdealClusterTask(), PndSdsIdealRecoTask::PndSdsIdealRecoTask(), PndSdsNoiseProducer::PndSdsNoiseProducer(), PndSdsPixelClusterTask::PndSdsPixelClusterTask(), PndSdsStripClusterTask::PndSdsStripClusterTask(), PndSdsStripHitProducer::PndSdsStripHitProducer(), PndSdsTimeWalkCorrTask::PndSdsTimeWalkCorrTask(), PndSttFindTracks::PndSttFindTracks(), PndSttHelixHitProducer::PndSttHelixHitProducer(), PndSttHitProducerIdeal::PndSttHitProducerIdeal(), PndSttHitProducerRealFast::PndSttHitProducerRealFast(), PndSttHitProducerRealFull::PndSttHitProducerRealFull(), PndSttMatchTracks(), PndSttMvdGemTracking::PndSttMvdGemTracking(), PndTrackSmearTask::PndTrackSmearTask(), PndTrkTracking2::PndTrkTracking2(), reco(), reco_complete(), reco_complete_gf2(), reco_complete_newSTT(), reco_complete_sec(), recoideal_complete(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndLmdPixelHitProducerFast::SetPersistance(), PndSdsHitProducerIdeal::SetPersistance(), PndSttMvdGemTracking::SetPersistenc(), PndMdtClusterTask::SetPersistence(), PndSttHelixHitProducer::SetPersistence(), PndMissingPzCleanerTask::SetPersistence(), PndFtsHitProducerRealFast::SetPersistence(), PndFtsHitProducerRealFull::SetPersistence(), PndSttHitProducerIdeal::SetPersistence(), PndSttHitProducerRealFull::SetPersistence(), PndSttHitProducerRealFast::SetPersistence(), PndFtsHitProducerIdeal::SetPersistence(), PndTrackSmearTask::SetPersistence(), PndSciTHitProducerIdeal::SetPersistence(), PndIdealTrackFinder::SetPersistence(), SetPersistence(), PndSttFindTracks::SetPersistence(), PndFtsHoughTrackerTask::SetPersistence(), PndTrkTracking2::SetPersistence(), PndEmcMakeRecoHit::SetStorageOfData(), PndEmcMakeClusterOnline::SetStorageOfData(), PndEmcFWEndcapDigi::SetStorageOfData(), PndEmcFWEndcapTimebasedWaveforms::SetStorageOfData(), PndEmcMakeDigi::SetStorageOfData(), PndMdtPointsToWaveform::SetStorageOfData(), PndEmc2DLocMaxFinder::SetStorageOfData(), PndEmcMakeCluster::SetStorageOfData(), PndEmcHitsToWaveform::SetStorageOfData(), PndEmcMakeBump::SetStorageOfData(), PndEmcTmpWaveformToDigi::SetStorageOfData(), PndEmcWaveformToDigi::SetStorageOfData(), PndEmcWaveformToCalibratedDigi::SetStorageOfData(), PndEmcMultiWaveformToCalibratedDigi::SetStorageOfData(), PndEmcExpClusterSplitter::SetStorageOfData(), PndEmcPhiBumpSplitter::SetStorageOfData(), standard_tracking(), and PndEmcFullDigiTask::StoreDigi().

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

Member Data Documentation

Bool_t PndSttMatchTracks::fCollectionsComplete
private

Definition at line 89 of file PndSttMatchTracks.h.

Referenced by AddAllCollections(), and PndSttMatchTracks().

TList PndSttMatchTracks::fHitCollectionList
private
std::vector<std::string> PndSttMatchTracks::fHitCollectionNames
private

Definition at line 91 of file PndSttMatchTracks.h.

Referenced by AddAllCollections(), AddHitCollectionName(), and ~PndSttMatchTracks().

TClonesArray* PndSttMatchTracks::fMatches
private

Definition at line 81 of file PndSttMatchTracks.h.

Referenced by Exec(), Init(), and PndSttMatchTracks().

std::map<Int_t, Int_t> PndSttMatchTracks::fMatchMap
private

Map from MCTrackID to number of common hits

Definition at line 84 of file PndSttMatchTracks.h.

Referenced by Exec().

TList PndSttMatchTracks::fPointCollectionList
private

Definition at line 94 of file PndSttMatchTracks.h.

Referenced by AddHitCollection(), Exec(), and GetPointFromCollections().

std::vector<std::string> PndSttMatchTracks::fPointCollectionNames
private

Definition at line 92 of file PndSttMatchTracks.h.

Referenced by AddAllCollections(), AddHitCollectionName(), and ~PndSttMatchTracks().

TClonesArray* PndSttMatchTracks::fTrackCandidates
private

Definition at line 80 of file PndSttMatchTracks.h.

Referenced by Exec(), and Init().

Int_t PndSttMatchTracks::fVerbose
private

Verbosity level

Definition at line 87 of file PndSttMatchTracks.h.

Referenced by Exec(), and PndSttMatchTracks().


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