]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hdisp.C
Updated to display cosmic runs (L. Molnar)
[u/mrichter/AliRoot.git] / HMPID / Hdisp.C
1 TCanvas *pAll=0;
2 AliRunLoader *gAL=0; AliLoader *gHL=0; AliESD *gEsd=0; TTree *gEsdTr=0; AliHMPID *gH=0;
3 Int_t gEvt=0; Int_t gMaxEvt=0;
4 TObjArray *pNmean;
5
6 TChain *pCosCh=new TChain("cosmic");                                                               //clm: Define TChain for cosmic
7 TObjArray *pCosDigAll=0;                                                                           //clm: Define global Digits
8 TObjArray *pCosCluAll=0;                                                                           //clm: Define global Clusters
9 Int_t gCosRun=0;                                                                                   //clm: global cosmic event number
10 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
11 void Hdisp(Int_t cosRun=44)                                                    //clm: Select cosmic file for display
12 {//display events from files if any in current directory or simulated events
13   pAll=new TCanvas("all","",1300,900); pAll->Divide(3,3,0,0);
14 //  pAll->ToggleEditor();
15   pAll->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",0,"","DoZoom(Int_t,Int_t,Int_t,TObject*)");
16
17   OpenCalib();  
18      if(gSystem->IsFileInIncludePath("galice.root")){// tries to open session
19       if(gAlice) delete gAlice;                                               //in case we execute this in aliroot delete default AliRun object 
20       gAL=AliRunLoader::Open();                                                                    //try to open galice.root from current dir 
21       gAL->LoadgAlice();                                                                           //take new AliRun object from galice.root   
22       gHL=gAL->GetDetectorLoader("HMPID");  gH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");  //get HMPID object from galice.root
23       gMaxEvt=gAL->GetNumberOfEvents()-1;
24       gHL->LoadHits(); gHL->LoadSDigits(); gHL->LoadDigits(); gHL->LoadRecPoints();
25
26       TFile *pEsdFl=TFile::Open("AliESDs.root"); gEsdTr=(TTree*) pEsdFl->Get("esdTree"); gEsdTr->SetBranchAddress("ESD", &gEsd);
27       pAll->cd(7); TButton *pBtn=new TButton("Next","ReadEvt()",0,0,0.2,0.1);   pBtn->Draw();
28                    TButton *pHitBtn=new TButton("Print hits","PrintHits()",0,0.2,0.3,0.3);   pHitBtn->Draw();
29                    TButton *pSdiBtn=new TButton("Print sdis","PrintSdis()",0,0.4,0.3,0.5);   pSdiBtn->Draw();
30                    TButton *pDigBtn=new TButton("Print digs","PrintDigs()",0,0.6,0.3,0.7);   pDigBtn->Draw();
31                    TButton *pCluBtn=new TButton("Print clus","PrintClus()",0,0.8,0.3,0.9);   pCluBtn->Draw();  
32       ReadEvt();
33     }
34     else if ( gSystem->IsFileInIncludePath(Form("cosmic%d.root",cosRun))){                          //clm: Check if cosmic file is in the folder
35     gCosRun=cosRun;
36     pCosCh->Add(Form("cosmic%d.root",gCosRun));                                                      //clm: Add cosmic file to chain
37     pCosCh->SetBranchAddress("Digs",&pCosDigAll);                                                   //clm: Set digit branch address
38     pCosCh->SetBranchAddress("Clus",&pCosCluAll);                                                   //clm: Set cluster branch address    
39     gMaxEvt=pCosCh->GetEntries()-1;                                                                 //clm: Get number of events from the cosmic chain
40     pAll->cd(7); TButton *pCosBtn=new TButton("Next Cosmic","ReadCosEvt()",0,0,0.3,0.1);   pCosBtn->Draw();   //clm: define next button
41     ReadCosEvt();                                                                                   //clm: Read first cosmic event  
42     }            
43     else{
44       pAll->cd(7); TButton *pBtn=new TButton("Next","SimEvt()",0,0,0.2,0.1);   pBtn->Draw(); 
45       SimEvt();
46           }      
47 }
48
49 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
50 void ReadCosEvt()
51 {// Read curent cosmic event and display it assumes that session is alredy opened
52   if(gEvt>gMaxEvt) gEvt=0; if(gEvt<0) gEvt=gMaxEvt;                                     //clm: set event limits
53   pCosCh->GetEntry(gEvt);                                                               //clm: read event from chain
54   DrawCosEvt(pCosDigAll,pCosCluAll);                                                    //clm: draw cosmic event
55   gEvt++;                                                                             
56 }
57 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
58 void ReadEvt()
59 {// Read curent event and display it assumes that session is alredy opened
60   TClonesArray hits("AliHMPIDHit");
61   if(gEvt>gMaxEvt) gEvt=0; if(gEvt<0) gEvt=gMaxEvt;
62   
63   gEsdTr->GetEntry(gEvt); gAL->GetEvent(gEvt); 
64   ReadHits(&hits); gHL->TreeS()->GetEntry(0); gHL->TreeD()->GetEntry(0); gHL->TreeR()->GetEntry(0);
65     
66   DrawEvt(&hits,gH->DigLst(),gH->CluLst(),gEsd);
67   gEvt++;
68 }
69 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
70 void SimEvt()
71 {
72   TClonesArray hits("AliHMPIDHit"); 
73   TClonesArray sdig("AliHMPIDDigit"); 
74   TObjArray    digs(7); for(Int_t i=0;i<7;i++) digs.AddAt(new TClonesArray("AliHMPIDDigit"),i);
75   TObjArray    clus(7); for(Int_t i=0;i<7;i++) clus.AddAt(new TClonesArray("AliHMPIDCluster"),i);
76   AliESD esd;
77   AliHMPIDDigit::fSigmas=4;
78   AliHMPIDDigitizer::DoNoise(kFALSE);
79   SimEsd(&esd);
80   SimHits(&esd,&hits);
81              AliHMPIDv1::Hit2Sdi(&hits,&sdig);                               
82       AliHMPIDDigitizer::Sdi2Dig(&sdig,&digs);     
83       AliHMPIDReconstructor::Dig2Clu(&digs,&clus);
84       AliHMPIDTracker::Recon(&esd,&clus,pNmean);
85   
86   pAll->cd(3);  gPad->Clear(); TLatex txt;txt.DrawLatex(0.2,0.2,Form("Simulated event %i",gEvt));
87   DrawEvt(&hits,&digs,&clus,&esd);  
88   gEvt++;
89 }//SimEvt()
90 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
91 void SimEsd(AliESD *pEsd)
92 {
93   TParticle part; TLorentzVector mom;
94   for(Int_t iTrk=0;iTrk<100;iTrk++){//stack loop
95     part.SetPdgCode(kProton);
96     part.SetProductionVertex(0,0,0,0);
97     Double_t eta= -0.2+gRandom->Rndm()*0.4; //rapidity is random [-0.2,+0.2]
98     Double_t phi= gRandom->Rndm()*60.*TMath::DegToRad();   //phi is random      [ 0  , 60 ] degrees    
99     mom.SetPtEtaPhiM(5,eta,phi,part.GetMass());
100     part.SetMomentum(mom);
101     AliESDtrack trk(&part);
102     pEsd->AddTrack(&trk);
103   }//stack loop  
104 }//EsdFromStack()
105 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
106 void SimHits(AliESD *pEsd, TClonesArray *pHits)
107 {
108   const Int_t kCerenkov=50000050;
109   const Int_t kFeedback=50000051;
110   
111   AliHMPIDRecon rec;
112   
113   Int_t hc=0; 
114   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop
115     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
116     Float_t xRa,yRa;
117     Int_t ch=AliHMPIDTracker::IntTrkCha(pTrk,xRa,yRa);
118     if(ch<0) continue; //this track does not hit HMPID
119     Float_t beta = pTrk->GetP()/(TMath::Sqrt(pTrk->GetP()*pTrk->GetP()+0.938*0.938));
120     Float_t ckov=TMath::ACos(1./(beta*1.292));
121
122     Float_t theta,phi,xPc,yPc,; 
123     pTrk->GetHMPIDtrk(xPc,yPc,theta,phi); 
124     rec.SetTrack(xRa,yRa,theta,phi); 
125     
126     if(!AliHMPIDDigit::IsInDead(xPc,yPc)) new((*pHits)[hc++]) AliHMPIDHit(ch,200e-9,kProton  ,iTrk,xPc,yPc);                 //mip hit
127     Int_t nPhots = (Int_t)(20.*TMath::Power(TMath::Sin(ckov),2)/TMath::Power(TMath::Sin(TMath::ACos(1./1.292)),2));
128     for(int i=0;i<nPhots;i++){
129       TVector2 pos;
130       pos=rec.TracePhot(ckov,gRandom->Rndm()*TMath::TwoPi());
131       if(!AliHMPIDDigit::IsInDead(pos.X(),pos.Y())) new((*pHits)[hc++]) AliHMPIDHit(ch,7.5e-9,kCerenkov,iTrk,pos.X(),pos.Y());
132     }                      //photon hits  
133   }//tracks loop    
134 }//SimHits()
135 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
136 void ReadHits(TClonesArray *pHitLst)
137 {
138   pHitLst->Delete();  Int_t cnt=0;
139   for(Int_t iEnt=0;iEnt<gHL->TreeH()->GetEntries();iEnt++){       //TreeH loop
140     gHL->TreeH()->GetEntry(iEnt);                                 //get current entry (prim)                
141     for(Int_t iHit=0;iHit<gH->Hits()->GetEntries();iHit++){       //hits loop
142       AliHMPIDHit *pHit = (AliHMPIDHit*)gH->Hits()->At(iHit);     //get current hit        
143       new((*pHitLst)[cnt++]) AliHMPIDHit(*pHit);
144     }//hits loop for this entry
145   }//tree entries loop
146   
147 }//ReadHits()
148 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
149 void DrawCh(Int_t iCh) 
150
151   gPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5); 
152   if(iCh>=0){TLatex txt; txt.SetTextSize(0.1); txt.DrawLatex(-5,-5,Form("%i",iCh));}
153   
154   for(Int_t iPc=AliHMPIDDigit::kMinPc;iPc<=AliHMPIDDigit::kMaxPc;iPc++){
155     TBox *pBox=new TBox(AliHMPIDDigit::fMinPcX[iPc],AliHMPIDDigit::fMinPcY[iPc],
156                         AliHMPIDDigit::fMaxPcX[iPc],AliHMPIDDigit::fMaxPcY[iPc]);
157     
158     pBox->SetFillStyle(0);
159     pBox->Draw();
160     
161   }//PC loop      
162 }//DrawCh()
163 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 void DrawEvt(TClonesArray *pHitLst,TObjArray *pDigLst,TObjArray *pCluLst,AliESD *pEsd)
165 {//draws all the objects of current event in given canvas
166   AliHMPIDRecon rec;  
167   TPolyMarker *pTxC[7], *pRin[7]; TMarker *pMip,*pCko,*pFee,*pDig,*pClu;
168   pMip=new TMarker; pMip->SetMarkerColor(kRed);  pMip->SetMarkerStyle(kOpenTriangleUp);
169   pCko=new TMarker; pCko->SetMarkerColor(kRed);  pCko->SetMarkerStyle(kOpenCircle);
170   pFee=new TMarker; pFee->SetMarkerColor(kRed);  pFee->SetMarkerStyle(kOpenDiamond);
171   pDig=new TMarker; pDig->SetMarkerColor(kGreen);pDig->SetMarkerStyle(kOpenSquare);
172   pClu=new TMarker; pClu->SetMarkerColor(kBlue); pClu->SetMarkerStyle(kStar);
173   for(Int_t ch=0;ch<7;ch++){
174     pTxC[ch]=new TPolyMarker; pTxC[ch]->SetMarkerStyle(2); pTxC[ch]->SetMarkerColor(kRed); pTxC[ch]->SetMarkerSize(3);
175     pRin[ch]=new TPolyMarker; pRin[ch]->SetMarkerStyle(6); pRin[ch]->SetMarkerColor(kMagenta);
176   }
177   
178   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop to collect cerenkov rings and intersection points
179     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
180     Int_t ch=pTrk->GetHMPIDcluIdx();
181     if(ch<0) continue; //this track does not hit HMPID
182     ch/=1000000; 
183     Float_t th,ph,xRad,yRad; pTrk->GetHMPIDtrk(xRad,yRad,th,ph);//get info on current track
184 //    pTxC[ch]->SetNextPoint(xPc,yPc);                          //add new intersection point  TEMPORARLY DISABLED...no more available in ESD!
185     
186     Float_t ckov=pTrk->GetHMPIDsignal();  Float_t err=TMath::Sqrt(pTrk->GetHMPIDchi2());
187     if(ckov>0){
188       Printf("theta %f phi %f ckov %f",th*TMath::RadToDeg(),ph*TMath::RadToDeg(),ckov);
189       rec.SetTrack(xRad,yRad,th,ph);
190       for(Int_t j=0;j<100;j++){ 
191         TVector2 pos;
192         pos=rec.TracePhot(ckov,j*0.0628);
193        if(!AliHMPIDDigit::IsInDead(pos.X(),pos.Y())) pRin[ch]->SetNextPoint(pos.X(),pos.Y());
194       }      
195     }
196   }//tracks loop
197       
198   Int_t totHit=pHitLst->GetEntriesFast(),totDig=0,totClu=0,totTxC=0;
199   Int_t totCkov=0, totFeed=0,totMip=0;
200   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop    
201     totTxC+=pTxC[iCh]->GetN();
202     totDig+=((TClonesArray*)pDigLst->At(iCh))->GetEntriesFast();
203     totClu+=((TClonesArray*)pCluLst->At(iCh))->GetEntriesFast();
204     
205     switch(iCh){
206       case 6: pAll->cd(1); break; case 5: pAll->cd(2); break;
207       case 4: pAll->cd(4); break; case 3: pAll->cd(5); break; case 2: pAll->cd(6); break;
208                                   case 1: pAll->cd(8); break; case 0: pAll->cd(9); break;
209     }
210     gPad->SetEditable(kTRUE); gPad->Clear(); 
211     DrawCh(iCh);
212                            
213     ((TClonesArray*)pDigLst->At(iCh))->Draw();               //draw digits
214     
215     totCkov=0;totFeed=0;totMip=0;
216     for(Int_t iHit=0;iHit<pHitLst->GetEntriesFast();iHit++){ // Draw hits
217       AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);
218       if(pHit->Ch()==iCh) pHit->Draw();
219       if(pHit->Pid()==50000050) totCkov++;
220       else if(pHit->Pid()==50000051) totFeed++;
221       else totMip++;
222     }    
223            
224     ((TClonesArray*)pCluLst->At(iCh))->Draw();              //draw clusters
225                             pTxC[iCh]->Draw();              //draw intersections
226                             pRin[iCh]->Draw("CLP");         //draw rings
227     gPad->SetEditable(kFALSE);
228   }//chambers loop
229   
230   pAll->cd(3);  gPad->Clear(); TLegend *pLeg=new TLegend(0.2,0.2,0.8,0.8);
231                                         pLeg->SetHeader(Form("Event %i Total %i",gEvt,gMaxEvt+1));
232                                         pLeg->AddEntry(pTxC[0],Form("TRKxPC %i"     ,totTxC),"p");
233                                         pLeg->AddEntry(pMip   ,Form("Mip hits %i"   ,totMip),"p");    
234                                         pLeg->AddEntry(pCko   ,Form("Ckov hits %i"  ,totCkov),"p");    
235                                         pLeg->AddEntry(pFee   ,Form("Feed hits %i"  ,totFeed),"p");    
236                                         pLeg->AddEntry(pDig   ,Form("Digs %i"       ,totDig),"p");    
237                                         pLeg->AddEntry(pClu   ,Form("Clus %i"       ,totClu),"p");    
238                                         pLeg->Draw();
239 }//DrawEvt()
240 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
241 void DrawCosEvt(TObjArray *pCosDigLst,TObjArray *pCosCluLst)                              //clm: new method to read and display cosmic events 
242 {//draws all the objects of current event in given canvas
243   TMarker *pDig,*pClu;
244   pDig=new TMarker; pDig->SetMarkerColor(kGreen);pDig->SetMarkerStyle(kOpenSquare);
245   pClu=new TMarker; pClu->SetMarkerColor(kBlue); pClu->SetMarkerStyle(kStar);
246   
247       
248   Int_t totDig=0,totClu=0;
249   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop                                            //clm: chambers loop is now not needed since iCh=0 but kept for provision
250     totDig+=((TClonesArray*)pCosDigLst->At(iCh))->GetEntriesFast();
251     totClu+=((TClonesArray*)pCosCluLst->At(iCh))->GetEntriesFast();
252     
253     switch(iCh){
254       case 6: pAll->cd(1); break; case 5: pAll->cd(2); break;
255       case 4: pAll->cd(4); break; case 3: pAll->cd(5); break; case 2: pAll->cd(6); break;
256                                   case 1: pAll->cd(8); break; case 0: pAll->cd(9); break;
257     }
258    gPad->SetEditable(kTRUE); gPad->Clear(); 
259    DrawCh(iCh);
260   if(gCosRun < 500 && iCh == 6) {                           //clm: hard coded selection since raw data does not contain ch info which is set to 0
261   ((TClonesArray*)pCosDigLst->At(0))->Draw();               //draw digits
262   ((TClonesArray*)pCosCluLst->At(0))->Draw();              //draw clusters
263   }
264   if (gCosRun >= 500 && iCh == 5) { 
265   ((TClonesArray*)pCosDigLst->At(0))->Draw();               //draw digits
266   ((TClonesArray*)pCosCluLst->At(0))->Draw();              //draw clusters
267   }
268    
269    gPad->SetEditable(kFALSE);
270   }//chambers loop
271   pAll->cd(3);  gPad->Clear(); TLegend *pLeg=new TLegend(0.2,0.2,0.8,0.8);
272                                         pLeg->SetHeader(Form("Cosmic Run %i Event %i Total %i",gCosRun,gEvt,gMaxEvt+1));
273                                         pLeg->AddEntry(pDig   ,Form("Digs %i"       ,totDig),"p");    
274                                         pLeg->AddEntry(pClu   ,Form("Clus %i"       ,totClu),"p");    
275                                         pLeg->Draw();
276 }//DrawCosEvt()
277 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
278 void DoZoom(Int_t evt, Int_t px, Int_t py, TObject *obj)
279 {
280   if(evt!=5 && evt!=6) return; //5- zoom in 6-zoom out
281   const Int_t minZoom=64;
282   const Int_t maxZoom=2;
283   static Int_t zoom=minZoom; //zoom level
284   if(evt==5&&zoom==maxZoom) return; 
285   if(evt==6&&zoom==minZoom) return; 
286   
287  // if(!obj->IsA()->InheritsFrom("TPad")) return;  //current object is not pad
288   TVirtualPad *pPad=gPad->GetSelectedPad();
289   if(pPad->GetNumber()==3 || pPad->GetNumber()==7) return; //current pad is wrong
290
291  // Printf("evt=%i (%i,%i) %s",evt,px,py,obj->GetName());
292     
293   Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py); 
294  
295   if(evt==5){ zoom=zoom/2;     pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom in
296   else      { zoom=zoom*2;     pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom out 
297   if(zoom==minZoom) pPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5);
298   ((TCanvas *)gTQSender)->SetTitle(Form("zoom x%i",minZoom/zoom));
299   pPad->Modified();
300   pPad->Update();                                              
301 }
302 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
303 void PrintHits()
304 {
305 //Prints a list of HMPID hits for a given event. Default is event number 0.
306   Printf("List of HMPID hits for event %i",gEvt);
307   
308   Int_t iTotHits=0;
309   for(Int_t iPrim=0;iPrim<gHL->TreeH()->GetEntries();iPrim++){//prims loop
310     gHL->TreeH()->GetEntry(iPrim);      
311     gH->Hits()->Print();
312     iTotHits+=gH->Hits()->GetEntries();
313   }
314   Printf("totally %i hits for event %i",iTotHits,gEvt);
315 }
316 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
317 void PrintSdis()
318 {
319 //prints a list of HMPID sdigits  for a given event
320   Printf("List of HMPID sdigits for event %i",gEvt);
321   gH->SdiLst()->Print();
322   Printf("totally %i sdigits for event %i",gH->SdiLst()->GetEntries(),gEvt);
323 }
324 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
325 void PrintDigs()
326 {
327 //prints a list of HMPID digits
328   Printf("List of HMPID digits for event %i",gEvt);  
329   gH->DigLst()->Print();
330   Int_t totDigs=0;  for(Int_t i=0;i<7;i++) {totDigs+=gH->DigLst(i)->GetEntries();}
331   Printf("totally %i digits for event %i",totDigs,gEvt);
332 }
333 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
334 void PrintClus()
335 {//prints a list of HMPID clusters  for a given event
336   Printf("List of HMPID clusters for event %i",gEvt);
337   gH->CluLst()->Print();
338   
339   Int_t iCluCnt=0; for(Int_t iCh=0;iCh<7;iCh++) iCluCnt+=gH->CluLst(iCh)->GetEntries();
340   
341   Printf("totally %i clusters for event %i",iCluCnt,gEvt);
342 }
343 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
344 void OpenCalib()
345 {
346   AliCDBManager* pCDB = AliCDBManager::Instance();
347   pCDB->SetDefaultStorage("local://$HOME");
348   AliCDBEntry *pQthreEnt=pCDB->Get("HMPID/Calib/Qthre",0);
349   AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean",0);
350   
351   if(!pQthreEnt || ! pNmeanEnt) return;
352   
353   pNmean=(TObjArray*)pNmeanEnt->GetObject(); 
354 }