35e4d088ac4ec9fd4c609f7f2e4abbb49ee73008
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDReconstructor.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 // HMPID base class to reconstruct an event
17 //.
18 //.
19 //.
20 #include "AliHMPIDReconstructor.h" //class header
21 #include "AliHMPID.h"              //Reconstruct() 
22 #include "AliHMPIDCluster.h"       //Dig2Clu()
23 #include "AliHMPIDParam.h"         //FillEsd() 
24 #include <AliCDBEntry.h>           //ctor
25 #include <AliCDBManager.h>         //ctor
26 #include <AliESD.h>                //FillEsd()
27 #include <AliRawReader.h>          //Reconstruct() for raw digits
28 ClassImp(AliHMPIDReconstructor)
29
30 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 AliHMPIDReconstructor::AliHMPIDReconstructor():AliReconstructor(),fUserCut(0),fDaqSig(0),fDig(0),fClu(0)
32 {
33 //
34 //ctor
35 //
36
37   fClu=new TObjArray(AliHMPIDDigit::kMaxCh+1); fClu->SetOwner(kTRUE);
38   fDig=new TObjArray(AliHMPIDDigit::kMaxCh+1); fDig->SetOwner(kTRUE);
39
40   for(int i=AliHMPIDDigit::kMinCh;i<=AliHMPIDDigit::kMaxCh;i++){ 
41     fDig->AddAt(new TClonesArray("AliHMPIDDigit"),i);
42     fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
43   }
44   
45   AliCDBEntry *pUserCutEnt =AliCDBManager::Instance()->Get("HMPID/Calib/UserCut");    //contains TObjArray of 14 TObject with n. of sigmas to cut charge 
46   if(!pUserCutEnt) return;                                                            //No request from User to apply a more severe cut on pad charge  
47   
48   fUserCut = pUserCutEnt->GetObject()->GetUniqueID();
49
50   Printf(" usercut %i ",fUserCut);
51   
52   AliCDBEntry *pDaqCutEnt =AliCDBManager::Instance()->Get("HMPID/Calib/DaqSigCut");     //contains TObjArray of 14 TObject with n. of sigmas to cut charge 
53   if(!pDaqCutEnt) AliFatal("No pedestal sigmas cut!");
54   TObjArray *pDaqCut = (TObjArray*)pDaqCutEnt->GetObject();
55
56
57   for(Int_t ddl=0;ddl<14;ddl++){
58     if(fUserCut<=(Int_t)pDaqCut->At(ddl)->GetUniqueID()) continue;
59     AliCDBEntry *pDaqSigEnt =AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");  //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads 
60     if(!pDaqSigEnt) AliFatal("No pedestals from DAQ!");
61     TObjArray *pDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
62     SigConv(pDaqSig);
63     break;
64   }
65 }//AliHMPIDReconstructor
66 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
67 void AliHMPIDReconstructor::SigConv(TObjArray *pDaqSig)
68 {
69 // Conversion from the pedestal objects in to 7 TMatrixF(padX,padY)
70 //Arguments: pDaqSig pointer to the pedestal objects from OCDB
71 //   Returs: none
72 //
73   fDaqSig = new TObjArray(AliHMPIDDigit::kMaxCh+1);fDaqSig->SetOwner(kTRUE);
74   
75   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++) fDaqSig->AddAt(new TMatrixF(1,10,0,47),iCh);
76  
77   AliHMPIDDigit dig;
78    
79   for(Int_t ddl=0;ddl<=13;ddl++){
80     for(Int_t row=1;row<=24;row++)
81       for(Int_t dil=1;dil<=10;dil++)
82         for(Int_t adr=0;adr<=47;adr++){
83           TObjArray *pDdl = (TObjArray*)pDaqSig->At(ddl);
84           TObjArray *pRow = (TObjArray*)pDdl->At(row-1);
85           TMatrixF *pM = (TMatrixF*)pRow;
86           Float_t sigma = (*pM)(dil,adr);
87           dig.Raw(ddl,row,dil,adr);
88           TMatrixF* pMConv = (TMatrixF*)(fDaqSig->At(dig.Ch()));
89           (*pMConv)(dig.PadChX(),dig.PadChY()) = sigma;
90         }
91   }
92 }//SigConv()
93 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94 void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Bool_t isTryUnfold)
95 {
96 // Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
97 // Puts all found clusters in separate lists, one per clusters. 
98 // Arguments: pDigAll     - list of digits for all chambers 
99 //            pCluAll     - list of clusters for all chambers
100 //            isTryUnfold - flag to choose between CoG and Mathieson fitting  
101 //  Returns: none    
102   TMatrixF padMap(AliHMPIDDigit::kMinPx,AliHMPIDDigit::kMaxPcx,AliHMPIDDigit::kMinPy,AliHMPIDDigit::kMaxPcy);  //pads map for single chamber 0..159 x 0..143 
103   
104   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++){                  //chambers loop 
105     TClonesArray *pDigCur=(TClonesArray*)pDigAll->At(iCh);                                //get list of digits for current chamber
106     if(pDigCur->GetEntriesFast()==0) continue;                                            //no digits for current chamber
107   
108     padMap=(Float_t)-1;                                                                   //reset map to -1 (means no digit for this pad)  
109     TClonesArray *pCluCur=(TClonesArray*)pCluAll->At(iCh);                                //get list of clusters for current chamber
110     
111     for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){                              //digits loop to fill pads map
112       AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigCur->At(iDig);                             //get current digit
113       padMap( pDig->PadChX(), pDig->PadChY() )=iDig;                                      //fill the map, (padx,pady) cell takes digit index
114     }//digits loop to fill digits map 
115     
116     AliHMPIDCluster clu;                                                                  //tmp cluster to be used as current
117     for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){                              //digits loop for current chamber
118       AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCur->At(iDig);                              //take current digit
119       if(!(pDig=UseDig(pDig->PadChX(),pDig->PadChY(),pDigCur,&padMap))) continue;         //this digit is already taken in FormClu(), go after next digit
120       FormClu(&clu,pDig,pDigCur,&padMap);                                                 //form cluster starting from this digit by recursion
121       
122       clu.Solve(pCluCur,isTryUnfold);                                                     //solve this cluster and add all unfolded clusters to provided list  
123       clu.Reset();                                                                        //empty current cluster
124     }//digits loop for current chamber
125   }//chambers loop
126 }//Dig2Clu()
127 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
128 void  AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
129 {
130 // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly  for all neighbours.
131 // Arguments: pClu - pointer to cluster being formed
132 //   Returns: none   
133   pClu->DigAdd(pDig);//take this digit in cluster
134
135   Int_t cnt=0,cx[4],cy[4];
136   
137   if(pDig->PadPcX() != AliHMPIDDigit::kMinPx){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY()  ;cnt++;}       //left
138   if(pDig->PadPcX() != AliHMPIDDigit::kMaxPx){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY()  ;cnt++;}       //right
139   if(pDig->PadPcY() != AliHMPIDDigit::kMinPy){cx[cnt]=pDig->PadChX()  ; cy[cnt]=pDig->PadChY()-1;cnt++;}       //down
140   if(pDig->PadPcY() != AliHMPIDDigit::kMaxPy){cx[cnt]=pDig->PadChX()  ; cy[cnt]=pDig->PadChY()+1;cnt++;}       //up
141   
142   for (Int_t i=0;i<cnt;i++){//neighbours loop
143     if((pDig=UseDig(cx[i],cy[i],pDigLst,pDigMap))) FormClu(pClu,pDig,pDigLst,pDigMap);   //check if this neighbour pad fired and mark it as taken  
144   }//neighbours loop  
145 }//FormClu()
146 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
147 void AliHMPIDReconstructor::Reconstruct(TTree *pDigTree,TTree *pCluTree)const
148 {
149 //Invoked  by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
150 //Arguments: pDigTree - pointer to Digit tree
151 //           pCluTree - poitner to Cluster tree
152 //  Returns: none    
153   AliDebug(1,"Start.");
154   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++) {
155     pCluTree->Branch(Form("HMPID%d",iCh),&((*fClu)[iCh]),4000,0);
156     pDigTree->SetBranchAddress(Form("HMPID%d",iCh),&((*fDig)[iCh]));
157   }   
158   pDigTree->GetEntry(0);
159   Dig2Clu(fDig,fClu);              //cluster finder 
160   pCluTree->Fill();                //fill tree for current event
161   
162   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++){
163     fDig->At(iCh)->Clear();
164     fClu->At(iCh)->Clear();
165   }
166   
167   AliDebug(1,"Stop.");      
168 }//Reconstruct(for simulated digits)
169 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
170 void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const
171 {
172 //Invoked  by AliReconstruction to convert raw digits from DDL files to clusters
173 //Arguments: pAL - ALICE run loader pointer
174 //           pRR - ALICE raw reader pointer
175 //  Returns: none
176   AliLoader *pRL=pAL->GetDetectorLoader("HMPID");  AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
177
178   AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it
179
180   TObjArray digLst; Int_t iDigCnt[7]; for(Int_t i=0;i<7;i++){digLst.AddAt(new TClonesArray("AliHMPIDDigit"),i); iDigCnt[i]=0;} //tmp list of digits for allchambers
181
182   Int_t iEvtN=0;
183   while(pRR->NextEvent()){//events loop
184     pAL->GetEvent(iEvtN++);
185     pRL->MakeTree("R");  pRich->MakeBranch("R");
186
187     for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
188       pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber
189       UInt_t w32=0;
190       while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
191         UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
192         dig.Raw(w32,ddl);
193         AliDebug(1,Form("Ch=%i DDL=%i raw=0x%x digit=(%3i,%3i,%3i,%3i) Q=%5.2f",iCh,ddl,w32,dig.Ch(),dig.Pc(),dig.PadPcX(),dig.PadPcY(),dig.Q()));
194         if(!IsDigSurvive(&dig)) continue;                                             //sigma cut test
195         new((*((TClonesArray*)digLst.At(iCh)))[iDigCnt[iCh]++]) AliHMPIDDigit(dig); //add this digit to the tmp list
196       }//raw records loop
197       pRR->Reset();
198     }//chambers loop
199
200     Dig2Clu(&digLst,pRich->CluLst());//cluster finder for all chambers
201     for(Int_t i=0;i<7;i++){digLst.At(i)->Delete(); iDigCnt[i]=0;}                    //clean up list of digits for all chambers
202
203     pRL->TreeR()->Fill();            //fill tree for current event
204     pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
205     pRich->CluReset();
206   }//events loop
207   pRL->UnloadRecPoints();
208 }//Reconstruct raw data
209 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
210 void AliHMPIDReconstructor::ConvertDigits(AliRawReader *pRR,TTree *pDigTree)const
211 {
212 //Invoked  by AliReconstruction to convert raw digits from DDL files to digits
213 //Arguments: pRR - ALICE raw reader pointer
214 //           pDigTree - pointer to Digit tree
215 //  Returns: none    
216   AliDebug(1,"Start.");
217   AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it
218   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++) {
219     pDigTree->Branch(Form("HMPID%d",iCh),&((*fDig)[iCh]),4000,0);
220     pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber      
221     Int_t iDigCnt=0;
222     UInt_t w32=0;
223                 while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
224       UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
225       dig.Raw(w32,ddl);  
226       AliDebug(1,Form("Ch=%i DDL=%i raw=0x%x digit=(%3i,%3i,%3i,%3i) Q=%5.2f",iCh,ddl,w32,dig.Ch(),dig.Pc(),dig.PadPcX(),dig.PadPcY(),dig.Q()));
227       if(!IsDigSurvive(&dig)) continue;                                       //sigma cut test
228       new((*((TClonesArray*)fDig->At(iCh)))[iDigCnt++]) AliHMPIDDigit(dig); //add this digit to the tmp list
229     }//raw records loop
230     pRR->Reset();        
231   }//chambers loop
232   pDigTree->Fill();
233   for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++)fDig->At(iCh)->Clear();
234   AliDebug(1,"Stop.");
235 }//Reconstruct digits from raw digits
236 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
237 void AliHMPIDReconstructor::FillESD(AliRunLoader *, AliESD *pESD) const
238 {
239 // Calculates probability to be a electron-muon-pion-kaon-proton
240 // from the given Cerenkov angle and momentum assuming no initial particle composition
241 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
242
243   AliPID ppp; //needed
244   Double_t pid[AliPID::kSPECIES],h[AliPID::kSPECIES];
245    
246   for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
247     AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
248     if(pTrk->GetHMPIDsignal()<=0){//HMPID does not find anything reasonable for this track, assign 0.2 for all species
249       for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) pid[iPart]=1.0/AliPID::kSPECIES;
250       pTrk->SetHMPIDpid(pid);
251       continue;
252     } 
253     Double_t pmod = pTrk->GetP();
254     Double_t hTot=0;
255     for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
256       Double_t mass = AliPID::ParticleMass(iPart);
257       Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(AliHMPIDParam::Instance()->MeanIdxRad()*pmod);
258       if(cosThetaTh<1) //calculate the height of theoretical theta ckov on the gaus of experimental one
259         h[iPart] =TMath::Gaus(TMath::ACos(cosThetaTh),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);      
260       else             //beta < 1/ref. idx. => no light at all  
261         h[iPart] =0 ;       
262       hTot    +=h[iPart]; //total height of all theoretical heights for normalization
263     }//species loop
264      
265     Double_t hMin=TMath::Gaus(pTrk->GetHMPIDsignal()-4*TMath::Sqrt(pTrk->GetHMPIDchi2()),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);//5 sigma protection
266     
267     for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)//species loop to assign probabilities
268       if(hTot>hMin)  
269         pid[iPart]=h[iPart]/hTot;
270       else                               //all theoretical values are far away from experemental one
271         pid[iPart]=1.0/AliPID::kSPECIES; 
272     pTrk->SetHMPIDpid(pid);
273   }//ESD tracks loop
274 }//FillESD()
275 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++