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

#include <PndMvdEventAnaTask.h>

Inheritance diagram for PndMvdEventAnaTask:

Public Member Functions

 PndMvdEventAnaTask ()
 
virtual ~PndMvdEventAnaTask ()
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void DrawTracksPerEvent (TString opt="")
 
void DrawHitsPerTrack (TString opt="")
 
void DrawEnergyPerHit (TString opt="")
 
void DrawPointRes (TString opt="")
 
void DrawPointResS (TString opt="")
 
void DrawPointResD (TString opt="")
 
void DrawPointResM (TString opt="")
 
void DrawEnergyRes (TString opt="")
 
void DrawDigisPerCluster (TString opt="")
 
void DrawPointResStrip (TString opt="")
 
void DrawEnergyResStrip (TString opt="")
 
void DrawDigisPerClusterStrip (TString opt="")
 
void DrawPtRes (TString opt="")
 
void DrawPRes (TString opt="")
 
void DrawRiemannRes (TString opt="")
 
void DrawRiemannFakes (TString opt="")
 
void DrawRiemannTracksPerTrack (TString opt="")
 
void DrawRiemannTracksPerTrackAdd (TString opt="")
 
void DrawRiemannVertexResolutionX (TString opt="")
 
void DrawRiemannVertexResolutionY (TString opt="")
 
void DrawRiemannVertexResolutionZ (TString opt="")
 

Private Member Functions

bool MCHitBelongsToCluster (int HitIndex, PndSdsCluster *cluster, bool pixCluster)
 
void GetTrackCandsForMCTrack (std::vector< int > pixHitId, std::vector< int > stripHitId, std::vector< int > &matches, std::vector< int > &result)
 
int GetRecoHit (int clIndex, bool pixel) const
 
std::vector< int > GetClusters (int MCHit, bool pixel)
 
void PrintClusterDigiInfo (int clIndex, std::vector< Int_t > digiInd, bool pixel)
 
void PrintRecoHitInfo (int hitInd, int digiSize, TVector3 MCPos, double MCEnergy, bool pixel) const
 
void PrintTrackCand (GFTrackCand *cand) const
 
void Register ()
 
void Reset ()
 
void ProduceHits ()
 
FairHit * GetFairHit (Int_t detId, Int_t hitId)
 
std::map< int, std::vector< int > > AssignHitsToTracks ()
 
 ClassDef (PndMvdEventAnaTask, 1)
 

Private Attributes

TClonesArray * fMCTracks
 
TClonesArray * fMCHits
 
TClonesArray * fPixDigis
 
TClonesArray * fStripDigis
 
TClonesArray * fPixReco
 
TClonesArray * fStripReco
 
TClonesArray * fPixCluster
 
TClonesArray * fStripCluster
 
TClonesArray * fTrackCand
 
std::vector< int > fGhostCand
 
TH1 * fHTracksPerEvent
 
TH1 * fHHitsPerTrack
 
TH1 * fHEnergyPerHit
 
TH1 * fHPointRes
 
TH1 * fHPointResS
 
TH1 * fHPointResD
 
TH1 * fHPointResM
 
TH1 * fHDigisPerCluster
 
TH1 * fHEnergyRes
 
TH1 * fHPRes
 
TH1 * fHPtRes
 
TH1 * fHPointResStrip
 
TH1 * fHPointResSStrip
 
TH1 * fHPointResDStrip
 
TH1 * fHPointResMStrip
 
TH1 * fHDigisPerClusterStrip
 
TH1 * fHEnergyResStrip
 
TH1 * fHRiemannRes
 
TH1 * fHRiemannFakes
 
TH1 * fHRiemannTracksPerTrack
 
TH1 * fHRiemannTracksPerTrackAdd
 
TH1 * fHRiemannVertexResolutionX
 
TH1 * fHRiemannVertexResolutionY
 
TH1 * fHRiemannVertexResolutionZ
 
bool fPrintTrack
 
bool fPrintMCHit
 
bool fPrintCluster
 
bool fPrintPixDigis
 
bool fPrintPixHit
 
bool fPrintStripCluster
 
bool fPrintStripDigis
 
bool fPrintStripHit
 
bool fPrintTrackMatch
 
bool fPrintGhosts
 
int fNTracks
 
int fNPossibleTracks
 
int fNCompleteTracks
 
int fNPartTracks
 
int fNNotFoundPossibleTracks
 
int fNNotFoundTracks
 
int fNGhostTracks
 
int fEventNr
 
std::map< int, std::vector< int > > fTrackPixHitIdMap
 
std::map< int, std::vector< int > > fTrackStripHitIdMap
 

Detailed Description

Definition at line 30 of file PndMvdEventAnaTask.h.

Constructor & Destructor Documentation

PndMvdEventAnaTask::PndMvdEventAnaTask ( )

Default constructor

Definition at line 35 of file PndMvdEventAnaTask.cxx.

35  : FairTask("Event Display task for PANDA PndMvd"),
36  fPrintTrack(true), fPrintMCHit(true), fPrintCluster(true), fPrintPixDigis(true), fPrintPixHit(true),
39  fNGhostTracks(0), fEventNr(0)
40 {
41 }
PndMvdEventAnaTask::~PndMvdEventAnaTask ( )
virtual

Destructor

Definition at line 46 of file PndMvdEventAnaTask.cxx.

47 {
48 }

Member Function Documentation

std::map< int, std::vector< int > > PndMvdEventAnaTask::AssignHitsToTracks ( )
private

Definition at line 640 of file PndMvdEventAnaTask.cxx.

References fMCHits, and i.

Referenced by Exec().

641 {
642  std::map<int, std::vector<int> > result;
643  for (int i = 0; i < fMCHits->GetEntriesFast(); i++){ //get all MC Hits
644  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMCHits->At(i)); //sort MCHits with Tracks
645  //PndMCTrack* myTrack = (PndMCTrack*)(fMCTracks->At(myPoint->GetTrackID()));
646  result[myPoint->GetTrackID()].push_back(i);
647 
648  }
649  return result;
650 }
Int_t i
Definition: run_full.C:25
TClonesArray * fMCHits
PndMvdEventAnaTask::ClassDef ( PndMvdEventAnaTask  ,
 
)
private
void PndMvdEventAnaTask::DrawDigisPerCluster ( TString  opt = "")
inline

Definition at line 61 of file PndMvdEventAnaTask.h.

References fHDigisPerCluster.

61 {fHDigisPerCluster->DrawClone(opt);}
void PndMvdEventAnaTask::DrawDigisPerClusterStrip ( TString  opt = "")
inline

Definition at line 65 of file PndMvdEventAnaTask.h.

References fHDigisPerClusterStrip.

65 {fHDigisPerClusterStrip->DrawClone(opt);}
void PndMvdEventAnaTask::DrawEnergyPerHit ( TString  opt = "")
inline

Definition at line 54 of file PndMvdEventAnaTask.h.

References fHEnergyPerHit.

54 {fHEnergyPerHit->DrawClone(opt);}
void PndMvdEventAnaTask::DrawEnergyRes ( TString  opt = "")
inline

Definition at line 60 of file PndMvdEventAnaTask.h.

