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