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