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