]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHReconstructor.cxx
New trigger descriptor for proton-proton which includes TOF minimum bias
[u/mrichter/AliRoot.git] / RICH / AliRICHReconstructor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 #include "AliRICHReconstructor.h" //class header
17 #include "AliRICH.h"              //Reconstruct(...) 
18 #include <AliRawReader.h>         //Reconstruct(...)
19 #include <AliRun.h>               //ConvertDigits uses gAlice
20 #include <AliESD.h>               //RichAna()
21 #include <AliStack.h>             //RichAna()
22 #include <TFile.h>                //RichAna()
23 #include <TFile.h>                //RichAna()
24 #include <TParticle.h>            //RichAna()
25 #include <TH1F.h>                 //CluQA() 
26 #include <TH2F.h>                 //CluQA() 
27 #include <TCanvas.h>              //CluQA() 
28 #include <TNtupleD.h>             //CheckPR()
29 #include <TPolyMarker.h>          //Test()
30 ClassImp(AliRICHReconstructor)
31
32 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
33 void AliRICHReconstructor::CluQA(AliRunLoader *pAL)
34 {
35 // Quality assesment plots for clusters. 
36 // This methode takes list of digits and form list of clusters again in order to 
37 // calculate cluster shape and cluster particle mixture    
38   AliLoader *pRL=pAL->GetDetectorLoader("RICH");  AliRICH *pRich=(AliRICH*)pAL->GetAliRun()->GetDetector("RICH");//get pointers for RICH and RICH loader
39   Int_t iNevt=pAL->GetNumberOfEvents();  if(iNevt==0)             {AliInfoClass("No events");return;}   
40                                          if(pRL->LoadDigits())    {AliInfoClass("No digits file");return;}
41                                             pAL->LoadHeader();
42                                             pAL->LoadKinematics(); 
43 //                                            AliStack *pStack=pAL->Stack();
44   TH1::AddDirectory(kFALSE);
45   
46         
47   TH1F*    pQ=new TH1F("RichCluQdc"     ,"All QDC;q [QDC]"           ,4000 ,0  ,4000);// QDC hists
48   TH1F* pCerQ=new TH1F("RichCluCerQdc"  ,"Ckov QDC;q [QDC]"          ,4000 ,0  ,4000);
49   TH1F* pMipQ=new TH1F("RichCluMipQdc"  ,"MIP QDC;q [QDC]"           ,4000 ,0  ,4000);
50   
51   TH1F*    pS=new TH1F("RichCluSize"    ,"Cluster size;size"         ,100  ,0  ,100 );// size hists
52   TH1F* pCerS=new TH1F("RichCluCerSize" ,"Ckov size;size"            ,100  ,0  ,100 );
53   TH1F* pMipS=new TH1F("RichCluMipSize" ,"MIP size;size"             ,100  ,0  ,100 );
54   
55   TH2F*    pM=new TH2F("RichCluMap"     ,"Cluster map;x [cm];y [cm]" ,1000 ,0  ,AliRICHParam::PcSizeX(),1000,0,AliRICHParam::PcSizeY()); // maps
56   TH2F* pMipM=new TH2F("RichCluMipMap"  ,"MIP map;x [cm];y [cm]"     ,1000 ,0  ,AliRICHParam::PcSizeX(),1000,0,AliRICHParam::PcSizeY());
57   TH2F* pCerM=new TH2F("RichCluCerMap"  ,"Ckov map;x [cm];y [cm]"    ,1000 ,0  ,AliRICHParam::PcSizeX(),1000,0,AliRICHParam::PcSizeY());
58  
59   
60   TClonesArray *pCluLst=new TClonesArray("AliRICHCluster");//tmp list of clusters
61   
62   for(Int_t iEvtN=0; iEvtN<iNevt; iEvtN++){
63     pAL->GetEvent(iEvtN);               
64     pRL->TreeD()->GetEntry(0); 
65     for(Int_t iChN=1;iChN<=AliRICHParam::kNch;iChN++){//chambers loop
66       if(pRich->Digs(iChN)->GetEntriesFast()>0) Dig2Clu(pRich->Digs(iChN),pCluLst,kFALSE);//cluster finder for the current chamber if any digits present
67     }//chambers loop
68     
69     for(Int_t iCluN=0 ; iCluN < pCluLst->GetEntriesFast() ; iCluN++){
70       AliRICHCluster *pClu = (AliRICHCluster*)pCluLst->At(iCluN);
71       Int_t cfm=0; for(Int_t iDig=0;iDig<pClu->Size();iDig++)  cfm+=pClu->Dig(iDig)->Cfm(); //collect ckov-fee-mip structure of current cluster
72       Int_t iNckov=cfm/1000000;      Int_t iNfee =cfm%1000000/1000;      Int_t iNmip =cfm%1000000%1000; 
73
74                                              pQ   ->Fill(pClu->Q()) ; pS   ->Fill(pClu->Size()) ; pM    ->Fill(pClu->X(),pClu->Y()); //all clusters                                      
75       if(iNckov!=0 && iNfee==0 && iNmip==0) {pCerQ->Fill(pClu->Q()) ; pCerS->Fill(pClu->Size()) ; pCerM ->Fill(pClu->X(),pClu->Y());}//ckov only cluster
76       if(iNckov==0 && iNfee==0 && iNmip!=0) {pMipQ->Fill(pClu->Q()) ; pMipS->Fill(pClu->Size()) ; pMipM ->Fill(pClu->X(),pClu->Y());}//mip only cluster
77                                        
78     }//clusters loop   
79   }//events loop
80   
81   pRL->UnloadDigits(); pAL->UnloadKinematics(); pAL->UnloadHeader();
82   TCanvas *pC=new TCanvas("RichCluQA",Form("QA for cluster from %i events",iNevt),1000,900); pC->Divide(3,3);
83   pC->cd(1);    pM->Draw();          pC->cd(2);    pQ->Draw();       pC->cd(3);    pS->Draw();        
84   pC->cd(4); pMipM->Draw();          pC->cd(5); pMipQ->Draw();       pC->cd(6); pMipS->Draw();        
85   pC->cd(7); pCerM->Draw();          pC->cd(8); pCerQ->Draw();       pC->cd(9); pCerS->Draw();        
86 }//CluQA()
87 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88 void AliRICHReconstructor::CheckPR()
89 {
90 //Pattern recognition with stack particles
91   TFile *pFile = new TFile("$(HOME)/RPR.root","RECREATE","RICH Pattern Recognition");
92   TNtupleD *hn = new TNtupleD("hn","ntuple","Pmod:Charge:TrackTheta:TrackPhi:TrackX:TrackY:MinX:MinY:ChargeMIP:ThetaCerenkov:NPhotons:MipIndex:Chamber:Particle");
93 //  printf("\n\n");
94 //  printf("Pattern Recognition done for event %5i",0);
95   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
96   AliMagF * magf = gAlice->Field();
97   AliTracker::SetFieldMap(magf,kTRUE);
98   for(Int_t iEvtN=0;iEvtN<pRich->GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvtN++) {
99     pRich->GetLoader()->GetRunLoader()->GetEvent(iEvtN);
100     AliRICHTracker *tr = new AliRICHTracker();
101     tr->RecWithStack(hn);
102     AliInfoClass(Form("Pattern Recognition done for event %i \b",iEvtN));
103 //    printf("\b\b\b\b\b%5i",iEvtN+1);
104   }
105   printf("\n\n");
106   pFile->Write();pFile->Close();
107 }
108 //__________________________________________________________________________________________________
109 void AliRICHReconstructor::Dig2Clu(TClonesArray *pDigLst,TClonesArray *pCluLst,Bool_t isTryUnfold)
110 {
111 //Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
112 //If pStack not 0 then also finds Ckov-Fee-Mip composition for formed clusters.  
113 //Puts all found clusters in the given clusters list. 
114 //Arguments: pDigLst     - list of digits provided not empty
115 //           pCluLst     - list of clusters, provided empty     
116 //           isTryUnfold - flag to choose between CoG and Mathieson fitting  
117 //  Returns: none    
118   TMatrixF digMap(1,AliRICHParam::NpadsX(),1,AliRICHParam::NpadsY());  digMap=(Float_t)-1; //digit map for one chamber reseted to -1
119   for(Int_t iDigN=0 ; iDigN < pDigLst->GetEntriesFast() ; iDigN++){                        //digits loop to fill digits map
120     AliRICHDigit *pDig= (AliRICHDigit*)pDigLst->At(iDigN);                                   //get current digit
121     digMap( pDig->PadX(), pDig->PadY() )=iDigN;                                              //fill the map, (padx,pady) cell takes digit index
122   }                                                                                        //digits loop to fill digits map 
123   
124   AliRICHCluster clu;                                                                      //tmp cluster to be used as current
125   
126   for(Int_t iDigN=0;iDigN<pDigLst->GetEntriesFast();iDigN++){                              //digits loop to form clusters list
127     AliRICHDigit *pDig=(AliRICHDigit*)pDigLst->At(iDigN);                                  //take current digit
128     if(!(pDig=UseDig(pDig->PadX(),pDig->PadY(),pDigLst,&digMap))) continue;                //this digit is already taken in FormClu(), go after next digit
129     FormClu(&clu,pDig,pDigLst,&digMap);                                                    //form cluster starting from this digit by recursion
130     clu.Solve(pCluLst,isTryUnfold);                                                        //solve this cluster and add all unfolded clusters to provided list  
131     clu.Reset();                                                                           //empty current cluster
132   }                                                                                        //digits loop to form clusters list
133 }//Dig2Clu()
134 //__________________________________________________________________________________________________
135 void  AliRICHReconstructor::FormClu(AliRICHCluster *pClu,AliRICHDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
136 {
137 //Forms the initial cluster as a sum of all adjascent digits. Starts from the given digit
138 //then calls itself recursevly  for all neighbours.
139 //Arguments: pClu - pointer to cluster being formed
140 //  Returns: none
141   pClu->DigAdd(pDig);//take this digit in cluster
142
143   Int_t x[4],y[4];
144   
145   Int_t iNnei=AliRICHParam::PadNeighbours(pDig->PadX(),pDig->PadY(),x,y);//returns in x,y all possible neighbours of the given padx,pady
146   for (Int_t i=0;i<iNnei;i++)
147     if((pDig=UseDig(x[i],y[i],pDigLst,pDigMap))) FormClu(pClu,pDig,pDigLst,pDigMap);   //check if this neighbour hit and mark it as taken  
148 }//FormClu()
149 //__________________________________________________________________________________________________
150 void AliRICHReconstructor::Reconstruct(AliRunLoader *pAL)const
151 {
152 //Invoked  by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
153 //Arguments: pAL - ALICE run loader pointer
154 //  Returns: none    
155   AliDebug(1,"Start.");
156   AliLoader *pRL=pAL->GetDetectorLoader("RICH");
157   AliRICH *pRich=(AliRICH*)pAL->GetAliRun()->GetDetector("RICH");//get pointers for RICH and RICH loader
158   pRL->LoadDigits();   
159   pRL->LoadRecPoints("recreate");
160   
161   for(Int_t iEvtN=0;iEvtN<pAL->GetNumberOfEvents();iEvtN++){//events loop
162     pAL->GetEvent(iEvtN); AliDebug(1,Form("Processing event %i...",iEvtN)); //switch current directory to next event    
163     pRL->TreeD()->GetEntry(0);  pRL->MakeTree("R");  pRich->MakeBranch("R");  //load digits to memory  and create branches for clusters              
164     for(Int_t iChN=1;iChN<=7;iChN++){//chambers loop
165       if(pRich->Digs(iChN)->GetEntriesFast()>0) Dig2Clu(pRich->Digs(iChN),pRich->Clus(iChN));//cluster finder for the current chamber if any digits present
166     }//chambers loop
167     pRL->TreeR()->Fill();            //fill tree for current event
168     pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
169     pRich->DigReset(); pRich->CluReset();
170   }//events loop  
171
172   pRL->UnloadDigits(); 
173   pRL->UnloadRecPoints();  
174     
175   AliDebug(1,"Stop.");      
176 }//Reconstruct(for simulated digits)
177 //__________________________________________________________________________________________________
178 void AliRICHReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const
179 {
180 //Invoked  by AliReconstruction to convert raw digits from DDL files to clusters
181 //Arguments: pAL - ALICE run loader pointer
182 //           pRR - ALICE raw reader pointer  
183 //  Returns: none    
184   AliLoader *pRL=pAL->GetDetectorLoader("RICH");  AliRICH *pRich=(AliRICH*)pAL->GetAliRun()->GetDetector("RICH");//get pointers for RICH and RICH loader
185   
186   AliRICHDigit dig; //tmp digit, raw digit will be converted to it
187   TClonesArray *pDigList=new TClonesArray("AliRICHDigit"); Int_t iDigCnt=0; //tmp list of digits for single chamber only
188   
189   Int_t iEvtN=0;
190   while(pRR->NextEvent()){//events loop
191     pAL->GetEvent(iEvtN++);
192     pRL->MakeTree("R");  pRich->MakeBranch("R");
193     
194     for(Int_t iChN=1;iChN<=7;iChN++){//chambers loop
195       pRR->Select(AliRICHDigit::kRichRawId,2*iChN-2,2*iChN-1);//select only DDL files for the current chamber      
196       UInt_t w32=0;
197       while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
198         UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
199         dig.Raw2Dig(ddl,w32);  
200         AliDebug(1,Form("Ch=%i DDL=%i raw=0x%x digit=(%3i,%3i,%3i,%3i) Q=%5.2f",iChN,ddl,w32,dig.C(),dig.S(),dig.PadX(),dig.PadY(),dig.Qdc()));
201         new((*pDigList)[iDigCnt++]) AliRICHDigit(dig); //add this digit to the tmp list
202       }//raw records loop
203       if(iDigCnt) Dig2Clu(pDigList,pRich->Clus(iChN));//cluster finder for the current chamber if any digits present
204       pRR->Reset();        
205       pDigList->Delete();  iDigCnt=0;//clean up list of digits for the current chamber
206     }//chambers loop
207     pRL->TreeR()->Fill();            //fill tree for current event
208     pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
209     pRich->CluReset();
210   }//events loop  
211   pRL->UnloadRecPoints();  
212 }//Reconstruct raw data
213 //__________________________________________________________________________________________________
214 void AliRICHReconstructor::RichAna(Int_t iNevMin,Int_t iNevMax,Bool_t askPatRec)
215 {
216   TFile *pFile=TFile::Open("AliESDs.root","read");
217   if(!pFile || !pFile->IsOpen()) {AliInfoClass("ESD file not open.");return;}      //open AliESDs.root                                                                    
218   TTree *pTree = (TTree*) pFile->Get("esdTree");
219   if(!pTree){AliInfoClass("ESD not found.");return;}                               //get ESD tree  
220   AliInfoClass("ESD found. Go ahead!");
221   
222   AliRICH *pRich=((AliRICH*)gAlice->GetDetector("RICH"));
223
224   AliMagF * magf = gAlice->Field();
225   AliTracker::SetFieldMap(magf,kTRUE);
226   pRich->GetLoader()->GetRunLoader()->LoadHeader();
227   pRich->GetLoader()->GetRunLoader()->LoadKinematics();
228   TString var1 = "Pmod:Charge:TrackTheta:TrackPhi:MinX:MinY:ThetaCerenkov:NPhotons:";
229   TString var2 = "ChargeMIP:Chamber:TOF:LengthTOF:prob1:prob2:prob3:";
230   TString var3 = "ErrPar1:ErrPar2:ErrPar3:Th1:Th2:Th3:nPhotBKG:pdgCode";
231   TString varList = var1+var2+var3;
232   
233   Double_t hnvec[30];
234
235 //  TFile *pFileRA = new TFile("$(HOME)/RichAna.root","RECREATE","RICH Pattern Recognition");
236   TFile *pFileRA = new TFile("./RichAna.root","RECREATE","RICH Pattern Recognition");
237   TNtupleD *hn = new TNtupleD("hn","ntuple",varList);
238   if(iNevMin<0) iNevMin=0;
239   if(iNevMin>iNevMax) {iNevMin=0;iNevMax=0;}  
240   if(iNevMax==0) iNevMax=999999;
241   if(pRich->GetLoader()->GetRunLoader()->GetNumberOfEvents()<iNevMax) iNevMax = pRich->GetLoader()->GetRunLoader()->GetNumberOfEvents();
242   AliESD *pESD=new AliESD;  pTree->SetBranchAddress("ESD", &pESD);
243   for(Int_t iEvtN=iNevMin;iEvtN<iNevMax;iEvtN++) {
244     pTree->GetEvent(iEvtN);
245     pRich->GetLoader()->GetRunLoader()->GetEvent(iEvtN);
246     pRich->GetLoader()->LoadRecPoints();
247     pRich->GetLoader()->TreeR()->GetEntry(0);
248     AliStack *pStack = pRich->GetLoader()->GetRunLoader()->Stack();
249 //Pattern recognition started
250     if(pESD->GetNumberOfTracks()) {
251       Int_t iNtracks=pESD->GetNumberOfTracks();
252       AliInfoClass(Form("Start with %i tracks",iNtracks));
253       for(Int_t iTrackN=0;iTrackN<iNtracks;iTrackN++){//ESD tracks loop
254         if(iTrackN%100==0)AliInfoClass(Form("Track %i to be processed",iTrackN));
255         AliRICHTracker *pTrRich = new AliRICHTracker();
256         if(askPatRec==kTRUE) pTrRich->PropagateBack(pESD);
257         AliESDtrack *pTrack = pESD->GetTrack(iTrackN);// get next reconstructed track
258         
259         Int_t lab=TMath::Abs(pTrack->GetLabel());
260         TParticle *pPart=pStack->Particle(lab);
261         Int_t code=pPart->GetPdgCode();
262
263         hnvec[0]=pTrack->GetP();
264         hnvec[1]=pTrack->GetSign();
265         Float_t dx,dy,trkTheta,trkPhi; 
266         pTrack->GetRICHthetaPhi(trkTheta,trkPhi); hnvec[2]=trkTheta; hnvec[3]=trkPhi;
267         pTrack->GetRICHdxdy(dx,dy);               hnvec[4]=dx;       hnvec[5]=dy;
268         hnvec[6]=pTrack->GetRICHsignal();
269         hnvec[7]=pTrack->GetRICHnclusters();
270         hnvec[9]=pTrack->GetRICHcluster()/1000000;
271         hnvec[8]=pTrack->GetRICHcluster()-hnvec[9]*1000000;
272         hnvec[10]=pTrack->GetTOFsignal();
273         hnvec[11]=pTrack->GetIntegratedLength();
274         Double_t prob[5];   pTrack->GetRICHpid(prob);
275         hnvec[12]=prob[0]+prob[1]+prob[2];     hnvec[13]=prob[3];     hnvec[14]=prob[4];
276         hnvec[15]=pTrRich->fErrPar[2];   hnvec[16]=pTrRich->fErrPar[3];   hnvec[17]=pTrRich->fErrPar[4];
277         for(Int_t i=0;i<3;i++) {
278           Double_t mass = AliRICHParam::fgMass[i+2];
279           Double_t refIndex=AliRICHParam::Instance()->IdxC6F14(AliRICHParam::EckovMean());
280           Double_t cosThetaTh = TMath::Sqrt(mass*mass+pTrack->GetP()*pTrack->GetP())/(refIndex*pTrack->GetP());
281           hnvec[18+i]=0;
282           if(cosThetaTh>=1) continue;
283           hnvec[18+i]= TMath::ACos(cosThetaTh);
284         }
285 //        if(askPatRec==kTRUE) hnvec[21]=pTrRich->fnPhotBKG; else hnvec[21]=0;
286         hnvec[22]=code;
287         hn->Fill(hnvec);
288       }
289     }
290     AliInfoClass(Form("Pattern Recognition done for event %i",iEvtN));
291   }
292   pFileRA->Write();pFileRA->Close();// close RichAna.root
293   delete pESD;  pFile->Close();//close AliESDs.root
294 }//RichAna()
295 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
296 void AliRICHReconstructor::Test(Bool_t isTryUnfold)
297 {
298 // Test the cluster finding algorithm by providing predifined set of digits
299 // Arguments: none
300 //   Returns: none  
301   TClonesArray *pDigTst=new TClonesArray("AliRICHDigit");       TClonesArray *pCluTst=new TClonesArray("AliRICHCluster");     
302   Int_t iDigCnt=0;
303   Int_t c,padx,pady,qdc;
304 //ckov cluster  
305   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 89,pady=13),qdc= 10);
306   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 90,pady=13),qdc=  7);
307   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 90,pady=12),qdc=  6);
308   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 91,pady=12),qdc=  7);
309 //mip cluster  
310   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 99,pady=21),qdc=  9);
311   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx= 99,pady=22),qdc= 26);
312   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx=100,pady=21),qdc= 39);
313   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx=100,pady=22),qdc=109);
314   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx=100,pady=23),qdc=  7);
315   new((*pDigTst)[iDigCnt++]) AliRICHDigit(AliRICHDigit::P2A(c=1,padx=101,pady=22),qdc= 11);
316   
317   Printf("Initial digits (1 ckov cluster and 1 mip cluster):"); pDigTst->Print();
318   Dig2Clu(pDigTst,pCluTst,isTryUnfold);   
319   Printf("Resulting clusters (expecting to have 2):"); pCluTst->Print(); 
320   delete pDigTst; delete pCluTst;
321 }//Test()
322 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
323 void AliRICHReconstructor::Test(TClonesArray *pDigLst,Bool_t isTryUnfold)
324 {
325 // Test the cluster finding algorithm for given list of digits. Note that list of digits will not be deleted.
326 // Arguments: pDigLst- list of digits
327 //   Returns: none  
328   TClonesArray *pCluLst=new TClonesArray("AliRICHCluster");     
329   Dig2Clu(pDigLst,pCluLst,isTryUnfold);   
330   
331   Int_t iNdig=pDigLst->GetEntriesFast();
332   Int_t iNclu=pCluLst->GetEntriesFast();
333   
334   TH2F *pH2=new TH2F("RDH2",Form("Tst dig->clu Digs: %i Clus: %i;cm;cm",iNdig,iNclu),AliRICHParam::NpadsX(),0,AliRICHParam::PcSizeX(),AliRICHParam::NpadsY(),0,AliRICHParam::PcSizeY());
335   pH2->SetStats(kFALSE);
336   for(Int_t iDig=0;iDig < iNdig;iDig++) {//digits loop
337     AliRICHDigit *pDig = (AliRICHDigit*)pDigLst->At(iDig);
338     TVector2 x2=AliRICHParam::Pad2Loc(pDig->Pad());
339     pH2->Fill(x2.X(),x2.Y(),pDig->Qdc());
340   }//digits loop
341   
342   TPolyMarker *pCluMarker=new TPolyMarker(iNclu); pCluMarker->SetMarkerStyle(5); pCluMarker->SetMarkerColor(kBlue);
343   for(Int_t iClu=0;iClu < iNclu;iClu++) {//clusters loop
344     AliRICHCluster *pClu = (AliRICHCluster*)pCluLst->At(iClu);
345     pCluMarker->SetNextPoint(pClu->X(),pClu->Y());
346   }//digits loop
347   
348   pH2->Draw("col");
349   pCluMarker->Draw();  
350   AliRICHParam::DrawSectors();
351   delete pCluLst;
352 }//Test()
353 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++