1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // HMPID base class to reconstruct an event
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 "AliHMPIDRawStream.h" //ConvertDigits()
30 #include "AliHMPIDRecoParam.h" //ctor
31 ClassImp(AliHMPIDReconstructor)
33 AliHMPIDRecoParam* AliHMPIDReconstructor::fgkRecoParam =0; //
34 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
35 AliHMPIDReconstructor::AliHMPIDReconstructor():AliReconstructor(),fUserCut(0),fDaqSig(0),fDig(0),fClu(0)
40 AliHMPIDParam::Instance(); //geometry loaded for reconstruction
41 fUserCut = new Int_t[7];
42 fClu=new TObjArray(AliHMPIDParam::kMaxCh+1); fClu->SetOwner(kTRUE);
43 fDig=new TObjArray(AliHMPIDParam::kMaxCh+1); fDig->SetOwner(kTRUE);
45 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
46 fDig->AddAt(new TClonesArray("AliHMPIDDigit",24000),i);
47 TClonesArray *pClus = new TClonesArray("AliHMPIDCluster",24000);
48 pClus->SetUniqueID(i);
52 if(fgkRecoParam!=0x0 && fgkRecoParam->GetUserCutMode()==kFALSE)
54 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
55 fUserCut[iCh] = fgkRecoParam->GetUserCut(iCh);
56 AliDebug(1,Form("UserCut successfully loaded (from RecoParam) for chamber %i -> %i ",iCh,fUserCut[iCh]));
60 AliCDBEntry *pUserCutEnt =AliCDBManager::Instance()->Get("HMPID/Calib/UserCut"); //contains TObjArray of 14 TObject with n. of sigmas to cut charge
62 TObjArray *pUserCut = (TObjArray*)pUserCutEnt->GetObject();
63 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ //chambers loop
64 fUserCut[iCh] = pUserCut->At(iCh)->GetUniqueID();
65 AliDebug(1,Form("UserCut successfully loaded (from OCDB) for chamber %i -> %i ",iCh,fUserCut[iCh]));
70 AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 7 TF1
71 if(!pQthreEnt) AliFatal("No Qthre available");
72 TObjArray *pQthre = (TObjArray*)pQthreEnt->GetObject();
73 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
74 for(Int_t isec=0;isec<=5;isec++) {
75 TF1 *pfQthre = (TF1*)pQthre->At(6*iCh+isec);
77 pfQthre->GetRange(tMin,tMax);
78 Double_t qthre=pfQthre->Eval(tMin);
79 AliDebug(1,Form("Qthre successfully loaded for chamber %i sector %i -> %f ",iCh,isec,qthre));
84 AliCDBEntry *pDaqSigEnt =AliCDBManager::Instance()->Get("HMPID/Calib/DaqSig"); //contains TObjArray of TObjArray 14 TMatrixF sigmas values for pads
85 if(!pDaqSigEnt) AliFatal("No pedestals from DAQ!");
86 fDaqSig = (TObjArray*)pDaqSigEnt->GetObject();
87 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ //chambers loop
88 AliDebug(1,Form("DaqSigCut successfully loaded for chamber %i -> %i ",iCh,(Int_t)fDaqSig->At(iCh)->GetUniqueID()));
90 }//AliHMPIDReconstructor
91 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
92 void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Int_t *pUserCut,Bool_t isTryUnfold)
94 // Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
95 // Puts all found clusters in separate lists, one per clusters.
96 // Arguments: pDigAll - list of digits for all chambers
97 // pCluAll - list of clusters for all chambers
98 // isTryUnfold - flag to choose between CoG and Mathieson fitting
100 TMatrixF padMap(AliHMPIDParam::kMinPx,AliHMPIDParam::kMaxPcx,AliHMPIDParam::kMinPy,AliHMPIDParam::kMaxPcy); //pads map for single chamber 0..159 x 0..143
102 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){ //chambers loop
103 TClonesArray *pDigCur=(TClonesArray*)pDigAll->At(iCh); //get list of digits for current chamber
104 if(pDigCur->GetEntriesFast()==0) continue; //no digits for current chamber
106 padMap=(Float_t)-1; //reset map to -1 (means no digit for this pad)
107 TClonesArray *pCluCur=(TClonesArray*)pCluAll->At(iCh); //get list of clusters for current chamber
109 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop to fill pads map
110 AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigCur->At(iDig); //get current digit
111 padMap( pDig->PadChX(), pDig->PadChY() )=iDig; //fill the map, (padx,pady) cell takes digit index
112 }//digits loop to fill digits map
114 AliHMPIDCluster clu; //tmp cluster to be used as current
115 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop for current chamber
116 AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCur->At(iDig); //take current digit
117 if(!(pDig=UseDig(pDig->PadChX(),pDig->PadChY(),pDigCur,&padMap))) continue; //this digit is already taken in FormClu(), go after next digit
118 FormClu(&clu,pDig,pDigCur,&padMap); //form cluster starting from this digit by recursion
119 clu.Solve(pCluCur,pUserCut,isTryUnfold); //solve this cluster and add all unfolded clusters to provided list
120 clu.Reset(); //empty current cluster
121 }//digits loop for current chamber
124 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
125 void AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
127 // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly for all neighbours.
128 // Arguments: pClu - pointer to cluster being formed
130 pClu->DigAdd(pDig);//take this digit in cluster
132 Int_t cnt=0,cx[4],cy[4];
134 if(pDig->PadPcX() != AliHMPIDParam::kMinPx){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY() ;cnt++;} //left
135 if(pDig->PadPcX() != AliHMPIDParam::kMaxPx){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY() ;cnt++;} //right
136 if(pDig->PadPcY() != AliHMPIDParam::kMinPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()-1;cnt++;} //down
137 if(pDig->PadPcY() != AliHMPIDParam::kMaxPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()+1;cnt++;} //up
139 for (Int_t i=0;i<cnt;i++){//neighbours loop
140 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
143 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
144 void AliHMPIDReconstructor::Reconstruct(TTree *pDigTree,TTree *pCluTree)const
146 //Invoked by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
147 //Arguments: pDigTree - pointer to Digit tree
148 // pCluTree - poitner to Cluster tree
150 AliDebug(1,"Start.");
151 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
152 pCluTree->Branch(Form("HMPID%d",iCh),&((*fClu)[iCh]),7);
153 pDigTree->SetBranchAddress(Form("HMPID%d",iCh),&((*fDig)[iCh]));
155 pDigTree->GetEntry(0);
156 Dig2Clu(fDig,fClu,fUserCut); //cluster finder
157 pCluTree->Fill(); //fill tree for current event
159 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
160 fDig->At(iCh)->Clear();
161 fClu->At(iCh)->Clear();
165 }//Reconstruct(for simulated digits)
166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167 void AliHMPIDReconstructor::ConvertDigits(AliRawReader *pRR,TTree *pDigTree)const
169 //Invoked by AliReconstruction to convert raw digits from DDL files to digits
170 //Arguments: pRR - ALICE raw reader pointer
171 // pDigTree - pointer to Digit tree
173 AliDebug(1,"Start.");
176 Int_t iDigCnt[7]={0,0,0,0,0,0,0};
178 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
179 pDigTree->Branch(Form("HMPID%d",iCh),&((*fDig)[iCh]));
182 AliHMPIDRawStream stream(pRR);
186 Int_t ch = AliHMPIDParam::DDL2C(stream.GetDDLNumber());
187 for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
188 AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
189 if(!IsDigSurvive(&dig)) continue;
190 new((*((TClonesArray*)fDig->At(ch)))[iDigCnt[ch]++]) AliHMPIDDigit(dig); //add this digit to the tmp list
196 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++)fDig->At(iCh)->Clear();
199 }//Reconstruct digits from raw digits
200 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
201 void AliHMPIDReconstructor::FillESD(TTree */*digitsTree*/, TTree */*clustersTree*/, AliESDEvent *pESD) const
203 // Fill ESD with all the infos from HMPID
204 // Probability vector from AliHMPIDPid
207 Double_t prob[AliPID::kSPECIES];
209 for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
210 AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
211 pID.FindPid(pTrk,AliPID::kSPECIES,prob);
212 pTrk->SetHMPIDpid(prob);