d3da6dc4 |
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() |
a1d55ff3 |
18 | #include "AliHMPIDCluster.h" //Dig2Clu() |
d3da6dc4 |
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() |
d3da6dc4 |
24 | ClassImp(AliHMPIDReconstructor) |
25 | |
26 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
3c6274c1 |
27 | void AliHMPIDReconstructor::Dig2Clu(TObjArray *pDigAll,TObjArray *pCluAll,Bool_t isTryUnfold) |
d3da6dc4 |
28 | { |
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. |
3c6274c1 |
31 | // Arguments: pDigAll - list of digits for all chambers |
32 | // pCluAll - list of clusters for all chambers |
d3da6dc4 |
33 | // isTryUnfold - flag to choose between CoG and Mathieson fitting |
34 | // Returns: none |
a1d55ff3 |
35 | TMatrixF padMap(AliHMPIDDigit::kMinPx,AliHMPIDDigit::kMaxPcx,AliHMPIDDigit::kMinPy,AliHMPIDDigit::kMaxPcy); //pads map for single chamber 0..159 x 0..143 |
3c6274c1 |
36 | |
a1d55ff3 |
37 | for(Int_t iCh=AliHMPIDDigit::kMinCh;iCh<=AliHMPIDDigit::kMaxCh;iCh++){ //chambers loop |
3c6274c1 |
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 |
d3da6dc4 |
40 | |
3c6274c1 |
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 |
43 | |
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 |
48 | |
49 | AliHMPIDCluster clu; //tmp cluster to be used as current |
d3da6dc4 |
50 | |
3c6274c1 |
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 |
55 | |
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 |
59 | }//chambers loop |
d3da6dc4 |
60 | }//Dig2Clu() |
61 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
62 | void AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap) |
63 | { |
a1d55ff3 |
64 | // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit then calls itself recursevly for all neighbours. |
d3da6dc4 |
65 | // Arguments: pClu - pointer to cluster being formed |
a1d55ff3 |
66 | // Returns: none |
d3da6dc4 |
67 | pClu->DigAdd(pDig);//take this digit in cluster |
68 | |
3c6274c1 |
69 | Int_t cnt=0,cx[4],cy[4]; |
d3da6dc4 |
70 | |
a1d55ff3 |
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 |
d3da6dc4 |
75 | |
a1d55ff3 |
76 | for (Int_t i=0;i<cnt;i++){//neighbours loop |
3c6274c1 |
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 |
a1d55ff3 |
78 | }//neighbours loop |
d3da6dc4 |
79 | }//FormClu() |
80 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
81 | void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL)const |
82 | { |
83 | //Invoked by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data |
84 | //Arguments: pAL - ALICE run loader pointer |
85 | // Returns: none |
86 | AliDebug(1,"Start."); |
87 | AliLoader *pRL=pAL->GetDetectorLoader("HMPID"); |
88 | AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader |
89 | pRL->LoadDigits(); |
90 | pRL->LoadRecPoints("recreate"); |
91 | |
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 |
95 | |
3c6274c1 |
96 | Dig2Clu(pRich->DigLst(),pRich->CluLst());//cluster finder |
d3da6dc4 |
97 | |
98 | pRL->TreeR()->Fill(); //fill tree for current event |
99 | pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event |
100 | pRich->DigReset(); pRich->CluReset(); |
101 | }//events loop |
102 | |
103 | pRL->UnloadDigits(); |
104 | pRL->UnloadRecPoints(); |
105 | |
106 | AliDebug(1,"Stop."); |
107 | }//Reconstruct(for simulated digits) |
a1d55ff3 |
108 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
d3da6dc4 |
109 | void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const |
110 | { |
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 |
114 | // Returns: none |
115 | AliLoader *pRL=pAL->GetDetectorLoader("HMPID"); AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader |
116 | |
117 | AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it |
118 | |
3c6274c1 |
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 |
120 | |
d3da6dc4 |
121 | Int_t iEvtN=0; |
122 | while(pRR->NextEvent()){//events loop |
123 | pAL->GetEvent(iEvtN++); |
124 | pRL->MakeTree("R"); pRich->MakeBranch("R"); |
125 | |
126 | for(Int_t iCh=0;iCh<7;iCh++){//chambers loop |
d3da6dc4 |
127 | pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber |
128 | UInt_t w32=0; |
129 | while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files) |
130 | UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13 |
1d4857c5 |
131 | dig.Raw(w32,ddl); |
da08475b |
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())); |
3c6274c1 |
133 | new((*((TClonesArray*)digLst.At(iCh)))[iDigCnt[iCh]++]) AliHMPIDDigit(dig); //add this digit to the tmp list |
d3da6dc4 |
134 | }//raw records loop |
d3da6dc4 |
135 | pRR->Reset(); |
d3da6dc4 |
136 | }//chambers loop |
3c6274c1 |
137 | |
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 |
140 | |
d3da6dc4 |
141 | pRL->TreeR()->Fill(); //fill tree for current event |
142 | pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event |
143 | pRich->CluReset(); |
144 | }//events loop |
145 | pRL->UnloadRecPoints(); |
146 | }//Reconstruct raw data |
147 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
d3da6dc4 |
148 | void AliHMPIDReconstructor::FillESD(AliRunLoader *, AliESD *pESD) const |
149 | { |
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) |
153 | |
154 | AliPID ppp; //needed |
155 | Double_t pid[AliPID::kSPECIES],h[AliPID::kSPECIES]; |
156 | |
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); |
162 | continue; |
163 | } |
164 | Double_t pmod = pTrk->GetP(); |
165 | Double_t hTot=0; |
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); |
a1d55ff3 |
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); |
d3da6dc4 |
171 | else //beta < 1/ref. idx. => no light at all |
172 | h[iPart] =0 ; |
173 | hTot +=h[iPart]; //total height of all theoretical heights for normalization |
174 | }//species loop |
175 | |
176 | Double_t hMin=TMath::Gaus(pTrk->GetHMPIDsignal()-4*TMath::Sqrt(pTrk->GetHMPIDchi2()),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);//5 sigma protection |
177 | |
178 | for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)//species loop to assign probabilities |
179 | if(hTot>hMin) |
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); |
184 | }//ESD tracks loop |
a1d55ff3 |
185 | }//FillESD() |
d3da6dc4 |
186 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |