]>
Commit | Line | Data |
---|---|---|
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() | |
18 | #include "AliHMPIDCluster.h" //CluQA() | |
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 | |
3c6274c1 | 35 | TMatrixF padMap(AliHMPIDDigit::kPadAllX,AliHMPIDDigit::kPadAllY); //pads map for single chamber 0..159 x 0..143 |
36 | ||
37 | for(Int_t iCh=0;iCh<7;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 | |
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 | { | |
da08475b | 64 | // Forms the initial cluster as a combination of all adjascent digits. Starts from the given digit |
65 | // then calls itself recursevly for all possible neighbours. | |
d3da6dc4 | 66 | // Arguments: pClu - pointer to cluster being formed |
3c6274c1 | 67 | // Returns: none |
d3da6dc4 | 68 | pClu->DigAdd(pDig);//take this digit in cluster |
69 | ||
3c6274c1 | 70 | Int_t cnt=0,cx[4],cy[4]; |
d3da6dc4 | 71 | |
3c6274c1 | 72 | if(pDig->PadPcX() != 0 ){cx[cnt]=pDig->PadChX()-1; cy[cnt]=pDig->PadChY() ;cnt++;} //left |
73 | if(pDig->PadPcX() != AliHMPIDDigit::kPadPcX-1){cx[cnt]=pDig->PadChX()+1; cy[cnt]=pDig->PadChY() ;cnt++;} //right | |
74 | if(pDig->PadPcY() != 0 ){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()-1;cnt++;} //down | |
75 | if(pDig->PadPcY() != AliHMPIDDigit::kPadPcY-1){cx[cnt]=pDig->PadChX() ; cy[cnt]=pDig->PadChY()+1;cnt++;} //up | |
d3da6dc4 | 76 | |
da08475b | 77 | for (Int_t i=0;i<cnt;i++) |
3c6274c1 | 78 | 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 |
d3da6dc4 | 79 | |
d3da6dc4 | 80 | }//FormClu() |
81 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
82 | void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL)const | |
83 | { | |
84 | //Invoked by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data | |
85 | //Arguments: pAL - ALICE run loader pointer | |
86 | // Returns: none | |
87 | AliDebug(1,"Start."); | |
88 | AliLoader *pRL=pAL->GetDetectorLoader("HMPID"); | |
89 | AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader | |
90 | pRL->LoadDigits(); | |
91 | pRL->LoadRecPoints("recreate"); | |
92 | ||
93 | for(Int_t iEvtN=0;iEvtN<pAL->GetNumberOfEvents();iEvtN++){//events loop | |
94 | pAL->GetEvent(iEvtN); AliDebug(1,Form("Processing event %i...",iEvtN)); //switch current directory to next event | |
95 | pRL->TreeD()->GetEntry(0); pRL->MakeTree("R"); pRich->MakeBranch("R"); //load digits to memory and create branches for clusters | |
96 | ||
3c6274c1 | 97 | Dig2Clu(pRich->DigLst(),pRich->CluLst());//cluster finder |
d3da6dc4 | 98 | |
99 | pRL->TreeR()->Fill(); //fill tree for current event | |
100 | pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event | |
101 | pRich->DigReset(); pRich->CluReset(); | |
102 | }//events loop | |
103 | ||
104 | pRL->UnloadDigits(); | |
105 | pRL->UnloadRecPoints(); | |
106 | ||
107 | AliDebug(1,"Stop."); | |
108 | }//Reconstruct(for simulated digits) | |
109 | //__________________________________________________________________________________________________ | |
110 | void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const | |
111 | { | |
112 | //Invoked by AliReconstruction to convert raw digits from DDL files to clusters | |
113 | //Arguments: pAL - ALICE run loader pointer | |
114 | // pRR - ALICE raw reader pointer | |
115 | // Returns: none | |
116 | AliLoader *pRL=pAL->GetDetectorLoader("HMPID"); AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader | |
117 | ||
118 | AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it | |
119 | ||
3c6274c1 | 120 | 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 |
121 | ||
d3da6dc4 | 122 | Int_t iEvtN=0; |
123 | while(pRR->NextEvent()){//events loop | |
124 | pAL->GetEvent(iEvtN++); | |
125 | pRL->MakeTree("R"); pRich->MakeBranch("R"); | |
126 | ||
127 | for(Int_t iCh=0;iCh<7;iCh++){//chambers loop | |
d3da6dc4 | 128 | pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber |
129 | UInt_t w32=0; | |
130 | while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files) | |
131 | UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13 | |
da08475b | 132 | dig.Raw(ddl,w32); |
133 | 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 | 134 | new((*((TClonesArray*)digLst.At(iCh)))[iDigCnt[iCh]++]) AliHMPIDDigit(dig); //add this digit to the tmp list |
d3da6dc4 | 135 | }//raw records loop |
d3da6dc4 | 136 | pRR->Reset(); |
d3da6dc4 | 137 | }//chambers loop |
3c6274c1 | 138 | |
139 | Dig2Clu(&digLst,pRich->CluLst());//cluster finder for all chambers | |
140 | for(Int_t i=0;i<7;i++){digLst.At(i)->Delete(); iDigCnt[i]=0;} //clean up list of digits for all chambers | |
141 | ||
d3da6dc4 | 142 | pRL->TreeR()->Fill(); //fill tree for current event |
143 | pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event | |
144 | pRich->CluReset(); | |
145 | }//events loop | |
146 | pRL->UnloadRecPoints(); | |
147 | }//Reconstruct raw data | |
148 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
149 | ||
150 | void AliHMPIDReconstructor::FillESD(AliRunLoader *, AliESD *pESD) const | |
151 | { | |
152 | // Calculates probability to be a electron-muon-pion-kaon-proton | |
153 | // from the given Cerenkov angle and momentum assuming no initial particle composition | |
154 | // (i.e. apriory probability to be the particle of the given sort is the same for all sorts) | |
155 | ||
156 | AliPID ppp; //needed | |
157 | Double_t pid[AliPID::kSPECIES],h[AliPID::kSPECIES]; | |
158 | ||
159 | for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop | |
160 | AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track | |
161 | if(pTrk->GetHMPIDsignal()<=0){//HMPID does not find anything reasonable for this track, assign 0.2 for all species | |
162 | for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) pid[iPart]=1.0/AliPID::kSPECIES; | |
163 | pTrk->SetHMPIDpid(pid); | |
164 | continue; | |
165 | } | |
166 | Double_t pmod = pTrk->GetP(); | |
167 | Double_t hTot=0; | |
168 | for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){ | |
169 | Double_t mass = AliPID::ParticleMass(iPart); | |
170 | Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(AliHMPIDParam::Instance()->MeanIdxRad()*pmod); | |
171 | if(cosThetaTh<1) //calculate the height of theortical theta ckov on the gaus of experimental one | |
172 | h[iPart] =TMath::Gaus(TMath::ACos(cosThetaTh),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE); | |
173 | ||
174 | else //beta < 1/ref. idx. => no light at all | |
175 | h[iPart] =0 ; | |
176 | hTot +=h[iPart]; //total height of all theoretical heights for normalization | |
177 | }//species loop | |
178 | ||
179 | Double_t hMin=TMath::Gaus(pTrk->GetHMPIDsignal()-4*TMath::Sqrt(pTrk->GetHMPIDchi2()),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);//5 sigma protection | |
180 | ||
181 | for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)//species loop to assign probabilities | |
182 | if(hTot>hMin) | |
183 | pid[iPart]=h[iPart]/hTot; | |
184 | else //all theoretical values are far away from experemental one | |
185 | pid[iPart]=1.0/AliPID::kSPECIES; | |
186 | pTrk->SetHMPIDpid(pid); | |
187 | }//ESD tracks loop | |
188 | //last line is to check if the nearest thetacerenkov to the teorethical one is within 5 sigma, otherwise no response (equal prob to every particle | |
189 | }//FillESD | |
190 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |