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 #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 <AliRunLoader.h> //Reconstruct() for simulated digits
22 #include <AliRawReader.h> //Reconstruct() for raw digits
23 #include <AliRun.h> //Reconstruct()
24 ClassImp(AliHMPIDReconstructor)
26 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Bool_t isTryUnfold)
29 // Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
30 // Puts all found clusters in separate lists, one per clusters.
31 // Arguments: pDigAll - list of digits for all chambers
32 // pCluAll - list of clusters for all chambers
33 // isTryUnfold - flag to choose between CoG and Mathieson fitting
35 TMatrixF padMap(AliHMPIDDigit::kMinPx,AliHMPIDDigit::kMaxPcx,AliHMPIDDigit::kMinPy,AliHMPIDDigit::kMaxPcy); //pads map for single chamber 0..159 x 0..143
37 for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++){ //chambers loop
38 TClonesArray *pDigCur=(TClonesArray*)pDigAll->At(iCh); //get list of digits for current chamber
39 if(pDigCur->GetEntriesFast()==0) continue; //no digits for current chamber
41 padMap=(Float_t)-1; //reset map to -1 (means no digit for this pad)
42 TClonesArray *pCluCur=(TClonesArray*)pCluAll->At(iCh); //get list of clusters for current chamber
44 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop to fill pads map
45 AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigCur->At(iDig); //get current digit
46 padMap( pDig->PadChX(), pDig->PadChY() )=iDig; //fill the map, (padx,pady) cell takes digit index
47 }//digits loop to fill digits map
49 AliHMPIDCluster clu; //tmp cluster to be used as current
51 for(Int_t iDig=0;iDig<pDigCur->GetEntriesFast();iDig++){ //digits loop for current chamber
52 AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCur->At(iDig); //take current digit
53 if(!(pDig=UseDig(pDig->PadChX(),pDig->PadChY(),pDigCur,&padMap))) continue; //this digit is already taken in FormClu(), go after next digit
54 FormClu(&clu,pDig,pDigCur,&padMap); //form cluster starting from this digit by recursion
56 clu.Solve(pCluCur,isTryUnfold); //solve this cluster and add all unfolded clusters to provided list
57 clu.Reset(); //empty current cluster
58 }//digits loop for current chamber
61 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
62 void AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
64 // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly for all neighbours.
65 // Arguments: pClu - pointer to cluster being formed
67 pClu->DigAdd(pDig);//take this digit in cluster
69 Int_t cnt=0,cx[4],cy[4];
71 if(pDig->PadPcX() != AliHMPIDDigit::kMinPx){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY() ;cnt++;} //left
72 if(pDig->PadPcX() != AliHMPIDDigit::kMaxPx){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY() ;cnt++;} //right
73 if(pDig->PadPcY() != AliHMPIDDigit::kMinPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()-1;cnt++;} //down
74 if(pDig->PadPcY() != AliHMPIDDigit::kMaxPy){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()+1;cnt++;} //up
76 for (Int_t i=0;i<cnt;i++){//neighbours loop
77 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
80 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
81 void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL)const
83 //Invoked by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
84 //Arguments: pAL - ALICE run loader pointer
87 AliLoader *pRL=pAL->GetDetectorLoader("HMPID");
88 AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
90 pRL->LoadRecPoints("recreate");
92 for(Int_t iEvtN=0;iEvtN<pAL->GetNumberOfEvents();iEvtN++){//events loop
93 pAL->GetEvent(iEvtN); AliDebug(1,Form("Processing event %i...",iEvtN)); //switch current directory to next event
94 pRL->TreeD()->GetEntry(0); pRL->MakeTree("R"); pRich->MakeBranch("R"); //load digits to memory and create branches for clusters
96 Dig2Clu(pRich->DigLst(),pRich->CluLst());//cluster finder
98 pRL->TreeR()->Fill(); //fill tree for current event
99 pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
100 pRich->DigReset(); pRich->CluReset();
104 pRL->UnloadRecPoints();
107 }//Reconstruct(for simulated digits)
108 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
109 void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const
111 //Invoked by AliReconstruction to convert raw digits from DDL files to clusters
112 //Arguments: pAL - ALICE run loader pointer
113 // pRR - ALICE raw reader pointer
115 AliLoader *pRL=pAL->GetDetectorLoader("HMPID"); AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
117 AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it
119 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 all chambers
122 while(pRR->NextEvent()){//events loop
123 pAL->GetEvent(iEvtN++);
124 pRL->MakeTree("R"); pRich->MakeBranch("R");
126 for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
127 pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber
129 while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
130 UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
132 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()));
133 new((*((TClonesArray*)digLst.At(iCh)))[iDigCnt[iCh]++]) AliHMPIDDigit(dig); //add this digit to the tmp list
138 Dig2Clu(&digLst,pRich->CluLst());//cluster finder for all chambers
139 for(Int_t i=0;i<7;i++){digLst.At(i)->Delete(); iDigCnt[i]=0;} //clean up list of digits for all chambers
141 pRL->TreeR()->Fill(); //fill tree for current event
142 pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
145 pRL->UnloadRecPoints();
146 }//Reconstruct raw data
147 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
148 void AliHMPIDReconstructor::FillESD(AliRunLoader *, AliESD *pESD) const
150 // Calculates probability to be a electron-muon-pion-kaon-proton
151 // from the given Cerenkov angle and momentum assuming no initial particle composition
152 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
155 Double_t pid[AliPID::kSPECIES],h[AliPID::kSPECIES];
157 for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
158 AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
159 if(pTrk->GetHMPIDsignal()<=0){//HMPID does not find anything reasonable for this track, assign 0.2 for all species
160 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) pid[iPart]=1.0/AliPID::kSPECIES;
161 pTrk->SetHMPIDpid(pid);
164 Double_t pmod = pTrk->GetP();
166 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
167 Double_t mass = AliPID::ParticleMass(iPart);
168 Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(AliHMPIDParam::Instance()->MeanIdxRad()*pmod);
169 if(cosThetaTh<1) //calculate the height of theoretical theta ckov on the gaus of experimental one
170 h[iPart] =TMath::Gaus(TMath::ACos(cosThetaTh),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);
171 else //beta < 1/ref. idx. => no light at all
173 hTot +=h[iPart]; //total height of all theoretical heights for normalization
176 Double_t hMin=TMath::Gaus(pTrk->GetHMPIDsignal()-4*TMath::Sqrt(pTrk->GetHMPIDchi2()),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);//5 sigma protection
178 for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)//species loop to assign probabilities
180 pid[iPart]=h[iPart]/hTot;
181 else //all theoretical values are far away from experemental one
182 pid[iPart]=1.0/AliPID::kSPECIES;
183 pTrk->SetHMPIDpid(pid);
186 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++