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