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