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