]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hdisp.C
Bug fix: correct deletion of fGRPList in case of CDB cache
[u/mrichter/AliRoot.git] / HMPID / Hdisp.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3 #include <TSystem.h>
4 #include <TFile.h>
5 #include <TVirtualX.h>
6 #include <TTree.h>
7 #include <TButton.h>
8 #include <TCanvas.h>
9 #include <TCanvasImp.h>
10 #include <TStyle.h>
11 #include <TString.h>
12 #include <TClonesArray.h>
13 #include <TParticle.h>
14 #include <TRandom.h>
15 #include <TLatex.h>
16 #include <TPDGCode.h>
17 #include <TLegend.h>
18 #include <TPolyMarker.h>
19 #include <TBox.h>
20 #include <AliESDEvent.h>
21 #include <AliCDBManager.h>
22 #include <AliCDBEntry.h>
23 #include "AliHMPIDHit.h"
24 #include "AliHMPIDv2.h"
25 #include "AliHMPIDReconstructor.h"
26 #include "AliHMPIDRecon.h"
27 #include "AliHMPIDParam.h"
28 #include "AliHMPIDCluster.h"
29 #include "AliTracker.h"
30 #include "AliStack.h"
31
32 #endif
33
34 AliHMPIDParam *fParam;
35 TDatabasePDG *fPdg;
36
37 TCanvas *fCanvas=0; Int_t fType=3; Int_t fEvt=-1; Int_t fNevt=0;                      
38 TCanvasImp *fCanvasImp;
39 TFile *fHitFile; TTree *fHitTree; TClonesArray *fHitLst; TPolyMarker *fRenMip[7]; TPolyMarker *fRenCko[7]; TPolyMarker *fRenFee[7];
40                                   TClonesArray *fSdiLst; 
41 TFile *fDigFile; TTree *fDigTree; TObjArray    *fDigLst; TBox *fRenDig[7][160*144]; TBox *fBox[7][160*144];   
42 TFile *fCluFile; TTree *fCluTree; TObjArray    *fCluLst; TPolyMarker *fRenClu[7];
43 TFile *fEsdFile; TTree *fEsdTree; AliESDEvent  *fEsd;    TPolyMarker *fRenTxC[7]; TPolyMarker *fRenRin[7];  
44 TFile *fCosFile; TTree *fCosTree;
45
46 TButton *fHitMipBok,*fHitCkoBok,*fHitFeeBok,*fDigBok,*fCluBok,*fEsdBok;
47
48 TString fStHitMip = "ON";
49 TString fStHitCko = "ON";
50 TString fStHitFee = "ON";
51 TString fStDig    = "ON";
52 TString fStClu    = "ON";
53 TString fStEsd    = "ON";
54
55 Int_t fTimeArrow = 1;
56
57 Int_t fMaxCharge = 200;  //maximum charge to saturate color to red
58
59 Int_t fChamN[9]={6,5,-1,4,3,2,-1,1,0};
60
61 Int_t fTotPads[7],fTotClus[7];
62
63 AliRunLoader *gAL=0; 
64
65 Int_t nDigs[7];
66
67 enum EObjectType {kHitMip=kOpenTriangleUp,kHitCko=kOpenCircle,kHitFee=kOpenDiamond,kCluster=kStar,kTrack=kPlus,kRing=kFullDotSmall};
68   //  kOpenTriangleUp  Hit Mip
69   //  kOpenCircle      Hit Ckov
70   //  kOpenDiamond     Hit Feedback
71   //  kStar            Cluster 
72   //  kPlus            Track
73   //  kFullDotSmall    Ring
74   
75 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
76 void CreateContainers()
77 {//to create all containers
78   fHitLst=new TClonesArray("AliHMPIDHit");
79   fSdiLst=new TClonesArray("AliHMPIDDigit");
80   fDigLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fDigLst->AddAt(new TClonesArray("AliHMPIDDigit"),i);       fDigLst->SetOwner(kTRUE);
81   fCluLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fCluLst->AddAt(new TClonesArray("AliHMPIDCluster"),i);     fCluLst->SetOwner(kTRUE); 
82   fEsd   =new AliESDEvent;
83 }
84 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
85 void CreateRenders()
86 {
87   for(Int_t ch=0;ch<7;ch++){
88     fRenMip[ch]=new TPolyMarker; fRenMip[ch]->SetMarkerStyle(kOpenTriangleUp);fRenMip[ch]->SetMarkerColor(kBlack)  ;
89     fRenCko[ch]=new TPolyMarker; fRenCko[ch]->SetMarkerStyle(kOpenCircle);    fRenCko[ch]->SetMarkerColor(kBlack)  ;
90     fRenFee[ch]=new TPolyMarker; fRenFee[ch]->SetMarkerStyle(kOpenDiamond);   fRenFee[ch]->SetMarkerColor(kBlack)  ;
91     fRenClu[ch]=new TPolyMarker; fRenClu[ch]->SetMarkerStyle(kStar);          fRenClu[ch]->SetMarkerColor(kMagenta);
92     fRenTxC[ch]=new TPolyMarker; fRenTxC[ch]->SetMarkerStyle(kPlus);          fRenTxC[ch]->SetMarkerColor(kRed)    ;
93                                  fRenTxC[ch]->SetMarkerSize(3);
94     fRenRin[ch]=new TPolyMarker; fRenRin[ch]->SetMarkerStyle(kFullDotSmall);  fRenRin[ch]->SetMarkerColor(kMagenta);
95
96     for(Int_t iDig=0;iDig<160*144;iDig++) {
97       fRenDig[ch][iDig] = new TBox;
98       fRenDig[ch][iDig]->SetFillStyle(1);
99       fBox[ch][iDig] = new TBox;
100       fBox[ch][iDig]->SetFillStyle(0);
101     }
102   }
103 }//CreateRenders()
104 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105 void ClearRenders()
106 {
107   for(Int_t ch=0;ch<7;ch++){
108     fRenTxC[ch]->SetPolyMarker(0); 
109     fRenRin[ch]->SetPolyMarker(0); 
110     fRenMip[ch]->SetPolyMarker(0); 
111     fRenCko[ch]->SetPolyMarker(0); 
112     fRenFee[ch]->SetPolyMarker(0); 
113     nDigs[ch] = 0;
114     fRenClu[ch]->SetPolyMarker(0); 
115   }
116 }//ClearRenders()
117 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
118 void PrintHits()
119 {
120 //Prints a list of HMPID hits for a given event. Default is event number 0.
121   if(!fHitTree) return;
122   Printf("List of HMPID hits for event %i",fEvt);
123   Int_t iTot=0;
124   for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){//entries loop
125     fHitTree->GetEntry(iEnt);      
126     fHitLst->Print();
127     iTot+=fHitLst->GetEntries();
128   }
129   Printf("totally %i hits for event %i",iTot,fEvt);
130 }//PrintHits();
131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
132 void PrintSdis()
133 {//prints a list of HMPID sdigits  for a given event
134   Printf("List of HMPID sdigits for event %i",fEvt);
135   fSdiLst->Print();
136   Printf("totally %i sdigits for event %i",fSdiLst->GetEntries(),fEvt);
137 }//PrintSdis()
138 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
139 void PrintDigs()
140 {//prints a list of HMPID digits
141   Printf("List of HMPID digits for event %i",fEvt);  
142   fDigLst->Print();
143   Int_t iTot=0;  for(Int_t i=0;i<7;i++) {iTot+=((TClonesArray*)fDigLst->At(i))->GetEntries();}
144   Printf("totally %i digits for event %i",iTot,fEvt);
145 }
146 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
147 void PrintClus()
148 {//prints a list of HMPID clusters  for a given event
149   Printf("List of HMPID clusters for event %i",fEvt);
150   
151   fCluLst->Print();
152   Int_t iTot=0; for(Int_t iCh=0;iCh<7;iCh++) iTot+=((TClonesArray*)fCluLst->At(iCh))->GetEntries();
153   
154   Printf("totally %i clusters for event %i",iTot,fEvt);
155 }
156 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
157 void PrintEsd()
158 {//prints a list of HMPID Esd  for a given event
159   Printf("List of HMPID ESD summary for event %i",fEvt);
160   for(Int_t iTrk=0;iTrk<fEsd->GetNumberOfTracks();iTrk++){
161     AliESDtrack *pTrk = fEsd->GetTrack(iTrk);
162     Float_t x,y;Int_t q,nacc;   pTrk->GetHMPIDmip(x,y,q,nacc);
163     Printf("Track %02i with p %7.2f with ThetaCer %5.3f with %i photons",iTrk,pTrk->GetP(),pTrk->GetHMPIDsignal(),nacc);
164   }  
165 }//PrintEsd()
166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167 void DrawChamber(Int_t iCh) 
168 {//used by Draw() to Draw() chamber structure
169   gPad->Range(-10,-10,AliHMPIDParam::SizeAllX()+5,AliHMPIDParam::SizeAllY()+5);
170   if(iCh>=0){TLatex txt; txt.SetTextSize(0.06); txt.DrawLatex(55,-9,Form("RICH %i",iCh));}
171   
172   for(Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++){
173     TBox *pBox=new TBox(AliHMPIDParam::MinPcX(iPc),AliHMPIDParam::MinPcY(iPc),
174                         AliHMPIDParam::MaxPcX(iPc),AliHMPIDParam::MaxPcY(iPc));
175     pBox->SetFillStyle(0);  pBox->Draw();
176   }//PC loop      
177 }//DrawChamber()
178 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179 void DrawLegend()
180 {//used by Draw() to draw legend
181   Int_t nTxC=0,nMip=0,nCko=0,nFee=0,nDig=0,nClu=0;
182   for(Int_t ch=0;ch<7;ch++){
183     nTxC+=fRenTxC[ch]->Size();
184     nMip+=fRenMip[ch]->Size();
185     nCko+=fRenCko[ch]->Size();
186     nFee+=fRenFee[ch]->Size();
187     nClu+=fRenClu[ch]->Size();
188     nDig+=nDigs[ch];
189   }
190   TLegend *pLeg=new TLegend(0.15,0.3,0.85,0.98);
191   pLeg->SetHeader(Form("Event %i Total %i",fEvt,fNevt));
192   pLeg->AddEntry(fRenTxC[0],   Form("TRKxPC %i"  ,nTxC),"p");
193   pLeg->AddEntry(fRenMip[0],   Form("Mip hits %i"  ,nMip),"p");    
194   pLeg->AddEntry(fRenCko[0],   Form("Ckov hits %i"  ,nCko),"p");    
195   pLeg->AddEntry(fRenFee[0],   Form("Feed hits %i"  ,nFee),"p");    
196   pLeg->AddEntry(fRenDig[0][0],Form("Digs %i"  ,nDig),"p");    
197   pLeg->AddEntry(fRenClu[0],   Form("Clus %i"  ,nClu),"p");    
198   pLeg->Draw();
199  
200 /*  TLegend *pLeg2 = new TLegend(0.4,0.087,0.8,0.97);
201   pLeg2->SetHeader("");
202   pLeg2->AddEntry(fRenMip[0],Form("Mip hits %i"   ,nMip),"p");
203   pLeg2->Draw();*/
204
205  TBox *pBox = new TBox(0.012,0.03,0.97,0.28);
206  pBox->Draw("l");
207
208  TText *pText = new TText(0.35,0.21, "ddl = (0....13)");
209  pText->SetTextSize(0.07);
210  pText->Draw();
211
212  TText *pText2 = new TText(0.09,0.13,"RICH n ---> ddl 2n (left) & 2n+1 (right)");
213  pText2->SetTextSize(0.07);
214  pText2->Draw();
215
216  TText *pText3 = new TText(0.1,0.05, "phcat (0,2,4) left,   phcat (1,3,5) right");
217  pText3->SetTextSize(0.07);
218  pText3->Draw();
219
220 }//DrawLegend()
221 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
222 void Draw()
223 {//draws all the objects of current event in given canvas
224   Int_t iPadN[7]={9,8,6,5,4,2,1};
225   Int_t sampleCol = fMaxCharge/gStyle->GetNumberOfColors();
226   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop    
227     fCanvas->cd(iPadN[iCh]);
228     gPad->SetEditable(kTRUE); gPad->Clear(); DrawChamber(iCh);
229     
230     if(fDigFile){
231       if(fStDig   =="ON") {
232         for(Int_t iDig=0;iDig<nDigs[iCh];iDig++) {
233           Int_t charge = fRenDig[iCh][iDig]->GetUniqueID();
234           Int_t color = charge/sampleCol;if(color>=gStyle->GetNumberOfColors()) color = gStyle->GetNumberOfColors()-1;
235           fRenDig[iCh][iDig]->SetFillColor(gStyle->GetColorPalette(color));
236           fRenDig[iCh][iDig]->Draw();
237           fBox[iCh][iDig]->Draw();
238         }
239       }
240     }
241     
242     if(fEsdFile){if(fStEsd   =="ON") fRenTxC[iCh]->Draw();}
243     
244     if(fHitFile){
245       if(fStHitMip=="ON") fRenMip[iCh]->Draw();
246       if(fStHitFee=="ON") fRenFee[iCh]->Draw();
247       if(fStHitCko=="ON") fRenCko[iCh]->Draw();
248     }
249     
250     if(fEsdFile){if(fStEsd   =="ON") fRenRin[iCh]->Draw();}
251     if(fCluFile){if(fStClu   =="ON") fRenClu[iCh]->Draw();}
252     gPad->SetEditable(kFALSE);
253   }//chambers loop
254   
255   fCanvas->cd(3);  gPad->Clear(); DrawLegend();
256   
257   fCanvas->cd(7);
258   
259   if(fHitFile){
260     fHitMipBok->SetTitle(Form("Mips  %s",fStHitMip.Data())); fHitMipBok->Modified();
261     fHitCkoBok->SetTitle(Form("Ckov  %s",fStHitCko.Data())); fHitCkoBok->Modified();
262     fHitFeeBok->SetTitle(Form("Fdbk  %s",fStHitFee.Data())); fHitFeeBok->Modified();
263   }
264   if(fDigFile){fDigBok->SetTitle(fStDig); fDigBok->Modified();}  
265   if(fCluFile){fCluBok->SetTitle(fStClu); fCluBok->Modified();} 
266   if(fEsdFile){fEsdBok->SetTitle(fStEsd); fEsdBok->Modified();}
267   
268 }//Draw()
269 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
270 void RenderHit(TClonesArray *pHitLst)
271 {//used by ReadEvent() and SimulateEvent() to render hits to polymarker structures, one per chamber
272   for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){       //hits loop
273     AliHMPIDHit *pHit = (AliHMPIDHit*)pHitLst->At(iHit); Int_t ch=pHit->Ch(); Float_t x=pHit->LorsX(); Float_t y=pHit->LorsY();    //get current hit        
274     switch(pHit->Pid()){
275       case 50000050: fRenCko[ch]->SetNextPoint(x,y);break; 
276       case 50000051: fRenFee[ch]->SetNextPoint(x,y);break;
277       default:       fRenMip[ch]->SetNextPoint(x,y);break;
278     }//switch hit PID      
279 //    Printf("----------ihit %i ch %i chhit %i MIP %i CKO %i FEE %i",iHit,ch,pHit->Pid());
280   }//hits loop for this entry
281 }//RenderHits()
282 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
283 void RenderDig(TObjArray *pDigLst)
284 {//used by ReadEvent() and SimulateEvent() to render digs to Tbox structures, one per chamber
285   for(Int_t iCh=0;iCh<=6;iCh++){                                    //chambers loop   
286     TClonesArray *pDigCham=(TClonesArray*)pDigLst->At(iCh);         //get digs list for this chamber
287     nDigs[iCh] = 0;
288     fTotPads[iCh] = pDigCham->GetEntries();
289     for(Int_t iDig=0;iDig<pDigCham->GetEntries();iDig++){            //digs loop
290       AliHMPIDDigit *pDig = (AliHMPIDDigit*)pDigCham->At(iDig); Float_t x=pDig->LorsX(); Float_t y=pDig->LorsY();    //get current hit        
291       fRenDig[iCh][iDig]->SetX1(x-0.5*AliHMPIDParam::SizePadX());
292       fRenDig[iCh][iDig]->SetX2(x+0.5*AliHMPIDParam::SizePadX());
293       fRenDig[iCh][iDig]->SetY1(y-0.5*AliHMPIDParam::SizePadY());
294       fRenDig[iCh][iDig]->SetY2(y+0.5*AliHMPIDParam::SizePadY());
295       fRenDig[iCh][iDig]->SetUniqueID((Int_t)pDig->Q());
296       fBox[iCh][iDig]->SetX1(x-0.5*AliHMPIDParam::SizePadX());
297       fBox[iCh][iDig]->SetX2(x+0.5*AliHMPIDParam::SizePadX());
298       fBox[iCh][iDig]->SetY1(y-0.5*AliHMPIDParam::SizePadY());
299       fBox[iCh][iDig]->SetY2(y+0.5*AliHMPIDParam::SizePadY());
300       nDigs[iCh]++;
301     }//Digit loop
302   }//hits loop for this entry
303 }//RenderHits()
304 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
305 void RenderClu(TObjArray *pClus)
306 {//used by ReadEvent() and SimulateEvent() to render clusters to polymarker structures, one per chamber
307   for(Int_t iCh=0;iCh<=6;iCh++){                                     //chambers loop   
308     TClonesArray *pClusCham=(TClonesArray*)pClus->At(iCh);           //get clusters list for this chamber
309     fTotClus[iCh] = pClusCham->GetEntries();
310     for(Int_t iClu=0;iClu<pClusCham->GetEntries();iClu++){           //clusters loop
311       AliHMPIDCluster *pClu = (AliHMPIDCluster*)pClusCham->At(iClu); //get current cluster        
312       fRenClu[iCh]->SetNextPoint(pClu->X(),pClu->Y());
313 //      Printf("RenderClu: ch %i x %f y %f",iCh,pClu->X(),pClu->Y());
314     }//cluster loop
315   }//chamber loop for this entry
316 }//RenderClus()
317 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
318 void RenderEsd(AliESDEvent *pEsd)
319 {//used by ReadEvent() or SimulateEvent() to render ESD to polymarker structures for rings and intersections one per chamber
320   AliHMPIDRecon rec;
321   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop to collect cerenkov rings and intersection points
322     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);    Int_t ch=pTrk->GetHMPIDcluIdx(); //get track and chamber intersected by it
323     if(ch<0) continue;                                                          //this track does not intersect any chamber
324     Float_t thRa,phRa,xRa,yRa; pTrk->GetHMPIDtrk(xRa,yRa,thRa,phRa);            //get info on current track
325     ch/=1000000;                            
326     Float_t xPc=0,yPc=0; AliHMPIDTracker::IntTrkCha(pTrk,xPc,yPc);              //find again intersection of track with PC--> it is not stored in ESD!
327     Int_t npTrk = fRenTxC[ch]->SetNextPoint(xPc,yPc);                           //add this intersection point
328     Float_t ckov=pTrk->GetHMPIDsignal();                                        //get ckov angle stored for this track  
329     if(ckov>0){
330       rec.SetTrack(xRa,yRa,thRa,phRa);
331       fRenRin[ch]->SetUniqueID(npTrk);
332       for(Int_t j=0;j<500;j++){
333         TVector2 pos; pos=rec.TracePhot(ckov,j*6.28/500.);
334        if(!AliHMPIDParam::IsInDead(pos.X(),pos.Y())) fRenRin[ch]->SetNextPoint(pos.X(),pos.Y());
335       }      
336     }//if ckov is valid
337   }//tracks loop  
338 }//RenEsd()
339 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
340 void SimulateEsd(AliESDEvent *pEsd)
341 {
342   TParticle part; TLorentzVector mom;
343   for(Int_t iTrk=0;iTrk<100;iTrk++){//stack loop
344     Printf("SimulateEsd: kProton %i",kProton);
345     part.SetPdgCode(kProton);
346     part.SetProductionVertex(0,0,0,0);
347     Double_t eta= -0.2+gRandom->Rndm()*0.4; //rapidity is random [-0.2,+0.2]
348     Double_t phi= gRandom->Rndm()*60.*TMath::DegToRad();   //phi is random      [ 0  , 60 ] degrees    
349     mom.SetPtEtaPhiM(5,eta,phi,part.GetMass());
350     part.SetMomentum(mom);
351     AliESDtrack trk(&part);
352     pEsd->AddTrack(&trk);
353   }//stack loop  
354 }//EsdFromStack()
355 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
356 void SimulateHits(AliESDEvent *pEsd, TClonesArray *pHits)
357 {//used by SimulateEvent to simulate hits out from provided ESD
358   const Int_t kCerenkov=50000050;
359   const Int_t kFeedback=50000051;
360   
361   AliHMPIDRecon rec;
362   Float_t eMip=200e-9,ePho=7.5e-9; 
363   Int_t hc=0; 
364   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop
365     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
366     Float_t xRa,yRa;
367     Int_t ch=AliHMPIDTracker::IntTrkCha(pTrk,xRa,yRa);
368     if(ch<0) continue; //this track does not hit HMPID
369     Float_t beta = pTrk->GetP()/(TMath::Sqrt(pTrk->GetP()*pTrk->GetP()+0.938*0.938));
370     Float_t ckov=TMath::ACos(1./(beta*1.292));
371
372     Float_t theta,phi,xPc,yPc,; pTrk->GetHMPIDtrk(xPc,yPc,theta,phi); rec.SetTrack(xRa,yRa,theta,phi); 
373     
374     if(!AliHMPIDParam::IsInDead(xPc,yPc)) new((*pHits)[hc++]) AliHMPIDHit(ch,eMip,kProton  ,iTrk,xPc,yPc);                 //mip hit
375     Int_t nPhots = (Int_t)(20.*TMath::Power(TMath::Sin(ckov),2)/TMath::Power(TMath::Sin(TMath::ACos(1./1.292)),2));
376     for(int i=0;i<nPhots;i++){
377       TVector2 pos;
378       pos=rec.TracePhot(ckov,gRandom->Rndm()*TMath::TwoPi());
379       if(!AliHMPIDParam::IsInDead(pos.X(),pos.Y())) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kCerenkov,iTrk,pos.X(),pos.Y());
380     }                      //photon hits  
381     for(int i=0;i<3;i++){//feedback photons
382       Float_t x=gRandom->Rndm()*160; Float_t y=gRandom->Rndm()*150;
383       if(!AliHMPIDParam::IsInDead(x,y)) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kFeedback,iTrk,x,y);                 //feedback hits  
384     }//photon hits loop                      
385   }//tracks loop    
386 }//SimulateHits()
387 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
388 TString Stack(Int_t evt,Int_t tid)
389 {
390 // Prints some useful info from stack
391 // Arguments: evt - event number. if not -1 print info only for that event
392 //            tid - track id. if not -1 then print it and all it's mothers if any   
393 //   Returns: mother tid of the given tid if any
394   if(gAL->LoadHeader()) return -1;
395   if(gAL->LoadKinematics()) return -1;
396   
397   gAL->GetEvent(evt);    
398   AliStack *pStack=gAL->Stack();  
399   TParticle *pTrack=pStack->Particle(tid);
400 //  Int_t mtid=pTrack->GetFirstMother();
401   TString str;
402   str.Append(Form("with P = %7.4f (%7.4f,%7.4f,%7.4f) GeV/c ",pTrack->P(),pTrack->Px(),pTrack->Py(),pTrack->Pz()));
403   gAL->UnloadHeader();  gAL->UnloadKinematics();
404   return str;
405 }
406 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
407 void DispTB(TString t0,TString t1,TString t2,TString t3)
408 {
409   fCanvasImp->SetStatusText(t0,0);
410   fCanvasImp->SetStatusText(t1,1);
411   fCanvasImp->SetStatusText(t2,2);
412   fCanvasImp->SetStatusText(t3,3);
413 }
414 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
415 Bool_t FindPos(TVector3 xyz,Int_t &cham,Int_t &sect,Int_t &gap,Int_t &padX,Int_t &padY)
416 {
417 // It finds the chamber parameters for a given point 
418 // Arguments:  xyz  coordinates of a given point in MARS
419 //   Returns:  cham   chamber ID ( 0-6 )
420 //             sect   sector  ID ( 0-5 )
421 //             padX   X (integer) of the pad in X in the given sector (1-80)
422 //             padY   Y (integer) of the pad in Y in the given sector (1-48)
423   cham = -1;
424   TGeoNode *node = gGeoManager->FindNode(xyz.X(),xyz.Y(),xyz.Z());
425   if(!node) return kFALSE;
426   TGeoVolume *vol = node->GetVolume();
427   if(!vol) return kFALSE;
428   TGeoManager *geo = vol->GetGeoManager();
429   if(!geo) return kFALSE;
430   Int_t i=0;
431   while(geo->GetMother(i)) {
432     TString s = geo->GetMother(i)->GetName();
433 //    if(s.Contains("Hcel_")) padX=(Int_t)(s.Remove(0,5)).Atof();
434 //    if(s.Contains("Hrow_")) padY=(Int_t)(s.Remove(0,5)).Atof();
435 //    if(s.Contains("Hsec_")) sect=(Int_t)(s.Remove(0,5)).Atof();
436 //    if(s.Contains("Hmp_"))  cham=(Int_t)(s.Remove(0,4)).Atof();
437     if(s.Contains("Hcel_")) padX=geo->GetMother(i)->GetNumber();
438     if(s.Contains("Hrow_")) padY=geo->GetMother(i)->GetNumber();
439     if(s.Contains("Hsec_")) sect=geo->GetMother(i)->GetNumber();
440     if(s.Contains("Hmp_"))  cham=geo->GetMother(i)->GetNumber();
441     i++; 
442   }
443   return (cham==-1) ? kFALSE : kTRUE;
444 }
445 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
446 void DisplayInfo(Int_t evt,Int_t px, Int_t py)
447 {
448   if(!gPad) return;
449   TVirtualPad *pPad=gPad->GetSelectedPad();
450   Int_t padN = pPad->GetNumber()-1;
451   if(padN<0 || padN>8) return;
452   Int_t ch,pc,padX,padY;
453   ch = fChamN[padN];
454   if(ch<0) return;
455   TString text0,text1,text2,text3;
456   
457   Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py); 
458   fParam->Lors2Pad(x,y,pc,padX,padY);
459   TVector3 xyz=fParam->Lors2Mars(ch,x,y);
460   
461   if(padX>=0&&padY>=0) {
462     text0.Append(Form("Pad(%i,%i)-LORS(%6.2f,%6.2f)-MARS(%7.2f,%7.2f,%7.2f)",padX,padY,x,y,xyz.X(),xyz.Y(),xyz.Z()));
463     text1.Append(Form("Module %i Sector %i",ch,pc));
464     text2="";
465     text3.Append(Form("Pads = %i - Clusters = %i - Occupancy %5.2f%%",fTotPads[ch],fTotClus[ch],100.*fTotPads[ch]/(144.*160.)));
466   }
467   
468   TObject *obj = fCanvas->GetSelected();
469   TString name = obj->GetName();
470   
471 // Find object DIGIT
472
473   if(name=="") {DispTB(text0,text1,text2,text3);return;} // TLatex not interesting
474     
475   if(name.Contains("TBox")) {
476     TBox *box = (TBox*)obj;
477     if(box->GetUniqueID()==0) return; // just black frame hit!
478     text2.Append(Form("charge = %i ADC",box->GetUniqueID()));
479     DispTB(text0,text1,text2,text3);
480     return;
481   }
482 //
483 // Find object HITs CLUSTER TRACK and RING (based on TPolyMarker and symbol)
484 //
485   TPolyMarker *b = (TPolyMarker*)obj;
486
487   const Int_t big = 9999;
488
489   // check if point is near one of the points
490   Int_t distance = big;
491   Int_t index=0;
492
493   Double_t *xPol = b->GetX();
494   Double_t *yPol = b->GetY();
495
496   for(Int_t i=0;i<b->Size();i++) {
497     Int_t pxp = pPad->XtoAbsPixel(pPad->XtoPad(xPol[i]));
498     Int_t pyp = pPad->YtoAbsPixel(pPad->YtoPad(yPol[i]));
499     Int_t d   = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
500     if (d < distance) {distance = d; index=i;}
501   }
502
503   Int_t symbol = b->GetMarkerStyle();
504
505   //  kOpenTriangleUp  Hit Mip
506   //  kOpenCircle      Hit Ckov
507   //  kOpenDiamond     Hit Feedback
508   //  kStar            Cluster 
509   //  kPlus            Track
510   //  kFullDotSmall    Ring
511   
512   // Case Hit Mip of Hit Ckov or Hit Feedback
513   if(symbol==kHitMip || symbol==kHitCko || symbol==kHitFee) {
514     if(evt!=kButton1Down) return;
515     TString nameHit;
516     Int_t iCko[7]={0,0,0,0,0,0,0};
517     Int_t iFee[7]={0,0,0,0,0,0,0};
518     Int_t iMip[7]={0,0,0,0,0,0,0};
519     Int_t iHit;
520     Int_t chHit;
521     Int_t indHit=0;
522     
523     Int_t indTrk = b->GetUniqueID(); // track index
524     for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){//entries loop
525       fHitTree->GetEntry(iEnt);                              
526       for(iHit=0;iHit<fHitLst->GetEntries();iHit++) {       //hits loop
527         AliHMPIDHit *pHit = (AliHMPIDHit*)fHitLst->At(iHit);
528         chHit = pHit->Ch();
529         Int_t pid = pHit->Pid();
530         if(pid==50000050) iCko[chHit]++;
531         else if(pid==50000051) iFee[chHit]++;
532         else iMip[chHit]++;
533
534         if(symbol==kHitMip && ch==chHit && index==iMip[chHit]-1) {indTrk=iEnt;indHit = iHit;break;}
535         if(symbol==kHitCko && ch==chHit && index==iCko[chHit]-1) {indTrk=iEnt;indHit = iHit;break;}
536         if(symbol==kHitFee && ch==chHit && index==iFee[chHit]-1) {indTrk=iEnt;indHit = iHit;break;}
537       }
538     }
539     fHitTree->GetEntry(indTrk);
540     AliHMPIDHit *pHit = (AliHMPIDHit*)fHitLst->At(indHit);
541
542     Int_t tid = pHit->Tid();
543     
544     TString str = Stack(fEvt,tid);
545      
546     Float_t x=pHit->LorsX(); 
547     Float_t y=pHit->LorsY();
548     Float_t charge = pHit->Q();
549     
550     TVector3 v = fParam->Lors2Mars(pHit->Ch(),pHit->LorsX(),pHit->LorsY());
551     Int_t cham,sect,gap,padX,padY;
552     FindPos(v,cham,sect,gap,padX,padY);
553     
554     if(pHit->Pid()==50000050) {nameHit = " Cherenkov Photon ";}
555     else if(pHit->Pid()==50000051) {nameHit = " Feedback Photon ";}
556     else if(fPdg->GetParticle(pHit->Pid())) nameHit = fPdg->GetParticle(pHit->Pid())->GetName();
557     nameHit.Append(Form(" (tid %i)",tid));
558     text0="";text0.Append(Form("Hit(%5.2f,%6.2f) LORS - Pad Volume(%i,%i) - hit %i/%i",x,y,padX,padY,indHit+1,fHitLst->GetEntries()));
559     text2="";text2.Append(Form("Q = %7.2f ADC",charge));
560     text3="";text3="Particle: "+ nameHit +" "+ str;
561   } // Hit  
562   else if (symbol==kCluster) {
563     TClonesArray *pClusCham=(TClonesArray*)fCluLst->At(ch);         //get clusters list for this chamber
564     AliHMPIDCluster *pClu = (AliHMPIDCluster*)pClusCham->At(index); //get current cluster
565     text0="";text0.Append(Form("CLUSTER: x %6.2f y %6.2f",pClu->X(),pClu->Y()));
566     text2="";text2.Append(Form("charge = %i ADC",(Int_t)pClu->Q()));
567     }
568   else if (symbol==kTrack || symbol==kRing) {
569     if(symbol==kRing) index = b->GetUniqueID();
570     AliESDtrack *pTrk=fEsd->GetTrack(index);
571     Float_t thRa,phRa,xRa,yRa; pTrk->GetHMPIDtrk(xRa,yRa,thRa,phRa);
572     Float_t xPc,yPc; AliHMPIDTracker::IntTrkCha(pTrk,xPc,yPc);
573     text0="";text0.Append(Form("TRACK n.%d: x %6.2f y %6.2f at PC plane",index,xPc,yPc));
574     text2="";text2.Append(Form("p = %7.2f GeV/c",pTrk->GetP()));
575     Float_t ckov=pTrk->GetHMPIDsignal();                             
576     Double_t prob[5];
577     pTrk->GetHMPIDpid(prob);
578     if(ckov>0){
579       Float_t x,y;Int_t q,nacc;   pTrk->GetHMPIDmip(x,y,q,nacc);
580       text3="";text3.Append(Form("Theta Cherenkov %5.3f with %i photons |  e %5.1f%% | u %5.1f%% | K %5.1f%% | pi %5.1f%% | p %5.1f%%",
581           pTrk->GetHMPIDsignal(),nacc,prob[0]*100,prob[1]*100,prob[2]*100,prob[3]*100,prob[4]*100));
582     }      
583   }
584 //Update toolbar status  
585
586   DispTB(text0,text1,text2,text3);  
587 }
588 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
589 void CloseInfo()
590 {
591 }
592 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
593 void DoZoom(Int_t evt, Int_t px, Int_t py, TObject *)
594 {
595 //  Printf(" px %i py%i event %i",px,py,evt);
596   if(evt==kMouseMotion||evt==kButton1Down) DisplayInfo(evt,px,py);
597   if(evt==11)  CloseInfo();
598   if(evt!=5 && evt!=6) return; //5- zoom in 6-zoom out
599   const Int_t minZoom=64;
600   const Int_t maxZoom=2;
601   static Int_t zoom[7]={64,64,64,64,64,64,64}; //zoom level
602   
603  // if(!obj->IsA()->InheritsFrom("TPad")) return;  //current object is not pad
604   TVirtualPad *pPad=gPad->GetSelectedPad();
605   Int_t padN = pPad->GetNumber()-1;
606   if(padN<0 || padN>8) return;
607   
608   Int_t iCh = fChamN[padN];
609   if(iCh < 0) return;
610   
611   if(evt==5&&zoom[iCh]==maxZoom) return; 
612   if(evt==6&&zoom[iCh]==minZoom) return; 
613     
614   Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py); 
615  
616   if(evt==5){ zoom[iCh]=zoom[iCh]/2;     pPad->Range(x-zoom[iCh]*2,y-zoom[iCh]*2,x+zoom[iCh]*2,y+zoom[iCh]*2);} //zoom in
617   if(evt==6){ zoom[iCh]=zoom[iCh]*2;     pPad->Range(x-zoom[iCh]*2,y-zoom[iCh]*2,x+zoom[iCh]*2,y+zoom[iCh]*2);} //zoom out 
618   if(zoom[iCh]==minZoom) pPad->Range(-10,-10,AliHMPIDParam::SizeAllX()+5,AliHMPIDParam::SizeAllY()+5);
619   ((TCanvas *)gTQSender)->SetTitle(Form("zoom x%i",minZoom/zoom[iCh]));
620   pPad->Modified();
621   pPad->Update();                                              
622 }
623 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
624 void ReadEvent()
625 {//used by NextEvent()  to read curent event and construct all render elements
626   if(fNevt && fEvt>=fNevt) fEvt=0; //loop over max event
627   if(fNevt && fEvt<0) fEvt=fNevt-1; //loop over max event
628   
629   if(fHitFile){ 
630     if(fHitTree) delete fHitTree;   fHitTree=(TTree*)fHitFile->Get(Form("Event%i/TreeH",fEvt));
631     
632                           fHitTree->SetBranchAddress("HMPID",&fHitLst);
633     for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){    
634                           fHitTree->GetEntry(iEnt);
635       RenderHit(fHitLst);   
636     }//prim loop
637   }//if hits file
638
639   if(fDigFile){ 
640     if(fDigTree) delete fDigTree;   fDigTree=(TTree*)fDigFile->Get(Form("Event%i/TreeD",fEvt));
641     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) fDigTree->SetBranchAddress(Form("HMPID%i",iCh),&(*fDigLst)[iCh]);      
642     fDigTree->GetEntry(0);                              
643     RenderDig(fDigLst);   
644   }//if digs file
645     
646   if(fCluFile){ 
647     if(fCluTree) delete fCluTree; fCluTree=(TTree*)fCluFile->Get(Form("Event%i/TreeR",fEvt));
648     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) fCluTree->SetBranchAddress(Form("HMPID%i",iCh),&(*fCluLst)[iCh]);      
649     fCluTree->GetEntry(0);
650     RenderClu(fCluLst);   
651   }//if clus file
652   
653   if(fEsdFile){//if ESD file open
654     fEsdTree->GetEntry(fEvt);  
655     RenderEsd(fEsd);
656   }//if ESD file 
657 }//ReadEvent()
658 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
659 void SimulateEvent()
660 {// used by NextEvent() to simulate all info
661   AliCDBManager* pCDB = AliCDBManager::Instance();  pCDB->SetDefaultStorage("local://$HOME"); pCDB->SetRun(0);
662   AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean");
663
664   SimulateEsd(fEsd);
665   SimulateHits(fEsd,fHitLst);
666                  AliHMPIDv2::Hit2Sdi(fHitLst,fSdiLst);                               
667           AliHMPIDDigitizer::Sdi2Dig(fSdiLst,fDigLst);     
668       AliHMPIDReconstructor::Dig2Clu(fDigLst,fCluLst);
669             AliHMPIDTracker::Recon(fEsd,fCluLst,(TObjArray*)pNmeanEnt->GetObject());
670             
671   RenderHit(fHitLst);
672   RenderClu(fCluLst);
673   RenderEsd(fEsd);            
674 }//SimulateEvent()
675 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
676 void CheckStatus()
677 {
678   if(fHitFile){
679     if(fStHitMip=="OFF") {fHitMipBok->SetFillColor(kRed);} else {fHitMipBok->SetFillColor(18);}
680     if(fStHitCko=="OFF") {fHitCkoBok->SetFillColor(kRed);} else {fHitCkoBok->SetFillColor(18);}
681     if(fStHitFee=="OFF") {fHitFeeBok->SetFillColor(kRed);} else {fHitFeeBok->SetFillColor(18);}
682   }
683   if(fDigFile){if(fStDig   =="OFF") {fDigBok->SetFillColor(kRed);} else {fDigBok->SetFillColor(18);}}
684   if(fCluFile){if(fStClu   =="OFF") {fCluBok->SetFillColor(kRed);} else {fCluBok->SetFillColor(18);}}
685   if(fEsdFile){if(fStEsd   =="OFF") {fEsdBok->SetFillColor(kRed);} else {fEsdBok->SetFillColor(18);}}
686 }
687 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
688 void GetEvent()
689 {
690   ClearRenders();
691   CheckStatus();
692   switch(fType){
693     case 1: ReadEvent();break;
694 //    case 2: ReadCosmic(); break;
695     case 3: SimulateEvent();     break;
696     default:                 return;
697   }  
698   Draw();
699 }
700 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
701 void NextEvent()
702 {
703   fTimeArrow = 1;
704   fEvt+=fTimeArrow;
705   GetEvent();
706 }
707 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
708 void PrevEvent()
709 {
710   fTimeArrow = -1;
711   fEvt+=fTimeArrow;
712   GetEvent();
713 }
714 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
715 void End()
716 {
717   gSystem->Exit(0);
718 }
719 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
720 void Hdisp()                                  
721 {//display events from files if any in current directory or simulated events
722   
723   fParam=AliHMPIDParam::Instance();                          // first invocation of AliHMPIDParam to initialize geometry...
724   fPdg = TDatabasePDG::Instance();                           // first invocation of TDatabasePDG to retrieve particle infos...
725   
726   CreateContainers();
727   CreateRenders();
728   
729   TString title="Session with";
730   if(gSystem->IsFileInIncludePath("HMPID.Hits.root")){// tries to open hits
731     fHitFile=TFile::Open("HMPID.Hits.root");       fNevt=fHitFile->GetNkeys(); fType=1; title+=Form(" HITS-%i ",fNevt);
732   }
733   
734   if(gSystem->IsFileInIncludePath("HMPID.Digits.root")){// tries to open clusters
735      fDigFile=TFile::Open("HMPID.Digits.root");  fNevt=fDigFile->GetNkeys(); fType=1;  title+=Form(" DIGITS-%i ",fNevt);
736    }
737   
738   if(gSystem->IsFileInIncludePath("HMPID.RecPoints.root")){// tries to open clusters
739     fCluFile=TFile::Open("HMPID.RecPoints.root");  fNevt=fCluFile->GetNkeys(); fType=1;  title+=Form(" CLUSTERS-%i ",fNevt);
740    }
741   
742   if(gSystem->IsFileInIncludePath("AliESDs.root")){     
743     fEsdFile=TFile::Open("AliESDs.root");
744     if(fEsdFile) { 
745       fEsdTree=(TTree*)fEsdFile->Get("esdTree"); 
746       fEsd->ReadFromTree(fEsdTree); fEsd->GetStdContent();                       //clm: new ESD schema: see Task Force meeting 20th June, 2007
747       if(fEsdTree) fNevt=fEsdTree->GetEntries(); fType=1;  title+=Form(" ESD-%i ",fNevt); else {delete fEsdFile;fEsdFile=0x0;}
748     } else {delete fEsdFile; delete fEsdTree;}
749     //clm: we need to set the magnetic field  
750     if(gSystem->IsFileInIncludePath("galice.root")){
751     if(gAlice) delete gAlice;                       
752     gAL=AliRunLoader::Open();                       
753     gAL->LoadgAlice();           
754 //    if(gAL)AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE);                         
755    }else{
756           Printf("=============== NO galice file! Magnetic field for ESD tracking is: %f ===============",AliTracker::GetBz());
757      }  
758       
759   }
760
761   if(gSystem->IsFileInIncludePath("cosmic.root")){                          //clm: Check if cosmic file is in the folder
762     fCosFile=TFile::Open("cosmic.root");  fCosTree=(TTree*)fCosFile->Get("cosmic");  fNevt=fCosTree->GetEntries();  fType=2;
763     fCosTree->SetBranchAddress("Digs",&fDigLst);   fCosTree->SetBranchAddress("Clus",&fCluLst); 
764   }            
765
766   fCanvas=new TCanvas("all","",-1024,768);fCanvas->SetWindowSize(1024,768);
767   fCanvasImp = fCanvas->GetCanvasImp();
768   fCanvasImp->ShowStatusBar();
769   fCanvas->Divide(3,3,0,0);//  pAll->ToggleEditor();
770   fCanvas->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",0,0,"DoZoom(Int_t,Int_t,Int_t,TObject*)");
771   fCanvas->cd(7); 
772   gStyle->SetPalette(1,0);
773   switch(fType){
774     case 1: fCanvas->SetTitle(title.Data());                      
775                          TButton *pNxtBtn = new TButton("Next"      ,"NextEvent()"   ,0.0,0.0,0.3,0.1); pNxtBtn->Draw(); 
776                          TButton *pPrvBtn = new TButton("Previous"  ,"PrevEvent()"   ,0.3,0.0,0.6,0.1); pPrvBtn->Draw(); 
777             if(fHitFile){TButton *pHitBtn = new TButton("Print hits","PrintHits()"   ,0.0,0.2,0.3,0.3); pHitBtn->Draw();}  
778             if(fDigFile){TButton *pDigBtn = new TButton("Print digs","PrintDigs()"   ,0.0,0.4,0.3,0.5); pDigBtn->Draw();}  
779             if(fCluFile){TButton *pCluBtn = new TButton("Print clus","PrintClus()"   ,0.0,0.6,0.3,0.7); pCluBtn->Draw();} 
780             if(fEsdFile){TButton *pEsdBtn = new TButton("Print ESD ","PrintEsd()"    ,0.0,0.8,0.3,0.9); pEsdBtn->Draw();}
781             if(fHitFile){      fHitMipBok = new TButton(Form("Mip  %s",fStHitMip.Data()),"SwitchHitMip()",0.3,0.16,0.6,0.22); fHitMipBok->Draw();}
782             if(fHitFile){      fHitCkoBok = new TButton(Form("Ckov %s",fStHitCko.Data()),"SwitchHitCko()",0.3,0.22,0.6,0.28); fHitCkoBok->Draw();}  
783             if(fHitFile){      fHitFeeBok = new TButton(Form("Fdbk %s",fStHitFee.Data()),"SwitchHitFee()",0.3,0.28,0.6,0.34); fHitFeeBok->Draw();}  
784             if(fDigFile){         fDigBok = new TButton(fStDig      ,"SwitchDigs()"  ,0.3,0.4,0.6,0.5); fDigBok->Draw();}  
785             if(fCluFile){         fCluBok = new TButton(fStClu      ,"SwitchClus()"  ,0.3,0.6,0.6,0.7); fCluBok->Draw();} 
786             if(fEsdFile){         fEsdBok = new TButton(fStEsd      ,"SwitchEsd()"   ,0.3,0.8,0.6,0.9); fEsdBok->Draw();}
787
788             if(fHitFile){TButton *fHitMipOnly = new TButton(" Only ","OnlyHitMip()",0.6,0.16,0.9,0.22); fHitMipOnly->Draw();}  
789             if(fHitFile){TButton *fHitCkoOnly = new TButton(" Only ","OnlyHitCko()",0.6,0.22,0.9,0.28); fHitCkoOnly->Draw();}  
790             if(fHitFile){TButton *fHitFeeOnly = new TButton(" Only ","OnlyHitFee()",0.6,0.28,0.9,0.34); fHitFeeOnly->Draw();}  
791             if(fDigFile){TButton *fDigOnly    = new TButton(" Only ","OnlyDigs()"  ,0.6,0.40,0.9,0.50); fDigOnly->Draw();}  
792             if(fCluFile){TButton *fCluOnly    = new TButton(" Only ","OnlyClus()"  ,0.6,0.60,0.9,0.70); fCluOnly->Draw();} 
793             if(fEsdFile){TButton *fEsdOnly    = new TButton(" Only ","OnlyEsd()"   ,0.6,0.80,0.9,0.90); fEsdOnly->Draw();}
794                          TButton *pAll        = new TButton(" ALL  ","All()"       ,0.5,0.90,0.7,1.00); pAll->Draw();
795             
796                                      TButton *pEnd    = new TButton("Quit"      ,"End()"         ,0.6,0.0,0.9,0.1); pEnd->Draw();
797     break;
798     case 2: fCanvas->SetTitle("COSMIC");
799                          TButton *pCosBtn = new TButton("Next"      ,"NextEvent()",0.0,0.0,0.3,0.1); pCosBtn->Draw(); 
800              
801     break;
802     case 3: fCanvas->SetTitle("SIMULATION"); 
803             TButton *pSimBtn = new TButton("Simulate"  ,"NextEvent()",0.0,0.0,0.3,0.1); pSimBtn->Draw(); 
804             TButton *pHitBtn = new TButton("Print hits","PrintHits()",0.0,0.2,0.3,0.3); pHitBtn->Draw();  
805             TButton *pDigBtn = new TButton("Print digs","PrintDigs()",0.0,0.4,0.3,0.5); pDigBtn->Draw();
806             TButton *pCluBtn = new TButton("Print clus","PrintClus()",0.0,0.6,0.3,0.7); pCluBtn->Draw(); 
807             TButton *pEsdBtn = new TButton("Print ESD" ,"PrintEsd()" ,0.0,0.8,0.3,0.9); pEsdBtn->Draw();
808             if(fHitFile){      fHitMipBok = new TButton(Form("Mips %s",fStHitMip.Data()),"SwitchHitMip()",0.3,0.16,0.6,0.22); fHitMipBok->Draw();}  
809             if(fHitFile){      fHitCkoBok = new TButton(Form("Ckov %s",fStHitCko.Data()),"SwitchHitCko()",0.3,0.22,0.6,0.28); fHitCkoBok->Draw();}  
810             if(fHitFile){      fHitFeeBok = new TButton(Form("Fdbk %s",fStHitFee.Data()),"SwitchHitFee()",0.3,0.28,0.6,0.34); fHitFeeBok->Draw();}  
811             if(fDigFile){         fDigBok = new TButton(fStDig      ,"SwitchDigs()"  ,0.3,0.4,0.6,0.5); fDigBok->Draw();}  
812             if(fCluFile){         fCluBok = new TButton(fStClu      ,"SwitchClus()"  ,0.3,0.6,0.6,0.7); fCluBok->Draw();} 
813             if(fEsdFile){         fEsdBok = new TButton(fStEsd      ,"SwitchEsd()"   ,0.3,0.8,0.6,0.9); fEsdBok->Draw();}
814             
815             if(fHitFile){TButton *fHitMipOnly = new TButton(" Only ","OnlyHitMip()",0.6,0.16,0.9,0.22); fHitMipOnly->Draw();}  
816             if(fHitFile){TButton *fHitCkoOnly = new TButton(" Only ","OnlyHitCko()",0.6,0.22,0.9,0.28); fHitCkoOnly->Draw();}  
817             if(fHitFile){TButton *fHitFeeOnly = new TButton(" Only ","OnlyHitFee()",0.6,0.28,0.9,0.34); fHitFeeOnly->Draw();}  
818             if(fDigFile){TButton *fDigOnly    = new TButton(" Only ","OnlyDigs()"  ,0.6,0.40,0.9,0.50); fDigOnly->Draw();}  
819             if(fCluFile){TButton *fCluOnly    = new TButton(" Only ","OnlyClus()"  ,0.6,0.60,0.9,0.70); fCluOnly->Draw();} 
820             if(fEsdFile){TButton *fEsdOnly    = new TButton(" Only ","OnlyEsd()"   ,0.6,0.80,0.9,0.90); fEsdOnly->Draw();}
821                          TButton *pSimAll     = new TButton(" ALL  ","All()"       ,0.5,0.90,0.7,1.00); pSimAll->Draw();
822     break;
823   }
824     
825   NextEvent();           
826       
827 }//Hdisp()
828 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
829 void SwitchHitMip()
830 {
831   fStHitMip=(fStHitMip=="OFF")? "ON":"OFF";
832   GetEvent();
833 }
834 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
835 void SwitchHitCko()
836 {
837   fStHitCko=(fStHitCko=="OFF")? "ON":"OFF";
838   GetEvent();
839 }
840 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
841 void SwitchHitFee()
842 {
843   fStHitFee=(fStHitFee=="OFF")? "ON":"OFF";
844   GetEvent();
845 }
846 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
847 void SwitchDigs()
848 {
849   fStDig=(fStDig=="OFF")? "ON":"OFF";
850   GetEvent();
851 }
852 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
853 void SwitchClus()
854 {
855   fStClu=(fStClu=="OFF")? "ON":"OFF";
856   GetEvent();
857 }
858 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
859 void SwitchEsd()
860 {
861   fStEsd=(fStEsd=="OFF")? "ON":"OFF";
862    GetEvent();
863 }
864 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
865 void AllNo()
866 {
867   if(fHitFile) {
868     fStHitMip = "OFF";
869     fStHitCko = "OFF";
870     fStHitFee = "OFF";
871   }
872   if(fDigFile) fStDig    = "OFF";
873   if(fCluFile) fStClu    = "OFF";
874   if(fEsdFile) fStEsd    = "OFF";
875 }
876 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
877 void AllYes()
878 {
879   if(fHitFile) {
880     fStHitMip = "ON";
881     fStHitCko = "ON";
882     fStHitFee = "ON";
883   }
884   if(fDigFile) fStDig    = "ON";
885   if(fCluFile) fStClu    = "ON";
886   if(fEsdFile) fStEsd    = "ON";
887 }
888 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
889 void All()
890 {
891   AllYes();
892   GetEvent();
893 }
894 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
895 void OnlyHitMip()
896 {
897   AllNo();
898   fStHitMip = "ON";
899   GetEvent();
900 }
901 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
902 void OnlyHitCko()
903 {
904   AllNo();
905   fStHitCko = "ON";
906   GetEvent();
907 }
908 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
909 void OnlyHitFee()
910 {
911   AllNo();
912   fStHitFee = "ON";
913   GetEvent();
914 }
915 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
916 void OnlyDigs()
917 {
918   AllNo();
919   fStDig    = "ON";
920   GetEvent();
921 }
922 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
923 void OnlyClus()
924 {
925   AllNo();
926   fStClu    = "ON";
927   GetEvent();
928 }
929 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
930 void OnlyEsd()
931 {
932   AllNo();
933   fStEsd    = "ON";
934   GetEvent();
935 }
936 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++