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