]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/Hdisp.C
eff C++ warnings corrected.
[u/mrichter/AliRoot.git] / HMPID / Hdisp.C
CommitLineData
e9f11028 1#include <TSystem.h>
2#include <TFile.h>
3#include <TTree.h>
4#include <TButton.h>
5#include <TCanvas.h>
6#include <TClonesArray.h>
7#include <TParticle.h>
8#include <TRandom.h>
9#include <TLatex.h>
10#include <TLegend.h>
11#include <TPolyMarker.h>
12#include <TBox.h>
13#include <AliESD.h>
14#include <AliCDBManager.h>
15#include <AliCDBEntry.h>
16#include "AliHMPIDHit.h"
17#include "AliHMPIDv2.h"
18#include "AliHMPIDReconstructor.h"
19#include "AliHMPIDRecon.h"
20#include "AliHMPIDParam.h"
21#include "AliHMPIDCluster.h"
22#include <TChain.h>
126b3e0e 23
e9f11028 24TCanvas *fCanvas=0; Int_t fType=3; Int_t fEvt=0; Int_t fNevt=0;
25TFile *fHitFile; TTree *fHitTree; TClonesArray *fHitLst; TPolyMarker *fRenMip[7]; TPolyMarker *fRenCko[7]; TPolyMarker *fRenFee[7];
26 TClonesArray *fSdiLst;
27TFile *fDigFile; TTree *fDigTree; TObjArray *fDigLst; TPolyMarker *fRenDig[7];
28TFile *fCluFile; TTree *fCluTree; TObjArray *fCluLst; TPolyMarker *fRenClu[7];
29TFile *fEsdFile; TTree *fEsdTree; AliESD *fEsd; TPolyMarker *fRenTxC[7]; TPolyMarker *fRenRin[7];
30TFile *fCosFile; TTree *fCosTree;
8282f9d0 31
e9f11028 32//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
33void CreateContainers()
34{//to create all containers
35 fHitLst=new TClonesArray("AliHMPIDHit");
36 fSdiLst=new TClonesArray("AliHMPIDDigit");
37 fDigLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fDigLst->AddAt(new TClonesArray("AliHMPIDDigit"),i); fDigLst->SetOwner(kTRUE);
38 fCluLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fCluLst->AddAt(new TClonesArray("AliHMPIDCluster"),i); fCluLst->SetOwner(kTRUE);
39 fEsd =new AliESD;
126b3e0e 40}
126b3e0e 41//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 42void CreateRenders()
43{
44 for(Int_t ch=0;ch<7;ch++){
45 fRenMip[ch]=new TPolyMarker; fRenMip[ch]->SetMarkerStyle(kOpenTriangleUp); fRenMip[ch]->SetMarkerColor(kRed);
46 fRenCko[ch]=new TPolyMarker; fRenCko[ch]->SetMarkerStyle(kOpenCircle); fRenCko[ch]->SetMarkerColor(kRed);
47 fRenFee[ch]=new TPolyMarker; fRenFee[ch]->SetMarkerStyle(kOpenDiamond); fRenFee[ch]->SetMarkerColor(kRed);
48 fRenDig[ch]=new TPolyMarker; fRenDig[ch]->SetMarkerStyle(kOpenSquare); fRenDig[ch]->SetMarkerColor(kGreen);
49 fRenClu[ch]=new TPolyMarker; fRenClu[ch]->SetMarkerStyle(kStar); fRenClu[ch]->SetMarkerColor(kBlue);
50 fRenTxC[ch]=new TPolyMarker; fRenTxC[ch]->SetMarkerStyle(kPlus); fRenTxC[ch]->SetMarkerColor(kRed); fRenTxC[ch]->SetMarkerSize(3);
51 fRenRin[ch]=new TPolyMarker; fRenRin[ch]->SetMarkerStyle(kFullDotSmall); fRenRin[ch]->SetMarkerColor(kMagenta);
52 }
53}//CreateRenders()
54//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
55void ClearRenders()
56{
57 for(Int_t ch=0;ch<7;ch++){
58 fRenTxC[ch]->SetPolyMarker(0);
59 fRenRin[ch]->SetPolyMarker(0);
60 fRenMip[ch]->SetPolyMarker(0);
61 fRenCko[ch]->SetPolyMarker(0);
62 fRenFee[ch]->SetPolyMarker(0);
63 fRenDig[ch]->SetPolyMarker(0);
64 fRenClu[ch]->SetPolyMarker(0);
65 }
66}//ClearRenders()
67//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
68void PrintHits()
69{
70//Prints a list of HMPID hits for a given event. Default is event number 0.
71 if(!fHitTree) return;
72 Printf("List of HMPID hits for event %i",fEvt);
73 Int_t iTot=0;
74 for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){//entries loop
75 fHitTree->GetEntry(iEnt);
76 fHitLst->Print();
77 iTot+=fHitLst->GetEntries();
78 }
79 Printf("totally %i hits for event %i",iTot,fEvt);
80}//PrintHits();
81//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
82void PrintSdis()
83{//prints a list of HMPID sdigits for a given event
84 Printf("List of HMPID sdigits for event %i",fEvt);
85 fSdiLst->Print();
86 Printf("totally %i sdigits for event %i",fSdiLst->GetEntries(),fEvt);
87}//PrintSdis()
88//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
89void PrintDigs()
90{//prints a list of HMPID digits
91 Printf("List of HMPID digits for event %i",fEvt);
92 fDigLst->Print();
93 Int_t iTot=0; for(Int_t i=0;i<7;i++) {iTot+=((TClonesArray*)fDigLst->At(i))->GetEntries();}
94 Printf("totally %i digits for event %i",iTot,fEvt);
8282f9d0 95}
96//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 97void PrintClus()
98{//prints a list of HMPID clusters for a given event
99 Printf("List of HMPID clusters for event %i",fEvt);
8282f9d0 100
e9f11028 101 fCluLst->Print();
102 Int_t iTot=0; for(Int_t iCh=0;iCh<7;iCh++) iTot+=((TClonesArray*)fCluLst->At(iCh))->GetEntries();
103
104 Printf("totally %i clusters for event %i",iTot,fEvt);
8282f9d0 105}
106//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 107void PrintEsd()
108{//prints a list of HMPID Esd for a given event
109 Printf("List of HMPID ESD summary for event %i",fEvt);
110 for(Int_t iTrk=0;iTrk<fEsd->GetNumberOfTracks();iTrk++){
111 AliESDtrack *pTrk = fEsd->GetTrack(iTrk);
112 Float_t x,y;Int_t q,nacc; pTrk->GetHMPIDmip(x,y,q,nacc);
113 Printf("Track %02i with p %7.2f with ThetaCer %5.3f with %i photons",iTrk,pTrk->GetP(),pTrk->GetHMPIDsignal(),nacc);
114 }
115}//PrintEsd()
116//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117void DrawChamber(Int_t iCh)
118{//used by Draw() to Draw() chamber structure
119 gPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5);
120 if(iCh>=0){TLatex txt; txt.SetTextSize(0.1); txt.DrawLatex(-5,-5,Form("%i",iCh));}
121
122 for(Int_t iPc=AliHMPIDDigit::kMinPc;iPc<=AliHMPIDDigit::kMaxPc;iPc++){
123 TBox *pBox=new TBox(AliHMPIDDigit::MinPcX(iPc),AliHMPIDDigit::MinPcY(iPc),
124 AliHMPIDDigit::MaxPcX(iPc),AliHMPIDDigit::MaxPcY(iPc));
125 pBox->SetFillStyle(0); pBox->Draw();
126 }//PC loop
127}//DrawChamber()
128//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129void DrawLegend()
130{//used by Draw() to draw legend
131 Int_t nTxC=0,nMip=0,nCko=0,nFee=0,nDig=0,nClu=0;
132 for(Int_t ch=0;ch<7;ch++){
133 nTxC+=fRenTxC[ch]->GetN();
134 nMip+=fRenMip[ch]->GetN();
135 nCko+=fRenCko[ch]->GetN();
136 nFee+=fRenFee[ch]->GetN();
137 nDig+=fRenDig[ch]->GetN();
138 nClu+=fRenClu[ch]->GetN();
139 }
140 TLegend *pLeg=new TLegend(0.2,0.2,0.8,0.8);
141 pLeg->SetHeader(Form("Event %i Total %i",fEvt,fNevt));
142 pLeg->AddEntry(fRenTxC[0],Form("TRKxPC %i" ,nTxC),"p");
143 pLeg->AddEntry(fRenMip[0],Form("Mip hits %i" ,nMip),"p");
144 pLeg->AddEntry(fRenCko[0],Form("Ckov hits %i" ,nCko),"p");
145 pLeg->AddEntry(fRenFee[0],Form("Feed hits %i" ,nFee),"p");
146 pLeg->AddEntry(fRenDig[0],Form("Digs %i" ,nDig),"p");
147 pLeg->AddEntry(fRenClu[0],Form("Clus %i" ,nClu),"p");
148 pLeg->Draw();
149}//DrawLegend()
150//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
151void Draw()
152{//draws all the objects of current event in given canvas
153 for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
154 switch(iCh){
155 case 6: fCanvas->cd(1); break; case 5: fCanvas->cd(2); break;
156 case 4: fCanvas->cd(4); break; case 3: fCanvas->cd(5); break; case 2: fCanvas->cd(6); break;
157 case 1: fCanvas->cd(8); break; case 0: fCanvas->cd(9); break;
158 }
159 gPad->SetEditable(kTRUE); gPad->Clear(); DrawChamber(iCh);
160 fRenTxC[iCh]->Draw();
161 fRenMip[iCh]->Draw();
162 fRenFee[iCh]->Draw();
163 fRenCko[iCh]->Draw();
164 fRenRin[iCh]->Draw("CLP");
165 fRenDig[iCh]->Draw();
166 fRenClu[iCh]->Draw();
167 gPad->SetEditable(kFALSE);
168 }//chambers loop
8282f9d0 169
e9f11028 170 fCanvas->cd(3); gPad->Clear(); DrawLegend();
171}//Draw()
172//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
173void RenderHit(TClonesArray *pHitLst)
174{//used by ReadEvent() and SimulateEvent() to render hits to polymarker structures, one per chamber
175 for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){ //hits loop
176 AliHMPIDHit *pHit = (AliHMPIDHit*)pHitLst->At(iHit); Int_t ch=pHit->Ch(); Float_t x=pHit->LorsX(); Float_t y=pHit->LorsY(); //get current hit
177 switch(pHit->Pid()){
178 case 50000050: fRenCko[ch]->SetNextPoint(x,y);break;
179 case 50000051: fRenFee[ch]->SetNextPoint(x,y);break;
180 default: fRenMip[ch]->SetNextPoint(x,y);break;
181 }//switch hit PID
182 }//hits loop for this entry
183}//RenderHits()
184//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
185void RenderClu(TObjArray *pClus,TPolyMarker **pMark)
186{//used by ReadEvent() and SimulateEvent() to render clusters to polymarker structures, one per chamber
187 for(Int_t iCh=0;iCh<=6;iCh++){ //chambers loop
188 TClonesArray *pClusCham=(TClonesArray*)pClus->At(iCh); //get clusters list for this chamber
189 for(Int_t iClu=0;iClu<pClusCham->GetEntries();iClu++){ //clusters loop
190 AliHMPIDCluster *pClu = (AliHMPIDCluster*)pClusCham->At(iClu); //get current cluster
191 pMark[iCh]->SetNextPoint(pClu->X(),pClu->Y());
192 }//switch hit PID
193 }//hits loop for this entry
194}//RenderClus()
195//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
196void RenderEsd(AliESD *pEsd)
197{//used by ReadEvent() or SimulateEvent() to render ESD to polymarker structures for rings and intersections one per chamber
198 AliHMPIDRecon rec;
199 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop to collect cerenkov rings and intersection points
200 AliESDtrack *pTrk=pEsd->GetTrack(iTrk); Int_t ch=pTrk->GetHMPIDcluIdx(); //get track and chamber intersected by it
201 if(ch<0) continue; //this track does not intersect any chamber
202 Float_t thRa,phRa,xRa,yRa; pTrk->GetHMPIDtrk(xRa,yRa,thRa,phRa); //get info on current track
203 ch/=1000000; //actual chamber number
204 Float_t xPc=0,yPc=0; AliHMPIDTracker::IntTrkCha(pTrk,xPc,yPc); //find again intersection of track with PC--> it is not stored in ESD!
205 fRenTxC[ch]->SetNextPoint(xPc,yPc); //add this intersection point
206 Float_t ckov=pTrk->GetHMPIDsignal(); //get ckov angle stored for this track
207 if(ckov>0){
208 rec.SetTrack(xRa,yRa,thRa,phRa);
209 for(Int_t j=0;j<100;j++){
210 TVector2 pos; pos=rec.TracePhot(ckov,j*0.0628);
211 if(!AliHMPIDDigit::IsInDead(pos.X(),pos.Y())) fRenRin[ch]->SetNextPoint(pos.X(),pos.Y());
212 }
213 }//if ckov is valid
214 }//tracks loop
215}//RenEsd()
a1d55ff3 216//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 217void SimulateEsd(AliESD *pEsd)
a1d55ff3 218{
219 TParticle part; TLorentzVector mom;
50e69536 220 for(Int_t iTrk=0;iTrk<100;iTrk++){//stack loop
a1d55ff3 221 part.SetPdgCode(kProton);
50e69536 222 part.SetProductionVertex(0,0,0,0);
223 Double_t eta= -0.2+gRandom->Rndm()*0.4; //rapidity is random [-0.2,+0.2]
224 Double_t phi= gRandom->Rndm()*60.*TMath::DegToRad(); //phi is random [ 0 , 60 ] degrees
225 mom.SetPtEtaPhiM(5,eta,phi,part.GetMass());
a1d55ff3 226 part.SetMomentum(mom);
227 AliESDtrack trk(&part);
228 pEsd->AddTrack(&trk);
229 }//stack loop
230}//EsdFromStack()
231//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 232void SimulateHits(AliESD *pEsd, TClonesArray *pHits)
233{//used by SimulateEvent to simulate hits out from provided ESD
50e69536 234 const Int_t kCerenkov=50000050;
235 const Int_t kFeedback=50000051;
f3bae3e2 236
a1d55ff3 237 AliHMPIDRecon rec;
e9f11028 238 Float_t eMip=200e-9,ePho=7.5e-9;
a591e55f 239 Int_t hc=0;
a1d55ff3 240 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop
241 AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
e9f11028 242 Float_t xRa,yRa;
243 Int_t ch=AliHMPIDTracker::IntTrkCha(pTrk,xRa,yRa);
a1d55ff3 244 if(ch<0) continue; //this track does not hit HMPID
50e69536 245 Float_t beta = pTrk->GetP()/(TMath::Sqrt(pTrk->GetP()*pTrk->GetP()+0.938*0.938));
246 Float_t ckov=TMath::ACos(1./(beta*1.292));
a1d55ff3 247
e9f11028 248 Float_t theta,phi,xPc,yPc,; pTrk->GetHMPIDtrk(xPc,yPc,theta,phi); rec.SetTrack(xRa,yRa,theta,phi);
a1d55ff3 249
e9f11028 250 if(!AliHMPIDDigit::IsInDead(xPc,yPc)) new((*pHits)[hc++]) AliHMPIDHit(ch,eMip,kProton ,iTrk,xPc,yPc); //mip hit
50e69536 251 Int_t nPhots = (Int_t)(20.*TMath::Power(TMath::Sin(ckov),2)/TMath::Power(TMath::Sin(TMath::ACos(1./1.292)),2));
252 for(int i=0;i<nPhots;i++){
a591e55f 253 TVector2 pos;
254 pos=rec.TracePhot(ckov,gRandom->Rndm()*TMath::TwoPi());
e9f11028 255 if(!AliHMPIDDigit::IsInDead(pos.X(),pos.Y())) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kCerenkov,iTrk,pos.X(),pos.Y());
a1d55ff3 256 } //photon hits
e9f11028 257 for(int i=0;i<3;i++){//feedback photons
258 Float_t x=gRandom->Rndm()*160; Float_t y=gRandom->Rndm()*150;
259 if(!AliHMPIDDigit::IsInDead(x,y)) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kFeedback,iTrk,x,y); //feedback hits
260 }//photon hits loop
a1d55ff3 261 }//tracks loop
e9f11028 262}//SimulateHits()
a1d55ff3 263//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 264void DoZoom(Int_t evt, Int_t px, Int_t py, TObject *)
a1d55ff3 265{
8282f9d0 266 if(evt!=5 && evt!=6) return; //5- zoom in 6-zoom out
267 const Int_t minZoom=64;
268 const Int_t maxZoom=2;
269 static Int_t zoom=minZoom; //zoom level
270 if(evt==5&&zoom==maxZoom) return;
271 if(evt==6&&zoom==minZoom) return;
272
d73ece6a 273 // if(!obj->IsA()->InheritsFrom("TPad")) return; //current object is not pad
50e69536 274 TVirtualPad *pPad=gPad->GetSelectedPad();
8282f9d0 275 if(pPad->GetNumber()==3 || pPad->GetNumber()==7) return; //current pad is wrong
276
277 // Printf("evt=%i (%i,%i) %s",evt,px,py,obj->GetName());
278
279 Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py);
6ed6a312 280
281 if(evt==5){ zoom=zoom/2; pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom in
282 else { zoom=zoom*2; pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom out
283 if(zoom==minZoom) pPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5);
8282f9d0 284 ((TCanvas *)gTQSender)->SetTitle(Form("zoom x%i",minZoom/zoom));
285 pPad->Modified();
286 pPad->Update();
a1d55ff3 287}
8282f9d0 288//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 289void ReadEvent()
290{//used by NextEvent() to read curent event and construct all render elements
291 if(fNevt && fEvt>=fNevt) fEvt=0; //loop over max event
d73ece6a 292
e9f11028 293 if(fHitFile){
294 if(fHitTree) delete fHitTree; fHitTree=(TTree*)fHitFile->Get(Form("Event%i/TreeH",fEvt));
295
296 fHitTree->SetBranchAddress("HMPID",&fHitLst);
297 for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){
298 fHitTree->GetEntry(iEnt);
299 RenderHit(fHitLst);
300 }//prim loop
301 }//if hits file
302
303 if(fCluFile){
304 if(fCluTree) delete fCluTree; fCluTree=(TTree*)fCluFile->Get(Form("Event%i/TreeR",fEvt));
305 for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++) fCluTree->SetBranchAddress(Form("HMPID%i",iCh),&(*fCluLst)[iCh]);
306 fCluTree->GetEntry(0);
307 RenderClu(fCluLst,fRenClu);
308 }//if clus file
309
310 if(fEsdFile){//if ESD file open
311 fEsdTree->GetEntry(fEvt);
312 RenderEsd(fEsd);
313 }//if ESD file
314}//ReadEvent()
d73ece6a 315//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 316void SimulateEvent()
317{// used by NextEvent() to simulate all info
318 AliCDBManager* pCDB = AliCDBManager::Instance(); pCDB->SetDefaultStorage("local://$HOME"); pCDB->SetRun(0);
319 AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean");
320
321 SimulateEsd(fEsd);
322 SimulateHits(fEsd,fHitLst);
323 AliHMPIDv2::Hit2Sdi(fHitLst,fSdiLst);
324 AliHMPIDDigitizer::Sdi2Dig(fSdiLst,fDigLst);
325 AliHMPIDReconstructor::Dig2Clu(fDigLst,fCluLst);
326 AliHMPIDTracker::Recon(fEsd,fCluLst,(TObjArray*)pNmeanEnt->GetObject());
327
328 RenderHit(fHitLst);
329 RenderClu(fCluLst,fRenClu);
330 RenderEsd(fEsd);
331}//SimulateEvent()
d73ece6a 332//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 333void NextEvent()
d73ece6a 334{
e9f11028 335 fEvt++;
336 ClearRenders();
337 switch(fType){
338 case 1: ReadEvent();break;
339// case 2: ReadCosmic(); break;
340 case 3: SimulateEvent(); break;
341 default: return;
342 }
343 Draw();
d73ece6a 344}
345//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
e9f11028 346void Hdisp()
347{//display events from files if any in current directory or simulated events
348 CreateContainers();
349 CreateRenders();
d73ece6a 350
e9f11028 351 TString title="Session with";
352 if(gSystem->IsFileInIncludePath("HMPID.Hits.root")){// tries to open hits
353 fHitFile=TFile::Open("HMPID.Hits.root"); fNevt=fHitFile->GetNkeys(); fType=1; title+=Form(" HITS-%i ",fNevt);
354 }
6b155609 355
e9f11028 356 if(gSystem->IsFileInIncludePath("HMPID.Digits.root")){// tries to open clusters
357 fDigFile=TFile::Open("HMPID.Digits.root"); fNevt=fDigFile->GetNkeys(); fType=1; title+=Form(" DIGITS-%i ",fNevt);
6b155609 358 }
359
e9f11028 360 if(gSystem->IsFileInIncludePath("HMPID.RecPoints.root")){// tries to open clusters
361 fCluFile=TFile::Open("HMPID.RecPoints.root"); fNevt=fCluFile->GetNkeys(); fType=1; title+=Form(" CLUSTERS-%i ",fNevt);
362 }
50e69536 363
e9f11028 364 if(gSystem->IsFileInIncludePath("AliESDs.root")){
365 fEsdFile=TFile::Open("AliESDs.root"); fEsdTree=(TTree*)fEsdFile->Get("esdTree"); fNevt=fEsdTree->GetEntries(); fType=1; title+=Form(" ESD-%i ",fNevt);
366 fEsdTree->SetBranchAddress("ESD", &fEsd);
367 }
368
369 if(gSystem->IsFileInIncludePath("cosmic.root")){ //clm: Check if cosmic file is in the folder
370 fCosFile=TFile::Open("cosmic.root"); fCosTree=(TTree*)fCosFile->Get("cosmic"); fNevt=fCosTree->GetEntries(); fType=2;
371 fCosTree->SetBranchAddress("Digs",&fDigLst); fCosTree->SetBranchAddress("Clus",&fCluLst);
372 }
373
374 fCanvas=new TCanvas("all","",600,400); fCanvas->Divide(3,3,0,0);// pAll->ToggleEditor();
375 fCanvas->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",0,0,"DoZoom(Int_t,Int_t,Int_t,TObject*)");
376 fCanvas->cd(7);
377 switch(fType){
378 case 1: fCanvas->SetTitle(title.Data());
379 TButton *pNxtBtn=new TButton("Next" ,"NextEvent()",0.0,0.0,0.3,0.1); pNxtBtn->Draw();
380 if(fHitFile){TButton *pHitBtn=new TButton("Print hits","PrintHits()",0.0,0.2,0.3,0.3); pHitBtn->Draw();}
381 if(fDigFile){TButton *pDigBtn=new TButton("Print digs","PrintDigs()",0.0,0.4,0.3,0.5); pDigBtn->Draw();}
382 if(fCluFile){TButton *pCluBtn=new TButton("Print clus","PrintClus()",0.0,0.6,0.3,0.7); pCluBtn->Draw();}
383 if(fEsdFile){TButton *pEsdBtn=new TButton("Print ESD" ,"PrintEsd()" ,0.0,0.8,0.3,0.9); pEsdBtn->Draw();}
384
385 break;
386 case 2: fCanvas->SetTitle("COSMIC");
387 TButton *pCosBtn=new TButton("Next" ,"NextEvent()",0.0,0.0,0.3,0.1); pCosBtn->Draw();
388
389 break;
390 case 3: fCanvas->SetTitle("SIMULATION");
391 TButton *pSimBtn=new TButton("Simulate" ,"NextEvent()",0.0,0.0,0.3,0.1); pSimBtn->Draw();
392 TButton *pHitBtn=new TButton("Print hits","PrintHits()",0.0,0.2,0.3,0.3); pHitBtn->Draw();
393 TButton *pDigBtn=new TButton("Print digs","PrintDigs()",0.0,0.4,0.3,0.5); pDigBtn->Draw();
394 TButton *pCluBtn=new TButton("Print clus","PrintClus()",0.0,0.6,0.3,0.7); pCluBtn->Draw();
395 TButton *pEsdBtn=new TButton("Print ESD" ,"PrintEsd()" ,0.0,0.8,0.3,0.9); pEsdBtn->Draw();
396 break;
397 }
398
399
0112d781 400 NextEvent();
50e69536 401
e9f11028 402
403}//Hdisp()
404//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++