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