FairRoot/PandaRoot
Functions
EventListing.C File Reference

Go to the source code of this file.

Functions

int EventListing ()
 

Function Documentation

int EventListing ( )

Definition at line 1 of file EventListing.C.

References creator, PndSdsCluster::DigiBelongsToCluster(), DigiFile, PndSdsDigi::GetIndex(), PndSdsDigi::GetNIndices(), h1, i, MCFile, PndSdsDigiPixel::Print(), PndSdsDigiStrip::Print(), PndSdsMCPoint::Print(), PndSdsHit::Print(), RecoFile, t, and TString.

2 {
3  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
4  gSystem->Load("libriemann.C");
5 
6  TString MCFile = "Mvd_Test.root";
7  TH1D* h1 = new TH1D("h1","PixError",100,-0.01,0.01);
8 
9  PndFileNameCreator creator(MCFile.Data());
10  TString DigiFile = creator.GetDigiFileName(false).c_str();
11  TString RecoFile = creator.GetRecoFileName(false).c_str();
12  //TString TrackFFile = creator.GetTrackFindingFileName(false).c_str();
13 
14  TFile* fMC = new TFile(MCFile.Data());
15  TFile* fDigi = new TFile(DigiFile.Data());
16  TFile* fReco = new TFile(RecoFile.Data());
17  //TFile* fTrack = new TFile(TrackFFile.Data());
18 
19  TTree* t = (TTree*)(fMC->Get("pndsim"));
20 
21  t->AddFriend("pndsim", fDigi);
22  t->AddFriend("pndsim", fReco);
23  //t->AddFriend("pndsim", fTrack);
24 
25  t->StartViewer();
26 
27  TClonesArray* MCHits = new TClonesArray("PndSdsMCPoint");
28  TClonesArray* PixDigis = new TClonesArray("PndSdsDigiPixel");
29  TClonesArray* StripDigis = new TClonesArray("PndSdsDigiStrip");
30  TClonesArray* PixReco = new TClonesArray("PndSdsHit");
31  TClonesArray* StripReco = new TClonesArray("PndSdsHit");
32  TClonesArray* PixCluster = new TClonesArray("PndSdsClusterPixel");
33  TClonesArray* StripCluster = new TClonesArray("PndSdsClusterStrip");
34  //TClonesArray* TrackCand = new TClonesArray("TrackCand");
35 
36  t->SetBranchAddress("MVDPoint", &MCHits);
37  t->SetBranchAddress("MVDPixelDigis", &PixDigis);
38  t->SetBranchAddress("MVDStripDigis", &StripDigis);
39  t->SetBranchAddress("MVDHitsPixel", &PixReco);
40  t->SetBranchAddress("MVDHitsStrip", &StripReco);
41  t->SetBranchAddress("MVDPixelClusterCand", &PixCluster);
42  t->SetBranchAddress("MVDStripClusterCand", &StripCluster);
43 
44  t->GetEntry(0);
45 
46  for (int i = 0; i < MCHits->GetEntriesFast(); i++){ //get all MC Hits
47  PndSdsMCPoint* myPoint = (PndSdsMCPoint*)(MCHits->At(i));
48  std::cout << "<<<<<<<<<<< MCPoint >>>>>>>>>> " << std::endl;
49  myPoint->Print(); //write out MC info
50 
51  for (int j = 0; j < PixDigis->GetEntriesFast(); j++){ //get all Digis
52  PndSdsDigiPixel* myPixDigi = (PndSdsDigiPixel*)PixDigis->At(j);
53  bool dig = false;
54  for (int ind = 0; ind < myPixDigi->GetNIndices(); ind++) //test if digi belongs to MCHit
55  if (myPixDigi->GetIndex(ind) == i) dig = true;
56  if (dig){
57  std::cout << "PixDigi: "; //write out DigiInfo
58  myPixDigi->Print();
59  for (int k = 0; k < PixCluster->GetEntriesFast(); k++){ //get all clusters
60  PndSdsClusterPixel* myPixCluster = (PndSdsClusterPixel*)PixCluster->At(k);
61  if (myPixCluster->DigiBelongsToCluster(j)){ //test if digi belongs to cluster
62  std::cout << "Digi " << j << " belongs to cluster: " << k << std::endl; //write out cluster info
63  for (int l = 0; l < PixReco->GetEntriesFast(); l++){ //get all RecoHits
64  PndSdsHit* myPixHit = (PndSdsHit*)PixReco->At(l);
65  if (myPixHit->GetRefIndex() == k){ //test if RecoHit belongs to cluster
66  std::cout << "PixHit: " << l << std::endl; //write out RecoHit
67  myPixHit->Print();
68  }
69  }
70  }
71  }
72  }
73  }
74 
75  for (int j = 0; j < StripDigis->GetEntriesFast(); j++){
76  PndSdsDigiStrip* myStripDigi = (PndSdsDigiStrip*)StripDigis->At(j);
77  dig = false;
78  for (int ind = 0; ind < myStripDigi->GetNIndices(); ind++)
79  if (myStripDigi->GetIndex(ind) == i) dig = true;
80  if (dig){
81  std::cout << "StripDigi: "; // << myStripDigi;
82  myStripDigi->Print();
83  for (int k = 0; k < StripCluster->GetEntriesFast(); k++){
84  PndSdsClusterStrip* myStripCluster = (PndSdsClusterStrip*)StripCluster->At(k);
85  if (myStripCluster->DigiBelongsToCluster(j)){
86  std::cout << "Digi " << j << " belongs to cluster: " << k << std::endl;
87  for (int l = 0; l < StripReco->GetEntriesFast(); l++){
88  PndSdsHit* myStripHit = (PndSdsHit*)StripReco->At(l);
89  if (myStripHit->GetRefIndex() == k){
90  std::cout << "StripHit: " << l << std::endl;
91  myStripHit->Print();
92  }
93  }
94  }
95  }
96  }
97  }
98  }
99  return 0;
100 
101 }
std::ostream & Print(std::ostream &out=std::cout) const
virtual void Print(const Option_t *opt=0) const
Definition: PndSdsHit.cxx:66
TString RecoFile
Int_t i
Definition: run_full.C:25
Int_t GetIndex(int i=0) const
Definition: PndSdsDigi.h:63
Class for digitised strip hits.
TString DigiFile
A simple class which adds the corresponding file extensions to a given base class.
PndMvdCreateDefaultApvMap * creator
TString MCFile
TTree * t
Definition: bump_analys.C:13
Data class to store the digi output of a pixel module.
bool DigiBelongsToCluster(Int_t digiIndex)
Int_t GetNIndices() const
Definition: PndSdsDigi.h:64
virtual void Print(const Option_t *opt=0) const