References fHEnergyRes.

60 {fHEnergyRes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawEnergyResStrip ( TString  opt = "")
inline

Definition at line 64 of file PndMvdEventAnaTask.h.

References fHEnergyResStrip.

64 {fHEnergyResStrip->DrawClone(opt);}
void PndMvdEventAnaTask::DrawHitsPerTrack ( TString  opt = "")
inline

Definition at line 53 of file PndMvdEventAnaTask.h.

References fHHitsPerTrack.

53 {fHHitsPerTrack->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPointRes ( TString  opt = "")
inline

Definition at line 56 of file PndMvdEventAnaTask.h.

References fHPointRes.

56 {fHPointRes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPointResD ( TString  opt = "")
inline

Definition at line 58 of file PndMvdEventAnaTask.h.

References fHPointResD.

58 {fHPointResD->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPointResM ( TString  opt = "")
inline

Definition at line 59 of file PndMvdEventAnaTask.h.

References fHPointResM.

59 {fHPointResM->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPointResS ( TString  opt = "")
inline

Definition at line 57 of file PndMvdEventAnaTask.h.

References fHPointResS.

57 {fHPointResS->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPointResStrip ( TString  opt = "")
inline

Definition at line 63 of file PndMvdEventAnaTask.h.

References fHPointResStrip.

63 {fHPointResStrip->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPRes ( TString  opt = "")
inline

Definition at line 68 of file PndMvdEventAnaTask.h.

References fHPRes.

68 {fHPRes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawPtRes ( TString  opt = "")
inline

Definition at line 67 of file PndMvdEventAnaTask.h.

References fHPtRes.

67 {fHPtRes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannFakes ( TString  opt = "")
inline

Definition at line 71 of file PndMvdEventAnaTask.h.

References fHRiemannFakes.

71 {fHRiemannFakes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannRes ( TString  opt = "")
inline

Definition at line 70 of file PndMvdEventAnaTask.h.

References fHRiemannRes.

70 {fHRiemannRes->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannTracksPerTrack ( TString  opt = "")
inline

Definition at line 72 of file PndMvdEventAnaTask.h.

References fHRiemannTracksPerTrack.

72 {fHRiemannTracksPerTrack->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannTracksPerTrackAdd ( TString  opt = "")
inline

Definition at line 73 of file PndMvdEventAnaTask.h.

References fHRiemannTracksPerTrackAdd.

73 {fHRiemannTracksPerTrackAdd->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannVertexResolutionX ( TString  opt = "")
inline

Definition at line 75 of file PndMvdEventAnaTask.h.

References fHRiemannVertexResolutionX.

75 {fHRiemannVertexResolutionX->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannVertexResolutionY ( TString  opt = "")
inline

Definition at line 76 of file PndMvdEventAnaTask.h.

References fHRiemannVertexResolutionY.

76 {fHRiemannVertexResolutionY->DrawClone(opt);}
void PndMvdEventAnaTask::DrawRiemannVertexResolutionZ ( TString  opt = "")
inline

Definition at line 77 of file PndMvdEventAnaTask.h.

References fHRiemannVertexResolutionZ.

77 {fHRiemannVertexResolutionZ->DrawClone(opt);}
void PndMvdEventAnaTask::DrawTracksPerEvent ( TString  opt = "")
inline

Definition at line 52 of file PndMvdEventAnaTask.h.

References fHTracksPerEvent.

52 {fHTracksPerEvent->DrawClone(opt);}
void PndMvdEventAnaTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 168 of file PndMvdEventAnaTask.cxx.

References PndRiemannTrack::addHit(), AssignHitsToTracks(), fEventNr, fGhostCand, fHDigisPerCluster, fHDigisPerClusterStrip, fHEnergyPerHit, fHHitsPerTrack, fHPRes, fHPtRes, fHRiemannFakes, fHRiemannRes, fHRiemannTracksPerTrack, fHRiemannTracksPerTrackAdd, fHRiemannVertexResolutionX, fHRiemannVertexResolutionY, fHRiemannVertexResolutionZ, fHTracksPerEvent, fMCHits, fMCTracks, fNCompleteTracks, fNGhostTracks, fNNotFoundPossibleTracks, fNNotFoundTracks, fNPartTracks, fNPossibleTracks, fNTracks, fPixCluster, fPrintGhosts, fPrintMCHit, fPrintTrack, fPrintTrackMatch, fStripCluster, fTrackCand, fTrackPixHitIdMap, fTrackStripHitIdMap, PndSdsCluster::GetClusterList(), GetClusters(), GFTrackCand::getCurv(), GFTrackCand::getDip(), GetFairHit(), GFTrackCand::getHit(), PndMCTrack::GetMomentum(), PndMCTrack::GetMotherID(), GFTrackCand::getNHits(), PndSdsMCPoint::GetPosition(), PndSdsMCPoint::GetPositionOut(), GetRecoHit(), PndMCTrack::GetStartVertex(), GetTrackCandsForMCTrack(), hit(), i, map, n, PndRiemannTrack::n(), p, p1, p2, PndMCTrack::Print(), PndSdsMCPoint::Print(), PrintClusterDigiInfo(), PrintRecoHitInfo(), PrintTrackCand(), pt(), PndRiemannTrack::refit(), PndRiemannTrack::SetVerbose(), and PndRiemannTrack::szFit().

169 {
170  std::map<int, std::vector<int> > mcHitMap; //Track -> MCHits
171  std::map<int, std::vector<int> > trackToTrackCandMap; //Track -> TrackCand
172 
173  std::vector<PndRiemannTrack> riemannTracks;
174  std::vector<int> MCTrackOrderRiemann; //information which riemannTrack belongs to which MCTrack;
175 
176  TVector3 MCPos, RecoPos;
177  double MCEnergy;
178  double TrackP, TrackPt;
179 
180  mcHitMap = AssignHitsToTracks();
181  fGhostCand.clear();
182  fGhostCand.resize(fTrackCand->GetEntriesFast(), 0);
183  fTrackPixHitIdMap.clear();
184  fTrackStripHitIdMap.clear();
185 
186  std::cout << "------------Event " << fEventNr << "-------------" << std::endl;
187  fEventNr++;
188  fHTracksPerEvent->Fill(mcHitMap.size());
189  fNTracks += mcHitMap.size();
190 
191  for (std::map<int, std::vector<int> >::const_iterator kIt = mcHitMap.begin(); kIt!= mcHitMap.end(); kIt++){ //go through all tracks
192  PndMCTrack* myTrack = (PndMCTrack*)(fMCTracks->At(kIt->first));
193  std::vector<int> MChits = kIt->second;
194  TrackPt = myTrack->GetMomentum().Pt();
195  TrackP = myTrack->GetMomentum().Mag();
196  fHHitsPerTrack->Fill(MChits.size());
197 
198  if (fPrintTrack){
199  std::cout << "<<<<<<<<<<< MCTrack >>>>>>>>>> " << std::endl;
200  myTrack->Print(kIt->first);
201  TVector3 startVertex = myTrack->GetStartVertex();
202  std::cout << "Pt: " << TrackPt << " GeV/c; P: " << TrackP << " GeV/c" << std::endl;
203  std::cout << "StartVertex: " << startVertex.X() << " " << startVertex.Y() << " " << startVertex.Z() << std::endl;
204  }
205  for (unsigned int p = 0; p < MChits.size(); p++){ //go through all hits in track
206  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(fMCHits->At(MChits[p]));
207  if (fPrintMCHit){
208  std::cout << "-------------------------------" << std::endl;
209  myPoint->Print();
210  }
211  MCPos = (myPoint->GetPosition() + myPoint->GetPositionOut())*0.5;
212  MCEnergy = myPoint->GetEnergyLoss()*10E9;
213  fHEnergyPerHit->Fill(MCEnergy);
214 
215  std::vector<int> pixCluster = GetClusters(MChits[p], true); //Indices of all clusters belonging to one hit
216 
217  for (unsigned int clInd = 0; clInd < pixCluster.size(); clInd++) {
218 
219  PndSdsCluster* myCluster = (PndSdsCluster*)(fPixCluster->At(pixCluster[clInd]));
220  std::vector<Int_t> digiInd = myCluster->GetClusterList(); //gets the list of pixel Digis belonging to the cluster
221  int recoHit = GetRecoHit(pixCluster[clInd], true);
222  fHDigisPerCluster->Fill(digiInd.size());
223  if (recoHit > -1)
224  fTrackPixHitIdMap[kIt->first].push_back(recoHit);
225 
226  PrintClusterDigiInfo(pixCluster[clInd], digiInd, true);
227  if (recoHit < 0)
228  std::cout << "-W- No Reco Hit found!" << std::endl;
229  else
230  PrintRecoHitInfo(recoHit, digiInd.size(), MCPos, MCEnergy, true);
231  }
232 
233  std::vector<int> stripCluster = GetClusters(MChits[p], false);
234 
235  int hitCount = 0;
236  int oldRecoHit = -1;
237  for (unsigned int clInd = 0; clInd < stripCluster.size(); clInd++) {
238  int recoHit = -1;
239 
240  PndSdsCluster* myCluster = (PndSdsCluster*)(fStripCluster->At(stripCluster[clInd]));
241  std::vector<Int_t> digiInd = myCluster->GetClusterList(); //gets the list of pixel Digis belonging to the cluster
242  recoHit = GetRecoHit(stripCluster[clInd], false);
243 
244  fHDigisPerClusterStrip->Fill(digiInd.size());
245  std::cout << "RecoHit: " << recoHit << std::endl;
246  if (recoHit > -1){
247  hitCount++;
248  if (oldRecoHit != recoHit){
249  fTrackStripHitIdMap[kIt->first].push_back(recoHit);
250  oldRecoHit = recoHit;
251  }
252  }
253  PrintClusterDigiInfo(stripCluster[clInd], digiInd, false);
254  //std::cout << " hitCount: " << hitCount << " oldReco " << oldRecoHit << " recoHit: " << recoHit << std::endl;
255  if (recoHit > -1){
256  if (hitCount > 1){
257  PrintRecoHitInfo(recoHit, digiInd.size(), MCPos, MCEnergy, false);
258  if ((oldRecoHit > -1) && (oldRecoHit != recoHit))
259  PrintRecoHitInfo(oldRecoHit, digiInd.size(), MCPos, MCEnergy, false);
260  }
261  }
262  }
263  }
264  std::vector<int> pixHits = fTrackPixHitIdMap[kIt->first];
265  std::vector<int> stripHits = fTrackStripHitIdMap[kIt->first];
266 
267  std::cout << std::endl;
268  std::cout << "TrackID " << kIt->first << ": ";
269  for (unsigned int testInd = 0; testInd < pixHits.size(); testInd++)
270  std::cout << " 5/" << pixHits[testInd];
271 
272  for (unsigned int testInd = 0; testInd < stripHits.size(); testInd++)
273  std::cout << " 4/" << stripHits[testInd];
274  std::cout << std::endl;
275 
276  std::cout << "MCHitMap: ";
277  for (unsigned int testInd = 0; testInd < mcHitMap[kIt->first].size(); testInd++)
278  std::cout << mcHitMap[kIt->first].at(testInd) << " ";
279  std::cout << std::endl;
280 
281  std::vector<int> matches;
282  std::vector<int> candidates;
283  GetTrackCandsForMCTrack(pixHits, stripHits, matches, candidates);
284  std::cout << "TrackCands for MCTrack: ";
285  for (unsigned int i = 0; i < candidates.size(); i++) std::cout << candidates[i];
286  std::cout << std::endl;
287  trackToTrackCandMap[kIt->first] = candidates;
288 
289  int TrackMatch = 0;
290  int oldmatches = 0;
291  int highestMatch=-1;
292  for (unsigned int i = 0; i < matches.size(); i++){
293  if (oldmatches < matches[i]){
294  oldmatches = matches[i];
295  highestMatch = i;
296  }
297  GFTrackCand* myCand = (GFTrackCand*)fTrackCand->At(candidates[i]);
298  PndRiemannTrack myRiemannTrack;
299  MCTrackOrderRiemann.push_back(kIt->first);
300  myRiemannTrack.SetVerbose(1);
301  for (UInt_t j = 0; j < myCand->getNHits(); j++){
302  unsigned int detId, hitId;
303  myCand->getHit(j, detId, hitId);
304  PndRiemannHit hit(GetFairHit(detId, hitId));
305  myRiemannTrack.addHit(hit);
306  }
307 
308  myRiemannTrack.refit();
309  myRiemannTrack.szFit();
310  TVectorD n = myRiemannTrack.n();
311  //std::cout << "RiemannNormal: " << n[0] << " " << n[1] << " " << n[2] << std::endl;
312  riemannTracks.push_back(myRiemannTrack);
313 
314  double curv = myCand->getCurv();
315  double dip = myCand->getDip();
316  double pt = 1/curv;
317  pt*=2;
318  pt*=0.00299792;
319  unsigned int detId, hitId;
320  if (matches[i]>2){
321  TrackMatch++;
322  myCand->getHit(0,detId, hitId);
323  if (detId > 0)
324  fHRiemannFakes->Fill(myCand->getNHits()- matches[i]);
325  else fHRiemannFakes->Fill(myCand->getNHits()- matches[i] -1);
326  }
327  if (fPrintTrackMatch){
328  myCand->getHit(0,detId, hitId);
329  std::cout << std::endl;
330  std::cout << "HitMatch: " << matches[i] << std::endl;
331  std::cout << "Hits in Track: " << myCand->getNHits() << " Difference: ";
332  if (detId > 0) std::cout << myCand->getNHits() - matches[i] << std::endl;
333  else std::cout << myCand->getNHits() - matches [i] -1 << std::endl;
334 
335  PrintTrackCand(myCand);
336  std::cout << "Pt: " << pt << " GeV/c " << " P: " << pt/dip << " GeV/c" << std::endl;
337  std::cout << "Pt error: " << TrackPt - pt << " P error: " << TrackP - pt/dip << std::endl;
338  std::cout << "TrackID: " << kIt->first << " mcHitMap.size: " << mcHitMap[kIt->first].size() << " matches: " << matches[i] << std::endl;
339  std::cout << "Missing Hits: " << mcHitMap[kIt->first].size() - matches[i] << std::endl;
340  if (mcHitMap[kIt->first].size() < (myCand->getNHits()-1))
341  fHRiemannRes->Fill(-2); //more hits than expected
342  }
343  }
344 
345  if (mcHitMap[kIt->first].size() > 2){
346  if (mcHitMap[kIt->first].size() > 3){
348  if (oldmatches == 0)
350  }
351  else if (myTrack->GetMotherID() < 0){
353  if (oldmatches == 0)
355  }
356  }
357  if (oldmatches == 0)
359 
360 
361  if (matches.size() > 0){ // there are riemann tracks to a MC Track
362  GFTrackCand* myCand = (GFTrackCand*)fTrackCand->At(candidates[highestMatch]);
363  double curv = myCand->getCurv();
364  double dip = myCand->getDip();
365  double pt = 1/curv;
366  pt*=2;
367  pt*=0.00299792;
368  fHPtRes->Fill(TrackPt - pt);
369  fHPRes->Fill(TrackP - pt/dip);
370  fHRiemannRes->Fill(mcHitMap[kIt->first].size() - matches[highestMatch]);
371  if (mcHitMap[kIt->first].size() - matches[highestMatch] == 0)
373  else fNPartTracks++;
374  }
375  std::cout << "Found TracksPerTrack: " << TrackMatch;
376  if (mcHitMap[kIt->first].size() > 2) std::cout << " 3Hits+" << std::endl;
377  if (mcHitMap[kIt->first].size() > 2) fHRiemannTracksPerTrackAdd->Fill(TrackMatch);
378  fHRiemannTracksPerTrack->Fill(TrackMatch);
379  }
380  if (riemannTracks.size() > 1){
381  std::cout << riemannTracks.size() << " Riemann Tracks used!" << std::endl;
382  for (unsigned int i = 0; i < riemannTracks.size() - 1; i++){
383  for (unsigned int j = i + 1; j < riemannTracks.size(); j++) {
384  TVector3 p1, p2;
385  //riemannTracks[i].refit();
386  //riemannTracks[i].szFit();
387  //std::cout << "Track: " << i << std::endl;
388  //riemannTracks[i].PrintHits();
389 
390  //riemannTracks[j].refit();
391  //riemannTracks[j].szFit();
392  //std::cout << "Track: " << j << std::endl;
393  //riemannTracks[j].PrintHits();
394 
395  std::cout << "Vertex Test for: " << MCTrackOrderRiemann[i] << " and " << MCTrackOrderRiemann[j] << std::endl;
396  int result = riemannTracks[i].calcIntersection(riemannTracks[j], p1, p2);
397 
398  if (result == 1){
399  std::cout << "Vertex1: " << p1.X() << " " << p1.Y() << " " << p1.Z() << std::endl;
400  std::cout << "Vertex2: " << p2.X() << " " << p2.Y() << " " << p2.Z() << std::endl;
401  PndMCTrack* myMCTrack = (PndMCTrack*)(fMCTracks->At(MCTrackOrderRiemann[i]));
402  TVector3 MCVertex = myMCTrack->GetStartVertex();
403  std::cout << "MCVertex " << MCTrackOrderRiemann[i] << ": "<< MCVertex.X() << " " << MCVertex.Y() << " " << MCVertex.Z() << std::endl;
404  MCVertex -= (p1+p2)*0.5;
405  std::cout << "Difference: " << MCVertex.X() << " " << MCVertex.Y() << " " << MCVertex.Z() << std::endl;
406  fHRiemannVertexResolutionX->Fill(MCVertex.X());
407  fHRiemannVertexResolutionY->Fill(MCVertex.Y());
408  fHRiemannVertexResolutionZ->Fill(MCVertex.Z());
409 
410 
411  }
412  }
413  }
414  }
415  if (fPrintGhosts){
416  std::cout << std::endl;
417  std::cout << "----Ghost TrackCandidates:-------" << std::endl;
418  for (unsigned int tcIndex = 0; tcIndex < fGhostCand.size(); tcIndex++){
419  if (fGhostCand[tcIndex] < 3){
420  GFTrackCand* myTrackCand = (GFTrackCand*)fTrackCand->At(tcIndex);
421  PrintTrackCand(myTrackCand);
422  fHRiemannRes->Fill(-1); //Ghosts
423  fNGhostTracks++;
424  }
425  }
426  std::cout << std::endl;
427  }
428 }
std::vector< Int_t > GetClusterList() const
Definition: PndSdsCluster.h:38
void PrintRecoHitInfo(int hitInd, int digiSize, TVector3 MCPos, double MCEnergy, bool pixel) const
Int_t i
Definition: run_full.C:25
Class to store the Digis which belong to one cluster This class holds the information which Digi belo...
Definition: PndSdsCluster.h:19
double getDip() const
Definition: GFTrackCand.h:115
PndTransMap * map
Definition: sim_emc_apd.C:99
FairHit * GetFairHit(Int_t detId, Int_t hitId)
std::map< int, std::vector< int > > AssignHitsToTracks()
unsigned int getNHits() const
Definition: GFTrackCand.h:113
int n
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
TVector3 GetPositionOut() const
Definition: PndSdsMCPoint.h:91
double getCurv() const
Definition: GFTrackCand.h:114
Double_t p
Definition: anasim.C:58
void PrintTrackCand(GFTrackCand *cand) const
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
const TVectorD & n() const
std::vector< int > fGhostCand
std::map< int, std::vector< int > > fTrackPixHitIdMap
TClonesArray * fStripCluster
void getHit(unsigned int i, unsigned int &detId, unsigned int &hitId) const
Get detector ID and cluster index (hitId) for hit number i.
Definition: GFTrackCand.h:84
std::map< int, std::vector< int > > fTrackStripHitIdMap
void refit(bool withErrorCalc=true)
void Print(Int_t iTrack=0) const
Definition: PndMCTrack.cxx:95
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
void PrintClusterDigiInfo(int clIndex, std::vector< Int_t > digiInd, bool pixel)
void SetVerbose(int i)
TClonesArray * fTrackCand
std::vector< int > GetClusters(int MCHit, bool pixel)
TPad * p2
Definition: hist-t7.C:117
TVector3 GetPosition() const
Definition: PndSdsMCPoint.h:90
int hit(Int_t nEvents=0, TString inFile="sim.root", TString parFile="par.root", TString inDigi="digi.root", TString outFile="hit.root", Int_t timeBased=0)
Definition: hit.C:1
void GetTrackCandsForMCTrack(std::vector< int > pixHitId, std::vector< int > stripHitId, std::vector< int > &matches, std::vector< int > &result)
TPad * p1
Definition: hist-t7.C:116
void szFit(bool withErrorCalc=true)
void addHit(PndRiemannHit &hit)
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
int GetRecoHit(int clIndex, bool pixel) const
TClonesArray * fPixCluster
virtual void Print(const Option_t *opt=0) const
TClonesArray * fMCTracks
TClonesArray * fMCHits
void PndMvdEventAnaTask::Finish ( )
virtual

Definition at line 430 of file PndMvdEventAnaTask.cxx.

References fNCompleteTracks, fNGhostTracks, fNNotFoundPossibleTracks, fNNotFoundTracks, fNPartTracks, fNPossibleTracks, and fNTracks.

431 {
432  std::cout << std::endl;
433  std::cout << "--------- Summary Track Fitting -----------" << std::endl;
434  std::cout << "Number of Tracks: " << fNTracks << std::endl;
435  std::cout << "Number of Tracks with more than 3 hits: " << fNPossibleTracks << std::endl;
436  std::cout << (double)fNCompleteTracks / fNTracks * 100<< " % found completely" << std::endl;
437  std::cout << (double)fNCompleteTracks / fNPossibleTracks * 100<< " % of possible tracks found completely" << std::endl;
438  std::cout << (double)fNPartTracks / fNTracks * 100<< " % found partially" << std::endl;
439  std::cout << (double)fNPartTracks / fNPossibleTracks * 100<< " % of possible tracks found partially" << std::endl;
440  std::cout << (double)fNNotFoundPossibleTracks / fNPossibleTracks * 100 << " % of possible tracks not found" << std::endl;
441  std::cout << (double)fNNotFoundTracks / fNTracks * 100<< " % of Tracks not found in total" << std::endl;
442  std::cout << (double)fNGhostTracks / fNTracks * 100<< " % of Ghost Tracks" << std::endl;
443 }
std::vector< int > PndMvdEventAnaTask::GetClusters ( int  MCHit,
bool  pixel 
)
private

Definition at line 462 of file PndMvdEventAnaTask.cxx.

References fPixCluster, fStripCluster, and MCHitBelongsToCluster().

Referenced by Exec().

463 {
464  TClonesArray* cluster;
465 
466  std::vector<int> result;
467 
468  if (pixel == true){
469  cluster = fPixCluster;
470  }
471  else{
472  cluster = fStripCluster;
473  }
474 
475  if (cluster != 0) {
476  for (int clIndex = 0; clIndex < cluster->GetEntriesFast(); clIndex++){ //go through all pixel clusters
477  PndSdsCluster* myCluster = (PndSdsCluster*)(cluster->At(clIndex));
478  if (MCHitBelongsToCluster(MCHit, myCluster, pixel)){
479  result.push_back(clIndex);
480  }
481  }
482  }
483  return result;
484 }
Class to store the Digis which belong to one cluster This class holds the information which Digi belo...
Definition: PndSdsCluster.h:19
TClonesArray * fStripCluster
bool MCHitBelongsToCluster(int HitIndex, PndSdsCluster *cluster, bool pixCluster)
TClonesArray * fPixCluster
FairHit * PndMvdEventAnaTask::GetFairHit ( Int_t  detId,
Int_t  hitId 
)
private

Definition at line 446 of file PndMvdEventAnaTask.cxx.

References fPixReco, fStripReco, and p.

Referenced by Exec().

447 {
448  FairHit* p;
449 
450  if (detId == 5){
451  return (FairHit *) fPixReco->At(hitId);
452  }
453  else if (detId == 4){
454  return p=(FairHit *) fStripReco->At(hitId);
455  }
456 
457  std::cout << "-E- FairRiemannTrackCandDraw::GetFairHit : Unknown Detector with ID: " << detId << std::endl;
458  return 0;
459 }
Double_t p
Definition: anasim.C:58
TClonesArray * fStripReco
TClonesArray * fPixReco
int PndMvdEventAnaTask::GetRecoHit ( int  clIndex,
bool  pixel 
) const
private

Definition at line 486 of file PndMvdEventAnaTask.cxx.

References fPixReco, fStripReco, PndSdsHit::GetBotIndex(), PndSdsHit::GetClusterIndex(), and reco().

Referenced by Exec().

487 {
488  TClonesArray* reco;
489 
490  if (pixel == true){
491  reco = fPixReco;
492  }
493  else{
494  reco = fStripReco;
495  }
496  if (reco != 0){
497  for (int hitIndex = 0; hitIndex < reco->GetEntriesFast(); hitIndex++){
498  PndSdsHit* myHit = (PndSdsHit*)reco->At(hitIndex);
499  if ((myHit->GetClusterIndex() == clIndex) || myHit->GetBotIndex() == clIndex){ //test if RecoHit belongs to cluster
500  return hitIndex;
501  }
502  }
503  }
504  return -1;
505 }
TClonesArray * fStripReco
Int_t GetBotIndex() const
Definition: PndSdsHit.h:96
int reco()
Int_t GetClusterIndex() const
Definition: PndSdsHit.h:94
TClonesArray * fPixReco
void PndMvdEventAnaTask::GetTrackCandsForMCTrack ( std::vector< int >  pixHitId,
std::vector< int >  stripHitId,
std::vector< int > &  matches,
std::vector< int > &  result 
)
private

Definition at line 601 of file PndMvdEventAnaTask.cxx.

References fGhostCand, fTrackCand, GFTrackCand::getHit(), GFTrackCand::getNHits(), and i.

Referenced by Exec().

603 {
604  unsigned int detId, hitId;
605  int hitMatch = 0;
606 
607  if (fTrackCand != 0){
608  for (int i = 0; i < fTrackCand->GetEntriesFast(); i++){
609  GFTrackCand* myTrackCand = (GFTrackCand*)fTrackCand->At(i);
610  hitMatch = 0;
611  for (unsigned int j = 0; j < myTrackCand->getNHits(); j++){
612  myTrackCand->getHit(j, detId, hitId);
613  //std::cout << "DetId: " << detId << " HitId: " << hitId << std::endl;
614  if (detId == 5){
615  for (unsigned int k = 0; k < pixHitId.size(); k++){
616  //std::cout << "PixHitId: " << pixHitId[k] << std::endl;
617  if ((int)hitId == pixHitId[k])
618  hitMatch++;
619  }
620  }
621  else if (detId == 4){
622  for (unsigned int k = 0; k < stripHitId.size(); k++){
623  //std::cout << "StripHitId: " << stripHitId[k] << std::endl;
624  if ((int)hitId == stripHitId[k])
625  hitMatch++;
626  }
627  }
628  }
629  //std::cout << "internal HitMatch: " << hitMatch << std::endl;
630  if (hitMatch > 2){
631  matches.push_back(hitMatch);
632  result.push_back(i);
633  }
634  if (fGhostCand[i] < hitMatch)
635  fGhostCand[i] = hitMatch;
636  }
637  }
638 }
Int_t i
Definition: run_full.C:25
unsigned int getNHits() const
Definition: GFTrackCand.h:113
std::vector< int > fGhostCand
void getHit(unsigned int i, unsigned int &detId, unsigned int &hitId) const
Get detector ID and cluster index (hitId) for hit number i.
Definition: GFTrackCand.h:84
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
TClonesArray * fTrackCand
InitStatus PndMvdEventAnaTask::Init ( )
virtual

Definition at line 51 of file PndMvdEventAnaTask.cxx.

References fHDigisPerCluster, fHDigisPerClusterStrip, fHEnergyPerHit, fHEnergyRes, fHEnergyResStrip, fHHitsPerTrack, fHPointRes, fHPointResD, fHPointResDStrip, fHPointResM, fHPointResMStrip, fHPointResS, fHPointResSStrip, fHPointResStrip, fHPRes, fHPtRes, fHRiemannFakes, fHRiemannRes, fHRiemannTracksPerTrack, fHRiemannTracksPerTrackAdd, fHRiemannVertexResolutionX, fHRiemannVertexResolutionY, fHRiemannVertexResolutionZ, fHTracksPerEvent, fMCHits, fMCTracks, fPixCluster, fPixDigis, fPixReco, fStripCluster, fStripDigis, fStripReco, and fTrackCand.

52 {
53  // Get RootManager
54  FairRootManager* ioman = FairRootManager::Instance();
55  if ( !ioman){
56  std::cout << "-E- PndMvdEventAnaTask::Init: "<< "RootManager not instantiated!" << std::endl;
57  return kFATAL;
58  }
59 
60  // Get input arrays
61  fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
62  if (!fMCTracks){
63  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MCTrack" << " array!" << std::endl;
64  return kERROR;
65  }
66 
67  fMCHits = (TClonesArray*) ioman->GetObject("MVDPoint");
68  if ( !fMCHits) {
69  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDPoint"<<" array!" << std::endl;
70  return kERROR;
71  }
72 
73  fStripDigis = (TClonesArray*) ioman->GetObject("MVDStripDigis");
74  if ( !fStripDigis){
75  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDStripDigis"<<" array!" << std::endl;
76  }
77 
78  fPixDigis = (TClonesArray*) ioman->GetObject("MVDPixelDigis");
79  if ( !fPixDigis){
80  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDPixDigis"<<" array!" << std::endl;
81  }
82 
83  fPixReco = (TClonesArray*) ioman->GetObject("MVDHitsPixel");
84  if ( !fPixReco) {
85  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDHitsPixel"<<" array!" << std::endl;
86  }
87 
88  fStripReco = (TClonesArray*) ioman->GetObject("MVDHitsStrip");
89  if ( !fStripReco) {
90  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDHitsStrip" << " array!" << std::endl;
91  }
92 
93  fPixCluster = (TClonesArray*) ioman->GetObject("MVDClusterCand");
94  if ( !fPixCluster){
95  std::cout << "-W- PndMvdEventAnaTask::Init: "<< "No MVDClusterCand" <<" array!" << std::endl;
96  }
97 
98  fStripCluster = (TClonesArray*) ioman->GetObject("MVDStripClusterCand");
99  if ( !fStripCluster) {
100  std::cout << "-W- PndMvdEventAnaTask::Init: " << "No MVDStripClusterCand" <<" array!" << std::endl;
101  }
102 
103  fTrackCand = (TClonesArray*) ioman->GetObject("MVDRiemannTrackCand");
104  if ( !fTrackCand) {
105  std::cout << "-W- PndMvdEventAnaTask::Init: " << "No MVDRiemannTrackCand" <<" array!" << std::endl;
106  fTrackCand = (TClonesArray*) ioman->GetObject("MVDIdealTrackCand");
107  if (!fTrackCand)
108  std::cout << "-W- PndMvdEventAnaTask::Init: " << "No MVDIdealTrackCand" <<" array!" << std::endl;
109 
110  }
111 
112  fHTracksPerEvent = new TH1I("HTracksPerEvent", "Tracks per Event", 21, -0.5, 20.5);
113  fHHitsPerTrack = new TH1I("HHitsPerTrack", "Hits per Track", 21,0.5,20.5);
114  fHEnergyPerHit = new TH1D("HEnergyPerHit","Energy per Hit", 200, 0, 10000000);
115 
116  fHPointRes = new TH1D("HPointRes", "Point Resolution", 1000, 0,0.1);
117  fHPointRes->SetLineColor(1);
118  fHPointResS = new TH1D("HPointResS", "Point Resolution", 1000, 0,0.1);
119  fHPointResS->SetLineColor(3);
120  fHPointResD = new TH1D("HPointResD", "Point Resolution", 1000, 0,0.1);
121  fHPointResD->SetLineColor(4);
122  fHPointResM = new TH1D("HPointResM", "Point Resolution", 1000, 0,0.1);
123  fHPointResM->SetLineColor(5);
124 
125  fHDigisPerCluster = new TH1I("HDigisPerCluster","Digis per Cluster", 20, 0, 20);
126  fHEnergyRes = new TH1D("HEnergyRes", "Energy Resolution", 1000, -100000, 100000);
127  fHPRes = new TH1D("HPRes","Momentum Resolution", 1000, -5,5);
128  fHPtRes = new TH1D("HPtRes","Transversal Momentum Resolution", 1000, -5,5);
129 
130  fHPointResStrip = new TH1D("HPointResStrip", "Point Resolution Strips", 1000, 0,0.1);
131  fHPointResStrip->SetLineColor(6);
132  fHPointResSStrip = new TH1D("HPointResSStrip", "Point Resolution Strip", 1000, 0,0.1);
133  fHPointResSStrip->SetLineColor(7);
134  fHPointResDStrip = new TH1D("HPointResDStrip", "Point Resolution Strip", 1000, 0,0.1);
135  fHPointResDStrip->SetLineColor(8);
136  fHPointResMStrip = new TH1D("HPointResMStrip", "Point Resolution Strip", 1000, 0,0.1);
137  fHPointResMStrip->SetLineColor(9);
138 
139  fHDigisPerClusterStrip = new TH1I("HDigisPerClusterStrip","Digis per Cluster Strips", 20, 0, 20);
140  fHDigisPerClusterStrip->SetLineColor(2);
141  fHEnergyResStrip = new TH1D("HEnergyResStrip", "Energy Resolution Strips", 1000, -100000, 100000);
142  fHEnergyResStrip->SetLineColor(2);
143 
144  fHRiemannRes = new TH1I("HRiemannRes", "Quality of the Riemann Trackfinding", 16,-5.5,10.5);
145  fHRiemannFakes = new TH1I("HRiemannFakes", "Wrong assign hits in Track", 16, -5.5,10.5);
146  fHRiemannFakes->SetLineColor(2);
147  fHRiemannTracksPerTrack = new TH1I("HRiemannTracksPerTrack", "Found Tracks per MC Track", 11, -0.5, 10.5);
148  fHRiemannTracksPerTrackAdd = new TH1I("HRiemannTracksPerTrackAdd", "Found Tracks per MC Track", 11, -0.5, 10.5);
149  fHRiemannTracksPerTrackAdd->SetLineColor(2);
150 
151  fHRiemannVertexResolutionX = new TH1D("HRiemannVertexResolutionX","Difference between reco vertex and MC vertex", 1000, -1, 1);
152  fHRiemannVertexResolutionY = new TH1D("HRiemannVertexResolutionY","Difference between reco vertex and MC vertex", 1000, -1, 1);
153  fHRiemannVertexResolutionZ = new TH1D("HRiemannVertexResolutionZ","Difference between reco vertex and MC vertex", 1000, -1, 1);
154 
155  return kSUCCESS;
156 }
TClonesArray * fStripDigis
TClonesArray * fStripCluster
TClonesArray * fStripReco
TClonesArray * fTrackCand
TClonesArray * fPixReco
TClonesArray * fPixCluster
TClonesArray * fMCTracks
TClonesArray * fMCHits
TClonesArray * fPixDigis
bool PndMvdEventAnaTask::MCHitBelongsToCluster ( int  HitIndex,
PndSdsCluster cluster,
bool  pixCluster 
)
private

Definition at line 583 of file PndMvdEventAnaTask.cxx.

References fPixDigis, fStripDigis, PndSdsCluster::GetClusterList(), PndSdsDigi::GetIndex(), PndSdsDigi::GetNIndices(), and i.

Referenced by GetClusters().

584 {
585  bool result = false;
586  std::vector<Int_t> clusterList = cluster->GetClusterList();
587  for (unsigned int i = 0; i < clusterList.size() && result == false; i++){
588  PndSdsDigi* myDigi;
589  if (pixCluster)
590  myDigi = (PndSdsDigi*)(fPixDigis->At(clusterList[i]));
591  else
592  myDigi = (PndSdsDigi*)(fStripDigis->At(clusterList[i]));
593  for (int j = 0; j < myDigi->GetNIndices(); j++)
594  if (myDigi->GetIndex(j) == HitIndex)
595  result = true;
596  }
597  return result;
598 }
std::vector< Int_t > GetClusterList() const
Definition: PndSdsCluster.h:38
TClonesArray * fStripDigis
Base class for Digi information.
Definition: PndSdsDigi.h:29
Int_t i
Definition: run_full.C:25
Int_t GetIndex(int i=0) const
Definition: PndSdsDigi.h:63
Int_t GetNIndices() const
Definition: PndSdsDigi.h:64
TClonesArray * fPixDigis
void PndMvdEventAnaTask::PrintClusterDigiInfo ( int  clIndex,
std::vector< Int_t >  digiInd,
bool  pixel 
)
private

Definition at line 507 of file PndMvdEventAnaTask.cxx.

References fPixDigis, fPrintCluster, fPrintPixDigis, fStripDigis, PndSdsDigi::Print(), and s.

Referenced by Exec().

509 {
510  TClonesArray* digis;
511 
512  if (pixel){
513  digis = fPixDigis;
514  }
515  else{
516  digis = fStripDigis;
517  }
518 
519  if (digis == 0)
520  fPrintPixDigis = false;
521 
522  if (fPrintCluster)
523  std::cout << "Hit belongs to cluster " << clIndex << std::endl;
524  if (fPrintPixDigis) {
525  for (unsigned int s = 0; s < digiInd.size(); s++) {
526  PndSdsDigi* myDigi = (PndSdsDigi*)(digis->At(digiInd[s]));
527  if (pixel) ((PndSdsDigiPixel*)myDigi)->Print();
528  else ((PndSdsDigiStrip*)myDigi)->Print();
529  }
530  }
531 }
TClonesArray * fStripDigis
Base class for Digi information.
Definition: PndSdsDigi.h:29
TLorentzVector s
Definition: Pnd2DStar.C:50
Class for digitised strip hits.
virtual void Print()
Definition: PndSdsDigi.h:87
Data class to store the digi output of a pixel module.
TClonesArray * fPixDigis
void PndMvdEventAnaTask::PrintRecoHitInfo ( int  hitInd,
int  digiSize,
TVector3  MCPos,
double  MCEnergy,
bool  pixel 
) const
private

Definition at line 533 of file PndMvdEventAnaTask.cxx.

References fHEnergyRes, fHEnergyResStrip, fHPointRes, fHPointResD, fHPointResDStrip, fHPointResM, fHPointResMStrip, fHPointResS, fHPointResSStrip, fHPointResStrip, fPixReco, fPrintPixHit, fStripReco, PndSdsHit::GetEloss(), PndSdsHit::GetPosition(), PndSdsHit::Print(), and reco().

Referenced by Exec().

534 {
535  TClonesArray* reco;
536  TH1* hPointRes;
537  TH1* hPointResS;
538  TH1* hPointResD;
539  TH1* hPointResM;
540  TH1* hEnergyRes;
541 
542  if (pixel){
543  reco = fPixReco;
544 
545  hPointRes = fHPointRes;
546  hPointResS = fHPointResS;
547  hPointResD = fHPointResD;
548  hPointResM = fHPointResM;
549  hEnergyRes = fHEnergyRes;
550  }
551  else{
552  reco = fStripReco;
553 
554  hPointRes = fHPointResStrip;
555  hPointResS = fHPointResSStrip;
556  hPointResD = fHPointResDStrip;
557  hPointResM = fHPointResMStrip;
558  hEnergyRes = fHEnergyResStrip;
559  }
560 
561  PndSdsHit* myHit = (PndSdsHit*)reco->At(hitInd);
562  if (fPrintPixHit)
563  {
564  if (pixel) std::cout << "PixHit: ";
565  else std::cout << "StripHit: ";
566  std::cout << hitInd << std::endl; //write out RecoHit
567  myHit->Print();
568  }
569  TVector3 RecoPos = myHit->GetPosition();
570  double RecoEnergy = myHit->GetEloss()*10E9;
571  TVector3 result = MCPos - RecoPos;
572  ((TH1D*)hPointRes)->Fill(result.Mag());
573  if (digiSize == 1)
574  hPointResS->Fill(result.Mag());
575  if (digiSize == 2)
576  hPointResD->Fill(result.Mag());
577  if (digiSize > 2)
578  hPointResM->Fill(result.Mag());
579  ((TH1D*)hEnergyRes)->Fill(MCEnergy-RecoEnergy);
580 }
virtual void Print(const Option_t *opt=0) const
Definition: PndSdsHit.cxx:66
TVector3 GetPosition() const
Definition: PndSdsHit.h:93
Double_t GetEloss() const
Definition: PndSdsHit.h:97
TClonesArray * fStripReco
int reco()
TClonesArray * fPixReco
void PndMvdEventAnaTask::PrintTrackCand ( GFTrackCand cand) const
private

Definition at line 652 of file PndMvdEventAnaTask.cxx.

References GFTrackCand::getHit(), GFTrackCand::getNHits(), hit(), and i.

Referenced by Exec().

653 {
654  unsigned int det, hit;
655  //std::cout << "TrackCand: " << cand->getCurv() << " curv, " << cand->getDip() << " dip, " << (int)(cand->inverted()) << " inverted." << "\n";
656  for (unsigned int i = 0; i < cand->getNHits(); i++){
657  cand->getHit(i, det, hit);
658  std::cout << det << "/" << hit << " ";
659  }
660  std::cout << "\n";
661 }
Int_t i
Definition: run_full.C:25
unsigned int getNHits() const
Definition: GFTrackCand.h:113
void getHit(unsigned int i, unsigned int &detId, unsigned int &hitId) const
Get detector ID and cluster index (hitId) for hit number i.
Definition: GFTrackCand.h:84
int hit(Int_t nEvents=0, TString inFile="sim.root", TString parFile="par.root", TString inDigi="digi.root", TString outFile="hit.root", Int_t timeBased=0)
Definition: hit.C:1
void PndMvdEventAnaTask::ProduceHits ( )
private
void PndMvdEventAnaTask::Register ( )
private
void PndMvdEventAnaTask::Reset ( )
private
void PndMvdEventAnaTask::SetParContainers ( )
virtual

Virtual method Init

Definition at line 158 of file PndMvdEventAnaTask.cxx.

159 {
160  // Get Base Container
161 // FairRun* ana = FairRun::Instance();
162 // FairRuntimeDb* rtdb=ana->GetRuntimeDb();
163 
164 }

Member Data Documentation

int PndMvdEventAnaTask::fEventNr
private

Definition at line 145 of file PndMvdEventAnaTask.h.

Referenced by Exec().

std::vector<int> PndMvdEventAnaTask::fGhostCand
private

Definition at line 94 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and GetTrackCandsForMCTrack().

TH1* PndMvdEventAnaTask::fHDigisPerCluster
private

Definition at line 104 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHDigisPerClusterStrip
private

Definition at line 113 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHEnergyPerHit
private

Definition at line 98 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHEnergyRes
private

Definition at line 105 of file PndMvdEventAnaTask.h.

Referenced by DrawEnergyRes(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHEnergyResStrip
private

Definition at line 114 of file PndMvdEventAnaTask.h.

Referenced by DrawEnergyResStrip(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHHitsPerTrack
private

Definition at line 97 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHPointRes
private

Definition at line 100 of file PndMvdEventAnaTask.h.

Referenced by DrawPointRes(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResD
private

Definition at line 102 of file PndMvdEventAnaTask.h.

Referenced by DrawPointResD(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResDStrip
private

Definition at line 111 of file PndMvdEventAnaTask.h.

Referenced by Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResM
private

Definition at line 103 of file PndMvdEventAnaTask.h.

Referenced by DrawPointResM(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResMStrip
private

Definition at line 112 of file PndMvdEventAnaTask.h.

Referenced by Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResS
private

Definition at line 101 of file PndMvdEventAnaTask.h.

Referenced by DrawPointResS(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResSStrip
private

Definition at line 110 of file PndMvdEventAnaTask.h.

Referenced by Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPointResStrip
private

Definition at line 109 of file PndMvdEventAnaTask.h.

Referenced by DrawPointResStrip(), Init(), and PrintRecoHitInfo().

TH1* PndMvdEventAnaTask::fHPRes
private

Definition at line 106 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHPtRes
private

Definition at line 107 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannFakes
private

Definition at line 117 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannRes
private

Definition at line 116 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannTracksPerTrack
private

Definition at line 118 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannTracksPerTrackAdd
private

Definition at line 119 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannVertexResolutionX
private

Definition at line 121 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannVertexResolutionY
private

Definition at line 122 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHRiemannVertexResolutionZ
private

Definition at line 123 of file PndMvdEventAnaTask.h.

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

TH1* PndMvdEventAnaTask::fHTracksPerEvent
private

Definition at line 96 of file PndMvdEventAnaTask.h.

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

TClonesArray* PndMvdEventAnaTask::fMCHits
private

Definition at line 85 of file PndMvdEventAnaTask.h.

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

TClonesArray* PndMvdEventAnaTask::fMCTracks
private

Definition at line 84 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Init().

int PndMvdEventAnaTask::fNCompleteTracks
private

Definition at line 139 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNGhostTracks
private

Definition at line 143 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNNotFoundPossibleTracks
private

Definition at line 141 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNNotFoundTracks
private

Definition at line 142 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNPartTracks
private

Definition at line 140 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNPossibleTracks
private

Definition at line 138 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

int PndMvdEventAnaTask::fNTracks
private

Definition at line 137 of file PndMvdEventAnaTask.h.

Referenced by Exec(), and Finish().

TClonesArray* PndMvdEventAnaTask::fPixCluster
private

Definition at line 90 of file PndMvdEventAnaTask.h.

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

TClonesArray* PndMvdEventAnaTask::fPixDigis
private

Definition at line 86 of file PndMvdEventAnaTask.h.

Referenced by Init(), MCHitBelongsToCluster(), and PrintClusterDigiInfo().

TClonesArray* PndMvdEventAnaTask::fPixReco
private

Definition at line 88 of file PndMvdEventAnaTask.h.

Referenced by GetFairHit(), GetRecoHit(), Init(), and PrintRecoHitInfo().

bool PndMvdEventAnaTask::fPrintCluster
private

Definition at line 128 of file PndMvdEventAnaTask.h.

Referenced by PrintClusterDigiInfo().

bool PndMvdEventAnaTask::fPrintGhosts
private

Definition at line 135 of file PndMvdEventAnaTask.h.

Referenced by Exec().

bool PndMvdEventAnaTask::fPrintMCHit
private

Definition at line 127 of file PndMvdEventAnaTask.h.

Referenced by Exec().

bool PndMvdEventAnaTask::fPrintPixDigis
private

Definition at line 129 of file PndMvdEventAnaTask.h.

Referenced by PrintClusterDigiInfo().

bool PndMvdEventAnaTask::fPrintPixHit
private

Definition at line 130 of file PndMvdEventAnaTask.h.

Referenced by PrintRecoHitInfo().

bool PndMvdEventAnaTask::fPrintStripCluster
private

Definition at line 131 of file PndMvdEventAnaTask.h.

bool PndMvdEventAnaTask::fPrintStripDigis
private

Definition at line 132 of file PndMvdEventAnaTask.h.

bool PndMvdEventAnaTask::fPrintStripHit
private

Definition at line 133 of file PndMvdEventAnaTask.h.

bool PndMvdEventAnaTask::fPrintTrack
private

Definition at line 126 of file PndMvdEventAnaTask.h.

Referenced by Exec().

bool PndMvdEventAnaTask::fPrintTrackMatch
private

Definition at line 134 of file PndMvdEventAnaTask.h.

Referenced by Exec().

TClonesArray* PndMvdEventAnaTask::fStripCluster
private

Definition at line 91 of file PndMvdEventAnaTask.h.

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

TClonesArray* PndMvdEventAnaTask::fStripDigis
private

Definition at line 87 of file PndMvdEventAnaTask.h.

Referenced by Init(), MCHitBelongsToCluster(), and PrintClusterDigiInfo().

TClonesArray* PndMvdEventAnaTask::fStripReco
private

Definition at line 89 of file PndMvdEventAnaTask.h.

Referenced by GetFairHit(), GetRecoHit(), Init(), and PrintRecoHitInfo().

TClonesArray* PndMvdEventAnaTask::fTrackCand
private

Definition at line 92 of file PndMvdEventAnaTask.h.

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

std::map<int, std::vector<int> > PndMvdEventAnaTask::fTrackPixHitIdMap
private

Definition at line 147 of file PndMvdEventAnaTask.h.

Referenced by Exec().

std::map<int, std::vector<int> > PndMvdEventAnaTask::fTrackStripHitIdMap
private

Definition at line 148 of file PndMvdEventAnaTask.h.

Referenced by Exec().


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