17 #if !defined(__CINT__) || defined(__MAKECINT__) 20 #include <Riostream.h> 24 #include <TParticle.h> 25 #include <TRefArray.h> 58 TH1F * hEta = (TH1F*)
new TH1F(
"hEta",
"reco #eta",1000, -0.71,0.71);
59 TH1F * hPhi = (TH1F*)
new TH1F(
"hPhi",
"reco #phi",360, 0,360);
60 TH1F * hE = (TH1F*)
new TH1F(
"hE" ,
"reco e",300, 0,30);
61 TH1F * hTime = (TH1F*)
new TH1F(
"hTime" ,
"reco time",1000, 0,1000);
62 hEta ->SetXTitle(
"#eta");
63 hPhi ->SetXTitle(
"#phi (deg)");
64 hE ->SetXTitle(
"E (GeV)");
65 hTime->SetXTitle(
"time (ns)");
68 TH1F * hMCEta = (TH1F*)
new TH1F(
"hMCEta",
"MC #eta",1000, -0.71,0.71);
69 TH1F * hMCPhi = (TH1F*)
new TH1F(
"hMCPhi",
"MC #phi",360, 0,360);
70 TH1F * hMCE = (TH1F*)
new TH1F(
"hMCE" ,
"MC e",300, 0,30);
71 hMCEta->SetXTitle(
"#eta");
72 hMCPhi->SetXTitle(
"#phi (deg)");
73 hMCE ->SetXTitle(
"E (GeV)");
76 TH1F * hMCDEta = (TH1F*)
new TH1F(
"hMCDEta",
" #eta cluster - #eta MC",500, -0.05,0.05);
77 TH1F * hMCDPhi = (TH1F*)
new TH1F(
"hMCDPhi",
" #phi cluster - #phi MC",500, -0.05,0.05);
78 TH1F * hMCDE = (TH1F*)
new TH1F(
"hMCDE" ,
"e cluster - e MC",200, -10,10);
79 hMCDEta->SetXTitle(
"#eta_{reco}-#eta_{MC}");
80 hMCDPhi->SetXTitle(
"#phi_{reco}-#phi_{MC} (rad)");
81 hMCDE ->SetXTitle(
"E_{reco}-E_{MC} (GeV)");
84 TH1F * hTMDEta = (TH1F*)
new TH1F(
"hTMDEta",
" #eta cluster - eta track",500, -0.05,0.05);
85 TH1F * hTMDPhi = (TH1F*)
new TH1F(
"hTMDPhi",
" #phi cluster - phi track",500, -0.05,0.05);
86 TH1F * hTMDEtaOut= (TH1F*)
new TH1F(
"hTMDEtaOut",
" #eta cluster - eta track-outer",500, -0.05,0.05);
87 TH1F * hTMDPhiOut= (TH1F*)
new TH1F(
"hTMDPhiOut",
" #phi cluster - phi track-outer",500, -0.05,0.05);
88 TH1F * hTMResEta = (TH1F*)
new TH1F(
"hTMResEta",
" Residual #eta cluster - eta track",500, -0.05,0.05);
89 TH1F * hTMResPhi = (TH1F*)
new TH1F(
"hTMResPhi",
" Residual #phi cluster - phi track",500, -0.05,0.05);
90 hTMDEta ->SetXTitle(
"#eta_{cluster}-#eta_{track}");
91 hTMDPhi ->SetXTitle(
"#phi_{cluster}-#phi_{track} (rad)");
92 hTMDEtaOut->SetXTitle(
"#eta_{cluster}-#eta_{track-out}");
93 hTMDPhiOut->SetXTitle(
"#phi_{cluster}-#phi_{track-out} (rad)");
95 TH1F * hTMEOverP = (TH1F*)
new TH1F(
"hTMEOverP",
" E_{cluster}/ p_{Track}",100,0,2);
96 TH1F * hTMEOverPOut = (TH1F*)
new TH1F(
"hTMEOverPOut",
" E_{cluster}/ p_{Track-out}",100,0,2);
97 hTMEOverP ->SetXTitle(
"E_{cluster}/ p_{Track}");
98 hTMEOverPOut ->SetXTitle(
"E_{cluster}/ p_{Track-out}");
101 TH2I * hEGAPatch =
new TH2I(
"hEGAPatch",
"EGA trigger",51,-0.5,50.5,65,-0.5,64.5);
102 TH2I * hEJEPatch =
new TH2I(
"hEJEPatch",
"EJE trigger",51,-0.5,50.5,65,-0.5,64.5);
103 hEGAPatch->SetXTitle(
"column (#eta direction)");
104 hEGAPatch->SetYTitle(
"row (#phi direction)");
105 hEJEPatch->SetXTitle(
"column (#eta direction)");
106 hEJEPatch->SetYTitle(
"row (#phi direction)");
109 TFile*
f =
new TFile(
"AliESDs.root");
110 TTree* esdTree = (TTree*)f->Get(
"esdTree");
123 for(Int_t mod = 0; mod < (geom->
GetEMCGeometry())->GetNumberOfSuperModules(); mod++)
144 TRefArray* caloClusters =
new TRefArray();
146 Int_t nEvt = esdTree->GetEntries();
148 for(Int_t iev = 0; iev < nEvt; iev++)
150 cout <<
"<<<< Event: " << iev+1 <<
"/" << nEvt <<
" >>>>";
154 esdTree->GetEvent(iev);
181 for (Int_t icell= 0; icell < nTotalCells; icell++)
184 cout<<
"Cell : " <<icell<<
"/"<<nTotalCells <<
" - ID: "<<cells.
GetCellNumber(icell)<<
"; High Gain? "<<cells.
GetHighGain(icell);
209 if (col > -1 && row > -1)
218 Bool_t isEGA1 = ((bit>> bitEGA ) & 0x1);
219 Bool_t isEGA2 = ((bit>> (bitEGA+1)) & 0x1);
220 Bool_t isEJE1 = ((bit>> bitEJE ) & 0x1);
221 Bool_t isEJE2 = ((bit>> (bitEJE+1)) & 0x1);
225 printf(
"Bit %d (G1 %d,G2 %d,J1 %d,J2 %d), cell in patch position (ieta,iphi)=(%d,%d), L1 amplitude %d\n ",
226 bit, isEGA1, isEGA2, isEJE1, isEJE2, col, row, ts);
227 if(bit>10) hEJEPatch->Fill(col,row);
228 else hEGAPatch->Fill(col,row);
242 Int_t nclus = caloClusters->GetEntries();
243 for (Int_t icl = 0; icl < nclus; icl++)
246 Float_t energy = clus->
E();
248 TVector3 vpos(pos[0],pos[1],pos[2]);
256 Double_t cphi = vpos.Phi();
257 if(cphi < 0) cphi +=TMath::TwoPi();
258 Double_t ceta = vpos.Eta();
263 Int_t labelIndex = clus->
GetLabel();
269 hPhi ->Fill(cphi*TMath::RadToDeg());
271 hTime->Fill(clus->
GetTOF()*1e9);
276 cout <<
"Cluster: " << icl+1 <<
"/" << nclus <<
" Energy: " << energy <<
"; Phi: " 277 << cphi*TMath::RadToDeg() <<
"; Eta: " << ceta
278 <<
"; NCells: " << nCells <<
"; NLM: " << nlm
279 <<
"; #Labels: " << nLabels <<
" Index: " << labelIndex
280 <<
"; Time "<<clus->
GetTOF()*1e9<<
" ns "<<endl;
286 printf(
"\t N matches %d, Residual phi %2.4f, eta %2.4f\n",
295 if(labelIndex >= 0 && labelIndex < stack->GetNtrack())
297 TParticle * particle = stack->
Particle(labelIndex);
299 hMCEta ->Fill(particle->Eta());
300 hMCPhi ->Fill(particle->Phi()*TMath::RadToDeg());
301 hMCE ->Fill(particle->Energy());
302 hMCDEta ->Fill(ceta-particle->Eta());
303 hMCDPhi ->Fill(cphi-particle->Phi());
304 hMCDE ->Fill(energy-particle->Energy());
307 cout<<
" MC More contributing primary: "<<particle->GetName()<<
"; with kinematics: "<<endl;
308 cout<<
" \t Energy: "<<particle->Energy()<<
"; Phi: "<<particle->Phi()*TMath::RadToDeg()<<
"; Eta: "<<particle->Eta()<<endl;
309 cout<<
" \t dE: "<<energy-particle->Energy()<<
"; dPhi: "<<(cphi-particle->Phi())*TMath::RadToDeg()<<
"; dEta: "<<ceta-particle->Eta()<<endl;
314 for(Int_t i = 1; i < nLabels; i++)
320 cout<<
" Other contributing primary: "<<particle->GetName()<<
"; Energy "<<particle->Energy()<<
324 else if( labelIndex >= stack->
GetNtrack())
325 cout <<
"PROBLEM, label is too large : "<<labelIndex<<
" >= particles in stack "<< stack->
GetNtrack() <<endl;
327 cout<<
"Negative label!!! : "<<labelIndex<<endl;
336 Double_t tphi = track->
Phi();
337 Double_t teta = track->
Eta();
338 Double_t tmom = track->
P();
340 Double_t deta = teta - ceta;
341 Double_t dphi = tphi - cphi;
343 if(dphi > TMath::Pi()) dphi -= 2*TMath::Pi();
344 if(dphi < -TMath::Pi()) dphi += 2*TMath::Pi();
346 Double_t dR = sqrt(dphi*dphi + deta*deta);
348 Double_t pOverE = tmom/energy;
350 cout <<
"\t Track Momentum : " << tmom <<
" phi: " << tphi <<
" eta: " << teta << endl;
352 hTMEOverP->Fill(pOverE);
353 hTMDEta ->Fill(deta);
354 hTMDPhi ->Fill(dphi);
362 cout <<
"\t Track Momentum - Outer: " << tmom <<
" phi: " << tphi <<
" eta: " << teta << endl;
367 if(dphi > TMath::Pi()) dphi -= 2*TMath::Pi();
368 if(dphi < -TMath::Pi()) dphi += 2*TMath::Pi();
370 dR = sqrt(dphi*dphi + deta*deta);
372 pOverE = tmom/energy;
374 if(dR < 0.02 && pOverE < 1.8 && nCells > 1)
376 cout <<
"\t Excellent MATCH! dR = " << dR <<
" p/E = " << pOverE <<
"Residuals phi "<< dphi*TMath::RadToDeg()<<
" eta "<< deta << endl;
379 hTMEOverPOut->Fill(pOverE);
380 hTMDEtaOut->Fill(deta);
381 hTMDPhiOut->Fill(dphi);
383 else printf(
"!!! NO OUTER PARAM !!!\n");
385 else printf(
"!!! NO TRACK !!!\n");
392 printf(
"PID weights: ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f, hadrons: pion %0.2f, kaon %0.2f, proton %0.2f , neutron %0.2f, kaon %0.2f \n",
406 for(Int_t i = 0; i < nCells ; i++)
408 Int_t absId = index[i];
409 Double_t ampFract = fraction[i];
414 cout<<
"\t Cluster Cell: AbsID : "<< absId <<
" == "<<clus->
GetCellAbsId(i) <<
"; Amplitude "<< amp;
415 cout<<
"; Fraction "<<ampFract<<
"; Time " <<time*1e9<<
"; High Gain? "<<hg<<endl;
423 printf(
"\t \t cell MC label index %d, cell MC deposited energy map %d: p1 %2.2f, p2 %2.2f, p3 %2.2f p4 %2.2f\n",
424 cells.
GetCellMCLabel(absId),map,eDepFrac[0],eDepFrac[1],eDepFrac[2],eDepFrac[3]);
440 iIphi, iIeta,iphi,ieta);
445 cout<<
" SModule "<<iSupMod<<
"; Tower "<<iTower <<
"; Eta "<<iIeta
446 <<
"; Phi "<<iIphi<<
"; Index: Cell Eta "<<ieta<<
"; Cell Phi "<<iphi
447 <<
"; Global: Cell Eta "<<cellEta<<
"; Cell Phi "<<cellPhi*TMath::RadToDeg()
450 if ( i == 0 ) sm = iSupMod;
451 else if ( sm != iSupMod)
printf(
" ******CLUSTER SHARED BY 2 SuperModules!!!!\n");
461 TFile * fhisto = (TFile*)
new TFile(
"histos.root",
"recreate");
483 hTMEOverPOut->Write();
484 hTMDEtaOut ->Write();
485 hTMDPhiOut ->Write();
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
static const TString & GetDefaultEventFolderName()
Bool_t kPrintMisalMatrix
Print the alignment matrices stored in ESDs.
virtual void GetCellMCEdepFractionArray(Int_t, Float_t *) const
virtual UShort_t * GetClusterMCEdepFraction() const
virtual Double_t * GetCellsAmplitudeFraction()
virtual Double_t E() const
virtual Double_t GetCellTime(Short_t cellNumber)=0
virtual Short_t GetCellNumber(Short_t pos) const =0
const AliExternalTrackParam * GetOuterParam() const
virtual Short_t GetNumberOfCells() const =0
void GetTriggerBits(Int_t &bits) const
Get the trigger bits for a given fastor position.
virtual Int_t GetNTracksMatched() const
Bool_t kPrintCaloTrigger
Print trigger patches information.
Virtual class for calorimeter cell data handling.
void SetMisalMatrix(const TGeoHMatrix *m, Int_t smod) const
void EtaPhiFromIndex(Int_t absId, Double_t &eta, Double_t &phi) const
virtual Bool_t Next()
Forward to next trigger entry (fastor / L0 patch)
virtual Int_t GetCellMCLabel(Short_t cellNumber)=0
virtual Double_t GetTrackDz(void) const
virtual UShort_t * GetCellsAbsId()
Virtual class for calorimeter cluster data handling.
Bool_t kPrintClusterCells
Print cells in clusters information.
void GetPosition(Int_t &col, Int_t &row) const
Access to position of the current fastor channel.
static AliRunLoader * Open(const char *filename="galice.root", const char *eventfoldername=AliConfig::GetDefaultEventFolderName(), Option_t *option="READ")
void GetCellPhiEtaIndexInSModule(Int_t nSupMod, Int_t nModule, Int_t nIphi, Int_t nIeta, Int_t &iphi, Int_t &ieta) const
void Print(const char *method, TStopwatch &timer, Int_t n)
virtual Int_t GetMCLabel(Short_t pos) const =0
virtual Int_t GetNtrack() const
virtual Double_t GetAmplitude(Short_t pos) const =0
virtual Bool_t GetCellHighGain(Short_t cellNumber)=0
AliESDCaloTrigger * GetCaloTrigger(TString calo) const
Bool_t GetCellIndex(Int_t absId, Int_t &nSupMod, Int_t &nModule, Int_t &nIphi, Int_t &nIeta) const
AliESDCaloCells * GetEMCALCells() const
Bool_t kPrintCaloCluster
Print cluster parameters.
Container with calorimeter trigger information in the ESD event.
TString GetFiredTriggerClasses() const
Int_t GetEMCALClusters(TRefArray *clusters) const
virtual Int_t * GetLabels() const
Bool_t kPrintCaloCells
Print cells parameters.
virtual Double_t GetCellAmplitude(Short_t cellNumber)=0
virtual Double_t GetTime(Short_t pos) const =0
virtual const Double_t * GetPID() const
virtual Double_t GetTrackDx(void) const
virtual UInt_t GetNLabels() const
virtual Int_t GetLabel() const
Int_t GetEvent(Int_t evno)
virtual Double_t GetTOF() const
TParticle * Particle(Int_t id, Bool_t useInEmbedding=kFALSE)
virtual UChar_t GetNExMax() const
void ReadFromTree(TTree *tree, Option_t *opt="")
AliEMCALEMCGeometry * GetEMCGeometry() const
virtual Int_t GetTrackMatchedIndex(Int_t=0) const
Only for ESDs.
virtual Double_t GetEFraction(Short_t pos) const =0
virtual UInt_t * GetCellsMCEdepFractionMap() const
AliESDtrack * GetTrack(Int_t i) const
const TGeoHMatrix * GetEMCALMatrix(Int_t i) const
static AliEMCALGeometry * GetInstance()
virtual Int_t GetNCells() const
virtual Bool_t GetHighGain(Short_t pos) const =0
void GetL1TimeSum(Int_t ×um) const
Get the L1 time sums (L1 ADC values) for the current fastor.
Bool_t kPrintClusterPID
Print clusters PID (bayesian) weights.
EMCal geometry, singleton.
Bool_t kPrintTrackMatches
Print cluster-track matching information.
Bool_t kPrintKine
Print MC related information. Do not use for raw data.
virtual Int_t GetCellAbsId(Int_t) const
Int_t LoadKinematics(Option_t *option="READ")
virtual void GetPosition(Float_t *) const