]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDReconstructor.cxx
Coding rule violation corrected.
[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 "AliHMPIDPid.h"           //FillEsd() 
24 #include "AliHMPIDParam.h"         //FillEsd() 
25 #include <AliCDBEntry.h>           //ctor
26 #include <AliCDBManager.h>         //ctor
27 #include <AliESDEvent.h>           //FillEsd()
28 #include <AliRawReader.h>          //Reconstruct() for raw digits
29 #include <AliLog.h>                //
30 #include "AliHMPIDRawStream.h"     //ConvertDigits()
31 #include "AliHMPIDRecoParamV1.h"   //ctor
32
33 ClassImp(AliHMPIDReconstructor)
34
35 Int_t AliHMPIDReconstructor::fgStreamLevel = 0;        // stream (debug) level
36
37 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 AliHMPIDReconstructor::AliHMPIDReconstructor():AliReconstructor(),fDaqSig(0),fDig(0),fClu(0)
39 {
40 //
41 //ctor
42 //
43   AliHMPIDParam::Instance();                                                        //geometry loaded for reconstruction
44
45   fClu=new TObjArray(AliHMPIDParam::kMaxCh+1); fClu->SetOwner(kTRUE);
46   fDig=new TObjArray(AliHMPIDParam::kMaxCh+1); fDig->SetOwner(kTRUE);
47   
48   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){ 
49     fDig->AddAt(new TClonesArray("AliHMPIDDigit",24000),i);
50     TClonesArray *pClus = new TClonesArray("AliHMPIDCluster",24000);
51     pClus->SetUniqueID(i);
52     fClu->AddAt(pClus,i);
53   }
54  
55   AliCDBEntry *pDaqSigEnt =AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig");  //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads 
56   if(!pDaqSigEnt) AliFatal("No pedestals from DAQ!");
57   fDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
58   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){                  //chambers loop 
59     AliDebug(1,Form("DaqSigCut successfully loaded for chamber %i -> %i ",iCh,(Int_t)fDaqSig->At(iCh)->GetUniqueID()));
60   }
61 }//AliHMPIDReconstructor
62 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
63 void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Int_t *pUserCut,Bool_t isTryUnfold)
64 {
65 // Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
66 // Puts all found clusters in separate lists, one per clusters. 
67 // Arguments: pDigAll     - list of digits for all chambers 
68 //            pCluAll     - list of clusters for all chambers
69 //            isTryUnfold - flag to choose between CoG and Mathieson fitting  
70 //  Returns: none    
71   TMatrixF padMap(AliHMPIDParam::kMinPx,AliHMPIDParam::kMaxPcx,AliHMPIDParam::kMinPy,AliHMPIDParam::kMaxPcy);  //pads map for single chamber 0..159 x 0..143 
72   
73   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){                  //chambers loop 
74     TClonesArray *pDigCur=(TClonesArray*)pDigAll->At(iCh);                                //get list of digits for current chamber
75     if(pDigCur->GetEntriesFast()==0) continue;                                            //no digits for current chamber
76   
77     padMap=(Float_t)-1;                                                                   //reset map to -1 (means no digit for this pad)  
78     TClonesArray *pCluCur=(TClonesArray*)pCluAll->At(iCh);                                //get list of clusters for current chamber
79     
80     for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){                              //digits loop to fill pads map
81       AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigCur->At(iDig);                             //get current digit
82       padMap( pDig->PadChX(), pDig->PadChY() )=iDig;                                      //fill the map, (padx,pady) cell takes digit index
83     }//digits loop to fill digits map 
84     
85     AliHMPIDCluster clu;                                                                  //tmp cluster to be used as current
86     for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){                              //digits loop for current chamber
87       AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCur->At(iDig);                              //take current digit
88       if(!(pDig=UseDig(pDig->PadChX(),pDig->PadChY(),pDigCur,&padMap))) continue;         //this digit is already taken in FormClu(), go after next digit
89       FormClu(&clu,pDig,pDigCur,&padMap);                                                 //form cluster starting from this digit by recursion
90       clu.Solve(pCluCur,pUserCut,isTryUnfold);                                 //solve this cluster and add all unfolded clusters to provided list  
91       clu.Reset();                                                                        //empty current cluster
92     }//digits loop for current chamber
93   }//chambers loop
94 }//Dig2Clu()
95 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
96 void  AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
97 {
98 // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly  for all neighbours.
99 // Arguments: pClu - pointer to cluster being formed
100 //   Returns: none   
101   pClu->DigAdd(pDig);//take this digit in cluster
102
103   Int_t cnt=0,cx[4],cy[4];
104   
105   if(pDig->PadPcX() != AliHMPIDParam::kMinPx){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY()  ;cnt++;}       //left
106   if(pDig->PadPcX() != AliHMPIDParam::kMaxPx){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY()  ;cnt++;}       //right
107   if(pDig->PadPcY() != AliHMPIDParam::kMinPy){cx[cnt]=pDig->PadChX()  ; cy[cnt]=pDig->PadChY()-1;cnt++;}       //down
108   if(pDig->PadPcY() != AliHMPIDParam::kMaxPy){cx[cnt]=pDig->PadChX()  ; cy[cnt]=pDig->PadChY()+1;cnt++;}       //up
109   
110   for (Int_t i=0;i<cnt;i++){//neighbours loop
111     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  
112   }//neighbours loop  
113 }//FormClu()
114 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
115 void AliHMPIDReconstructor::Reconstruct(TTree *pDigTree,TTree *pCluTree)const
116 {
117 //Invoked  by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
118 //Arguments: pDigTree - pointer to Digit tree
119 //           pCluTree - poitner to Cluster tree
120 //  Returns: none    
121   AliDebug(1,"Start.");
122   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
123     pCluTree->Branch(Form("HMPID%d",iCh),&((*fClu)[iCh]),7);
124     pDigTree->SetBranchAddress(Form("HMPID%d",iCh),&((*fDig)[iCh]));
125   }   
126   
127   Int_t pUserCut[7];
128   
129    if(AliHMPIDReconstructor::GetRecoParam()) {
130     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
131       pUserCut[iCh] = AliHMPIDReconstructor::GetRecoParam()->GetHmpUserCut(iCh);
132       AliDebug(1,Form("UserCut successfully loaded (from AliHMPIDRecoParamV1) for chamber %i -> %i ",iCh,pUserCut[iCh]));
133     }
134   }
135   else {
136     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
137       pUserCut[iCh] = 4;                                                                             // minimal requirement for sigma cut
138       AliDebug(1,Form("UserCut loaded from defaults for chamber %i -> %i ",iCh,pUserCut[iCh]));
139       AliDebug(1,Form("Cannot get AliHMPIDRecoParamV1!"));
140       }
141     }  
142   
143   pDigTree->GetEntry(0);
144   Dig2Clu(fDig,fClu,pUserCut);     //cluster finder 
145   pCluTree->Fill();                //fill tree for current event
146   
147   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
148     fDig->At(iCh)->Clear();
149     fClu->At(iCh)->Clear();
150   }
151   
152   
153   AliDebug(1,"Stop.");      
154 }//Reconstruct(for simulated digits)
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 //  Int_t digcnt=0;
164   
165   Int_t iDigCnt[7]={0,0,0,0,0,0,0};
166
167   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
168     pDigTree->Branch(Form("HMPID%d",iCh),&((*fDig)[iCh]));
169   }
170     
171   Int_t pUserCut[7];
172   
173    if(AliHMPIDReconstructor::GetRecoParam()) {
174     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
175       pUserCut[iCh] = AliHMPIDReconstructor::GetRecoParam()->GetHmpUserCut(iCh);
176       AliDebug(1,Form("UserCut successfully loaded (from AliHMPIDRecoParamV1) for chamber %i -> %i ",iCh,pUserCut[iCh]));
177     }
178   }
179   else {
180     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
181       pUserCut[iCh] = 4;                                                                             // minimal requirement for sigma cut
182       AliDebug(1,Form("UserCut loaded from defaults for chamber %i -> %i ",iCh,pUserCut[iCh]));
183       AliDebug(1,Form("Cannot get AliHMPIDRecoParamV1!"));
184       }
185     }  
186       
187   AliHMPIDRawStream stream(pRR);    
188   
189   while(stream.Next())
190   {
191     Int_t ch = AliHMPIDParam::DDL2C(stream.GetDDLNumber());
192     for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
193       AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
194       if(!IsDigSurvive(pUserCut,&dig)) continue; 
195       new((*((TClonesArray*)fDig->At(ch)))[iDigCnt[ch]++]) AliHMPIDDigit(dig); //add this digit to the tmp list 
196     }
197   }
198   
199   pDigTree->Fill();
200   
201   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)fDig->At(iCh)->Clear();
202   
203   
204   AliDebug(1,"Stop.");
205 }//Reconstruct digits from raw digits
206 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
207 void AliHMPIDReconstructor::FillESD(TTree */*digitsTree*/, TTree */*clustersTree*/, AliESDEvent *pESD) const
208 {
209 // Fill ESD with all the infos from HMPID
210 // Probability vector from AliHMPIDPid
211 //...
212   AliHMPIDPid pID;
213   Double_t prob[AliPID::kSPECIES];
214   
215   for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
216     AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
217     pID.FindPid(pTrk,AliPID::kSPECIES,prob);
218     pTrk->SetHMPIDpid(prob);
219   }//ESD tracks loop
220 }//FillESD()