]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDReconstructor.cxx
Changed selector to fit with the new AlidNdEtaCorrection scheme.
[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 #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()
24 #include <TH1F.h>                 //CluQA() 
25 #include <TH2F.h>                 //CluQA() 
26 #include <TCanvas.h>              //CluQA() 
27 #include <TNtupleD.h>             //CheckPR()
28 ClassImp(AliHMPIDReconstructor)
29
30 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 void AliHMPIDReconstructor::CluQA(AliRunLoader *pAL)
32 {
33 // Quality assesment plots for clusters. 
34 // This methode takes list of digits and form list of clusters again in order to 
35 // calculate cluster shape and cluster particle mixture    
36   AliLoader *pRL=pAL->GetDetectorLoader("HMPID");  AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
37   Int_t iNevt=pAL->GetNumberOfEvents();  if(iNevt==0)             {AliInfoClass("No events");return;}   
38                                          if(pRL->LoadDigits())    {AliInfoClass("No digits file");return;}
39                                             pAL->LoadHeader();
40                                             pAL->LoadKinematics(); 
41 //                                            AliStack *pStack=pAL->Stack();
42   TH1::AddDirectory(kFALSE);
43   
44         
45   TH1F*    pQ=new TH1F("RiAllQ"  ,"Charge All"           ,4000 ,0  ,4000);// Q hists
46   TH1F* pCerQ=new TH1F("RiCerQ"  ,"Charge Ckov"          ,4000 ,0  ,4000);
47   TH1F* pMipQ=new TH1F("RiMipQ"  ,"Charge MIP"           ,4000 ,0  ,4000);
48   
49   TH1F*    pS=new TH1F("RichCluSize"    ,"Cluster size;size"         ,100  ,0  ,100 );// size hists
50   TH1F* pCerS=new TH1F("RichCluCerSize" ,"Ckov size;size"            ,100  ,0  ,100 );
51   TH1F* pMipS=new TH1F("RichCluMipSize" ,"MIP size;size"             ,100  ,0  ,100 );
52   
53   TH2F*    pM=new TH2F("RichCluMap"     ,"Cluster map;x [cm];y [cm]" ,1000 ,0  ,AliHMPIDDigit::SizePcX(),1000,0,AliHMPIDDigit::SizePcY()); // maps
54   TH2F* pMipM=new TH2F("RichCluMipMap"  ,"MIP map;x [cm];y [cm]"     ,1000 ,0  ,AliHMPIDDigit::SizePcX(),1000,0,AliHMPIDDigit::SizePcY());
55   TH2F* pCerM=new TH2F("RichCluCerMap"  ,"Ckov map;x [cm];y [cm]"    ,1000 ,0  ,AliHMPIDDigit::SizePcX(),1000,0,AliHMPIDDigit::SizePcY());
56  
57   
58   
59   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){
60     pAL->GetEvent(iEvt);               
61     pRL->TreeD()->GetEntry(0); 
62     TClonesArray *pCluLst=new TClonesArray("AliHMPIDCluster");//tmp list of clusters for this event
63     
64     for(Int_t iCh=0;iCh<7;iCh++) Dig2Clu(pRich->DigLst(iCh),pCluLst,kFALSE);//cluster finder for all chamber if any digits present
65     
66     for(Int_t iClu=0;iClu<pCluLst->GetEntriesFast();iClu++){
67       AliHMPIDCluster *pClu = (AliHMPIDCluster*)pCluLst->At(iClu);
68       Int_t cfm=0; for(Int_t iDig=0;iDig<pClu->Size();iDig++)  cfm+=pClu->Dig(iDig)->Ch(); //collect ckov-fee-mip structure of current cluster ?????
69       Int_t iNckov=cfm/1000000;      Int_t iNfee =cfm%1000000/1000;      Int_t iNmip =cfm%1000000%1000; 
70
71                                              pQ   ->Fill(pClu->Q()) ; pS   ->Fill(pClu->Size()) ; pM    ->Fill(pClu->X(),pClu->Y()); //all clusters                                      
72       if(iNckov!=0 && iNfee==0 && iNmip==0) {pCerQ->Fill(pClu->Q()) ; pCerS->Fill(pClu->Size()) ; pCerM ->Fill(pClu->X(),pClu->Y());}//ckov only cluster
73       if(iNckov==0 && iNfee==0 && iNmip!=0) {pMipQ->Fill(pClu->Q()) ; pMipS->Fill(pClu->Size()) ; pMipM ->Fill(pClu->X(),pClu->Y());}//mip only cluster
74                                        
75     }//clusters loop   
76     pCluLst->Clear();delete pCluLst;
77   }//events loop
78   
79   pRL->UnloadDigits(); pAL->UnloadKinematics(); pAL->UnloadHeader();
80   TCanvas *pC=new TCanvas("RichCluQA",Form("QA for cluster from %i events",iNevt),1000,900); pC->Divide(3,3);
81   pC->cd(1);    pM->Draw();          pC->cd(2);    pQ->Draw();       pC->cd(3);    pS->Draw();        
82   pC->cd(4); pMipM->Draw();          pC->cd(5); pMipQ->Draw();       pC->cd(6); pMipS->Draw();        
83   pC->cd(7); pCerM->Draw();          pC->cd(8); pCerQ->Draw();       pC->cd(9); pCerS->Draw();        
84 }//CluQA()
85 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
86 void AliHMPIDReconstructor::Dig2Clu(TClonesArray *pDigLst,TClonesArray *pCluLst,Bool_t isTryUnfold)
87 {
88 // Finds all clusters for a given digits list provided not empty. Currently digits list is a list of all digits for a single chamber.
89 // Puts all found clusters in separate lists, one per clusters. 
90 // Arguments: pDigLst     - list of digits provided not empty
91 //            pCluLst     - list of clusters, provided empty     
92 //            isTryUnfold - flag to choose between CoG and Mathieson fitting  
93 //  Returns: none    
94   TMatrixF digMap(AliHMPIDDigit::kPadAllX,AliHMPIDDigit::kPadAllY);  digMap=(Float_t)-1;     //digit map for single chamber reseted to -1
95   for(Int_t iDig=0;iDig<pDigLst->GetEntriesFast();iDig++){                                 //digits loop to fill digits map
96     AliHMPIDDigit *pDig= (AliHMPIDDigit*)pDigLst->At(iDig);                                  //get current digit
97     digMap( pDig->PadX(), pDig->PadY() )=iDig;                                             //fill the map, (padx,pady) cell takes digit index
98   }                                                                                        //digits loop to fill digits map 
99   
100   AliHMPIDCluster clu;                                                                      //tmp cluster to be used as current
101   
102   for(Int_t iDig=0;iDig<pDigLst->GetEntriesFast();iDig++){                                 //digits loop to form clusters list
103     AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigLst->At(iDig);                                   //take current digit
104     if(!(pDig=UseDig(pDig->PadX(),pDig->PadY(),pDigLst,&digMap))) continue;                //this digit is already taken in FormClu(), go after next digit
105     FormClu(&clu,pDig,pDigLst,&digMap);                                                    //form cluster starting from this digit by recursion
106     clu.Solve(pCluLst,isTryUnfold);                                                        //solve this cluster and add all unfolded clusters to provided list  
107     clu.Reset();                                                                           //empty current cluster
108   }                                                                                        //digits loop to form clusters list
109 }//Dig2Clu()
110 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111 void  AliHMPIDReconstructor::FormClu(AliHMPIDCluster *pClu,AliHMPIDDigit *pDig,TClonesArray *pDigLst,TMatrixF *pDigMap)
112 {
113 // Forms the initial cluster as a sum of all adjascent digits. Starts from the given digit
114 // then calls itself recursevly  for all neighbours.
115 // Arguments: pClu - pointer to cluster being formed
116 //  Returns: none
117   pClu->DigAdd(pDig);//take this digit in cluster
118
119   Int_t x[4],y[4];
120   
121   
122   Int_t iPadCnt=0;  Int_t iPadX=pDig->PadX(); Int_t iPadY=pDig->PadY();
123   if(iPadX != AliHMPIDDigit::kPad1)     {x[iPadCnt]=iPadX-1; y[iPadCnt]=iPadY;   iPadCnt++;}       //left
124   if(iPadX != AliHMPIDDigit::kPadPcX)   {x[iPadCnt]=iPadX+1; y[iPadCnt]=iPadY;   iPadCnt++;}       //right
125   if(iPadY != AliHMPIDDigit::kPad1)     {x[iPadCnt]=iPadX;   y[iPadCnt]=iPadY-1; iPadCnt++;}       //down
126   if(iPadY != AliHMPIDDigit::kPadPcY)   {x[iPadCnt]=iPadX;   y[iPadCnt]=iPadY+1; iPadCnt++;}       //up
127   
128   for (Int_t i=0;i<iPadCnt;i++)
129     if((pDig=UseDig(x[i],y[i],pDigLst,pDigMap))) FormClu(pClu,pDig,pDigLst,pDigMap);   //check if this neighbour hit and mark it as taken  
130 }//FormClu()
131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
132 void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL)const
133 {
134 //Invoked  by AliReconstruction to convert digits to clusters i.e. reconstruct simulated data
135 //Arguments: pAL - ALICE run loader pointer
136 //  Returns: none    
137   AliDebug(1,"Start.");
138   AliLoader *pRL=pAL->GetDetectorLoader("HMPID");
139   AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
140   pRL->LoadDigits();   
141   pRL->LoadRecPoints("recreate");
142   
143   for(Int_t iEvtN=0;iEvtN<pAL->GetNumberOfEvents();iEvtN++){//events loop
144     pAL->GetEvent(iEvtN); AliDebug(1,Form("Processing event %i...",iEvtN)); //switch current directory to next event    
145     pRL->TreeD()->GetEntry(0);  pRL->MakeTree("R");  pRich->MakeBranch("R");  //load digits to memory  and create branches for clusters              
146     
147     for(Int_t iCh=0;iCh<7;iCh++) Dig2Clu(pRich->DigLst(iCh),pRich->CluLst(iCh));//cluster finder 
148       
149     pRL->TreeR()->Fill();            //fill tree for current event
150     pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
151     pRich->DigReset(); pRich->CluReset();
152   }//events loop  
153
154   pRL->UnloadDigits(); 
155   pRL->UnloadRecPoints();  
156     
157   AliDebug(1,"Stop.");      
158 }//Reconstruct(for simulated digits)
159 //__________________________________________________________________________________________________
160 void AliHMPIDReconstructor::Reconstruct(AliRunLoader *pAL,AliRawReader* pRR)const
161 {
162 //Invoked  by AliReconstruction to convert raw digits from DDL files to clusters
163 //Arguments: pAL - ALICE run loader pointer
164 //           pRR - ALICE raw reader pointer  
165 //  Returns: none    
166   AliLoader *pRL=pAL->GetDetectorLoader("HMPID");  AliHMPID *pRich=(AliHMPID*)pAL->GetAliRun()->GetDetector("HMPID");//get pointers for HMPID and HMPID loader
167   
168   AliHMPIDDigit dig; //tmp digit, raw digit will be converted to it
169   
170   Int_t iEvtN=0;
171   while(pRR->NextEvent()){//events loop
172     pAL->GetEvent(iEvtN++);
173     pRL->MakeTree("R");  pRich->MakeBranch("R");
174     
175     for(Int_t iCh=0;iCh<7;iCh++){//chambers loop
176       TClonesArray *pDigLst=new TClonesArray("AliHMPIDDigit"); Int_t iDigCnt=0; //tmp list of digits for single chamber
177       pRR->Select("HMPID",2*iCh,2*iCh+1);//select only DDL files for the current chamber      
178       UInt_t w32=0;
179       while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
180         UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
181         dig.ReadRaw(ddl,w32);  
182         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.PadX(),dig.PadY(),dig.Q()));
183         new((*pDigLst)[iDigCnt++]) AliHMPIDDigit(dig); //add this digit to the tmp list
184       }//raw records loop
185       if(iDigCnt) Dig2Clu(pDigLst,pRich->CluLst(iCh));//cluster finder for the current chamber if any digits present
186       pRR->Reset();        
187       pDigLst->Delete();  iDigCnt=0;//clean up list of digits for the current chamber
188     }//chambers loop
189     pRL->TreeR()->Fill();            //fill tree for current event
190     pRL->WriteRecPoints("OVERWRITE");//write out clusters for current event
191     pRich->CluReset();
192   }//events loop  
193   pRL->UnloadRecPoints();  
194 }//Reconstruct raw data
195 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
196
197 void AliHMPIDReconstructor::FillESD(AliRunLoader *, AliESD *pESD) const
198 {
199 // Calculates probability to be a electron-muon-pion-kaon-proton
200 // from the given Cerenkov angle and momentum assuming no initial particle composition
201 // (i.e. apriory probability to be the particle of the given sort is the same for all sorts)
202
203   AliPID ppp; //needed
204   Double_t pid[AliPID::kSPECIES],h[AliPID::kSPECIES];
205    
206   for(Int_t iTrk=0;iTrk<pESD->GetNumberOfTracks();iTrk++){//ESD tracks loop
207     AliESDtrack *pTrk = pESD->GetTrack(iTrk);// get next reconstructed track
208     if(pTrk->GetHMPIDsignal()<=0){//HMPID does not find anything reasonable for this track, assign 0.2 for all species
209       for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++) pid[iPart]=1.0/AliPID::kSPECIES;
210       pTrk->SetHMPIDpid(pid);
211       continue;
212     } 
213     Double_t pmod = pTrk->GetP();
214     Double_t hTot=0;
215     for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++){
216       Double_t mass = AliPID::ParticleMass(iPart);
217       Double_t cosThetaTh = TMath::Sqrt(mass*mass+pmod*pmod)/(AliHMPIDParam::Instance()->MeanIdxRad()*pmod);
218       if(cosThetaTh<1) //calculate the height of theortical theta ckov on the gaus of experimental one
219         h[iPart] =TMath::Gaus(TMath::ACos(cosThetaTh),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);
220       
221       else             //beta < 1/ref. idx. => no light at all  
222         h[iPart] =0 ;       
223       hTot    +=h[iPart]; //total height of all theoretical heights for normalization
224     }//species loop
225      
226     Double_t hMin=TMath::Gaus(pTrk->GetHMPIDsignal()-4*TMath::Sqrt(pTrk->GetHMPIDchi2()),pTrk->GetHMPIDsignal(),TMath::Sqrt(pTrk->GetHMPIDchi2()),kTRUE);//5 sigma protection
227     
228     for(Int_t iPart=0;iPart<AliPID::kSPECIES;iPart++)//species loop to assign probabilities
229       if(hTot>hMin)  
230         pid[iPart]=h[iPart]/hTot;
231       else                               //all theoretical values are far away from experemental one
232         pid[iPart]=1.0/AliPID::kSPECIES; 
233     pTrk->SetHMPIDpid(pid);
234   }//ESD tracks loop
235   //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
236 }//FillESD
237 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++