3 #include <TCanvas.h> //Terminate()
5 #include <TBenchmark.h>
6 #include <TFile.h> //docosmic()
7 #include <AliSelector.h> //base class
9 #include <AliBitPacking.h> //HmpidPayload()
10 #include "AliHMPIDDigit.h"
11 #include "AliHMPIDCluster.h"
12 #include "AliHMPIDReconstructor.h" //docosmic()
14 class AliHMPIDSelector : public AliSelector {
16 AliHMPIDSelector():AliSelector(),fChain(0),fEsd(0),fCkovP(0),fMipXY(0),fDifX(0),fSigP(0) {for(Int_t i=0;i<5;i++) fProb[i]=0;}
17 virtual ~AliHMPIDSelector() {delete fEsd;}
20 virtual Int_t Version () const {return 2;}
21 virtual void Begin (TTree *) {}
22 virtual void SlaveBegin (TTree *tree);
23 virtual void Init (TTree *tree);
24 virtual Bool_t Notify () {return kTRUE;}
25 virtual Bool_t Process (Long64_t entry);
26 virtual void SetOption (const char *option) { fOption = option; }
27 virtual void SetObject (TObject *obj) { fObject = obj; }
28 virtual void SetInputList (TList *input) {fInput = input;}
29 virtual TList *GetOutputList () const { return fOutput; }
30 virtual void SlaveTerminate ();
31 virtual void Terminate ();
34 TTree *fChain ; //!pointer to the analyzed TTree or TChain
37 TH2F *fCkovP,*fMipXY; //!
42 ClassDef(AliHMPIDSelector,0);
46 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 void AliHMPIDSelector::SlaveBegin(TTree *tree)
49 // The SlaveBegin() function is called after the Begin() function. When running with PROOF SlaveBegin() is called on each slave server.
50 // The tree argument is deprecated (on PROOF 0 is passed).
54 TString option = GetOption();
56 // create histograms on each slave server
57 fCkovP = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]", 150, 0, 7 ,500, 0, 1);
58 fSigP = new TH2F("SigP" ,"#sigma_{#theta_c}" , 150, 0, 7 ,100, 0, 1e20);
59 fMipXY = new TH2F("MipXY" ,"mip position" , 260, 0,130 ,252,0,126);
60 fDifX = new TH1F("DifX" ,"diff" , 200, -5, 5);
62 fProb[0] = new TH1F("PidE" ,"PID: e yellow #mu magenta" ,100,0,1); fProb[0]->SetLineColor(kYellow);
63 fProb[1] = new TH1F("PidMu","pid of #mu" ,100,0,1); fProb[1]->SetLineColor(kMagenta);
64 fProb[2] = new TH1F("PidPi","PID: #pi red K green p blue",100,0,1); fProb[2]->SetLineColor(kRed);
65 fProb[3] = new TH1F("PidK" ,"pid of K" ,100,0,1); fProb[3]->SetLineColor(kGreen);
66 fProb[4] = new TH1F("PidP" ,"pid of p" ,100,0,1); fProb[4]->SetLineColor(kBlue);
68 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
69 void AliHMPIDSelector::Init(TTree *pTr)
71 // Called when the selector needs to initialize a new tree of chain. Typically here the branch addresses of the tree is set
72 // Init() will be called many times when running with PROOF.
76 fChain->SetBranchAddress("ESD", &fEsd) ;
77 fChain->SetBranchStatus("*", 0);
78 fChain->SetBranchStatus("fTracks.*", 1);
80 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
81 Bool_t AliHMPIDSelector::Process(Long64_t entry)
83 fChain->GetTree()->GetEntry(entry);
85 for(Int_t iTrk=0;iTrk<fEsd->GetNumberOfTracks();iTrk++){
86 AliESDtrack *pTrk=fEsd->GetTrack(iTrk);
88 // if(pTrk->GetHMPIDsignal()==-1) continue;
90 fCkovP->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal()) ;
91 fSigP ->Fill(pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
93 Float_t xm,ym; Int_t q,np; pTrk->GetHMPIDmip(xm,ym,q,np); fMipXY->Fill(xm,ym); //mip info
94 Float_t xd,yd,th,ph; pTrk->GetHMPIDtrk(xd,yd,th,ph); fDifX->Fill(xd-xm); //track info
96 Double_t pid[5]; pTrk->GetHMPIDpid(pid); for(Int_t i =0;i<5;i++) fProb[i]->Fill(pid[i]);
101 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
102 void AliHMPIDSelector::SlaveTerminate()
104 // The SlaveTerminate() function is called after all entries or objects
105 // have been processed. When running with PROOF SlaveTerminate() is called
106 // on each slave server.
108 // Add the histograms to the output on each slave server
110 fOutput->Add(fCkovP);
112 fOutput->Add(fMipXY);
115 for(Int_t i=0;i<5;i++) fOutput->Add(fProb[i]);
117 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
118 void AliHMPIDSelector::Terminate()
120 // The last function to be called. It always runs on the client, it can be used to present
121 // the results graphically or save the results to file.
122 fCkovP = dynamic_cast<TH2F*>(fOutput->FindObject("CkovP")) ;
123 fSigP = dynamic_cast<TH2F*>(fOutput->FindObject("SigP")) ;
124 fMipXY = dynamic_cast<TH2F*>(fOutput->FindObject("MipXY")) ;
125 fDifX = dynamic_cast<TH1F*>(fOutput->FindObject("DifX")) ;
127 fProb[0] = dynamic_cast<TH1F*>(fOutput->FindObject("PidE")) ;
128 fProb[1] = dynamic_cast<TH1F*>(fOutput->FindObject("PidMu")) ;
129 fProb[2] = dynamic_cast<TH1F*>(fOutput->FindObject("PidPi")) ;
130 fProb[3] = dynamic_cast<TH1F*>(fOutput->FindObject("PidK")) ;
131 fProb[4] = dynamic_cast<TH1F*>(fOutput->FindObject("PidP")) ;
133 Float_t n=1.292; //mean freon ref idx
134 TF1 *pPi=new TF1("RiPiTheo","acos(sqrt(x*x+[0]*[0])/(x*[1]))",1.2,7); pPi->SetLineWidth(1); pPi->SetParameter(1,n);
135 AliPID ppp; pPi->SetLineColor(kRed); pPi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion)); //mass
136 TF1 *pK=(TF1*)pPi->Clone(); pK ->SetLineColor(kGreen); pK ->SetParameter(0,AliPID::ParticleMass(AliPID::kKaon));
137 TF1 *pP=(TF1*)pPi->Clone(); pP ->SetLineColor(kBlue); pP ->SetParameter(0,AliPID::ParticleMass(AliPID::kProton));
139 TCanvas *pC=new TCanvas("c1","ESD QA");pC->SetFillColor(10); pC->SetHighLightColor(10); pC->Divide(3,2);
140 pC->cd(1); fCkovP->Draw(); pPi->Draw("same"); pK->Draw("same"); pP->Draw("same"); pC->cd(2); fMipXY->Draw(); pC->cd(3); fProb[0]->Draw(); fProb[1]->Draw("same");
141 pC->cd(4); fSigP ->Draw(); pC->cd(5); fDifX->Draw(); pC->cd(6); fProb[2]->Draw(); fProb[3]->Draw("same"); fProb[4]->Draw("same");
143 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
146 TChain* pChain =new TChain("esdTree"); pChain->Add("AliESDs.root");
148 pChain->Process("AliHMPIDSelector.C+");
150 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
151 Int_t DateHeader(ifstream *pFile,Bool_t isPrint=0)
154 pFile->read((char*)&iSize,4);
156 pFile->seekg(16*4,ios::cur);
160 Printf("Event size %i bytes",iSize); //1 common DATE header 17 words
161 pFile->read((char*)&w32,4); Printf("Event magic 0x%x" ,w32); //2
162 pFile->read((char*)&w32,4); Printf("Event head size %i bytes",w32); //3
163 pFile->read((char*)&w32,4); Printf("Event version 0x%x" ,w32); //4
164 pFile->read((char*)&w32,4); Printf("Event type %i (%s)" ,w32,(w32==7)? "physics":"SOR"); //5
165 pFile->read((char*)&w32,4); Printf("Run number %i" ,w32); //6
167 pFile->read((char*)&w32,4); Printf("Event ID 1 %i" ,w32); //7
168 pFile->read((char*)&w32,4); Printf("Event ID 2 %i" ,w32); //8
170 pFile->read((char*)&w32,4); Printf("Trigger pattern 1 %i" ,w32); //9
171 pFile->read((char*)&w32,4); Printf("Trigger pattern 2 %i" ,w32); //10
173 pFile->read((char*)&w32,4); Printf("Detector pattern %i" ,w32); //11
175 pFile->read((char*)&w32,4); Printf("Type attribute 1 %i" ,w32); //12
176 pFile->read((char*)&w32,4); Printf("Type attribute 2 %i" ,w32); //13
177 pFile->read((char*)&w32,4); Printf("Type attribute 3 %i" ,w32); //14
179 pFile->read((char*)&w32,4); Printf("LDC ID %i" ,w32); //15
180 pFile->read((char*)&w32,4); Printf("GDC ID %i" ,w32); //16
181 pFile->read((char*)&w32,4); TDatime time(w32); time.Print(); //17
187 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
188 void HmpidHeader(ifstream *pFile,Bool_t isPrint=kFALSE)
190 // Prints hmpid trailer and returns number of words for this trailer
191 if(!isPrint) {pFile->seekg(15*4,ios::cur);return;}
193 Printf("\nHMPID Header:");//private HEADER is 15 words
194 for(Int_t i=1;i<=11;i++) { pFile->read((char*)&w32,4); Printf("Word #%2i=%12i meaningless",i,w32);}
195 pFile->read((char*)&w32,4); Printf("Word #12=%12i event counter",w32);
196 for(Int_t i=13;i<=15;i++){ pFile->read((char*)&w32,4); Printf("Word #%2i=%12i empty",i,w32);}
198 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
199 void HmpidPayload(ifstream *pFile,Int_t iDdl,TObjArray *pDigAll)
201 // payload is 8 sequences with structure: WC36A8 then WC number of w32
203 TClonesArray *pDig1=(TClonesArray*)pDigAll->At(iDdl/2); //get list of digits for requested chamber
206 Int_t iDigCnt=pDig1->GetEntriesFast();
207 for(Int_t row=1;row<=8;row++){
208 pFile->read((char*)&w32,4); Int_t wc=AliBitPacking::UnpackWord(w32,16,31); Int_t ma=AliBitPacking::UnpackWord(w32, 0,15);
209 if(ma!=0x36a8) Printf("ERROR ERROR ERROR WRONG Marker=0x%x ERROR ERROR ERROR",ma);
210 for(Int_t i=1;i<=wc;i++){//words loop
211 pFile->read((char*)&w32,4);
212 if(w32&0x08000000) continue; //it's DILOGIC CW
213 AliHMPIDDigit *pDig=new AliHMPIDDigit;
215 new ((*pDig1)[iDigCnt++]) AliHMPIDDigit(*pDig);
219 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
220 void docosmic(const char* name,Int_t iMaxEvt=9999999)
222 gBenchmark->Start("Cosmic converter");
226 Bool_t isPrint=kFALSE;
228 TString rooName=name; rooName.Replace(1+rooName.First('.'),4,"root");
231 TFile *pOut=new TFile(rooName,"recreate");
232 TTree *pTr=new TTree("cosmic","some for time being");
235 TObjArray *pDigAll=new TObjArray(7); pDigAll->SetOwner(); for(Int_t i=0;i<7;i++) pDigAll->AddAt(new TClonesArray("AliHMPIDDigit") ,i); pTr->Branch("Digs",&pDigAll,64000,0);
236 TObjArray *pCluAll=new TObjArray(7); pCluAll->SetOwner(); for(Int_t i=0;i<7;i++) pCluAll->AddAt(new TClonesArray("AliHMPIDCluster"),i); pTr->Branch("Clus",&pCluAll,64000,0);
240 Int_t iSize=DateHeader(&in, isPrint); if(iSize==68) continue; //DATE header
242 HmpidHeader (&in, isPrint); //HMPID header
243 HmpidPayload(&in,ddl+1,pDigAll); //HMPID payload
244 HmpidHeader (&in, isPrint); //HMPID header
245 HmpidPayload(&in,ddl ,pDigAll); //HMPID payload
247 AliHMPIDReconstructor::Dig2Clu(pDigAll,pCluAll);
249 for(Int_t i=0;i<7;i++){
250 ((TClonesArray*)pDigAll->At(i))->Clear();
251 ((TClonesArray*)pCluAll->At(i))->Clear();
254 if(!(iEvt%200)) Printf("Event %i processed",iEvt);
256 if(iEvt==iMaxEvt) break;
259 pTr->Write(); pOut->Close();
260 pDigAll->Delete(); pCluAll->Delete();
263 Printf("Total %i events processed",iEvt);
264 gBenchmark->Show("Cosmic converter");
266 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
269 TChain *pCh=new TChain("cosmic"); pCh->Add("test.root");
272 TH1F *pDigQ=new TH1F("digQ","Digits QDC",500,0,4100);
273 TH1F *pDigO=new TH1F("digO","Digits # per event",500,0,8000);
274 TH2F *pDigM=new TH2F("digM","Digits map",500,0,131,500,0,127);
276 TH1F *pCluQ =new TH1F("cluQ","Clusters QDC",500,0,12100);
277 TH1F *pCluQmax=new TH1F("cluQmax","Clusters Max QDC",500,0,12100); pCluQmax->SetLineColor(kRed);
278 TH1F *pCluO=new TH1F("cluO","Clusters # per event",500,0,5000);
279 TH2F *pCluM=new TH2F("cluM","Clusters map",500,0,131,500,0,127);
281 TObjArray *pDigAll=0; pCh->SetBranchAddress("Digs",&pDigAll);
282 TObjArray *pCluAll=0; pCh->SetBranchAddress("Clus",&pCluAll);
285 for(Int_t iEvt=0;iEvt<pCh->GetEntries();iEvt++){
288 TClonesArray *pDigCh=(TClonesArray*)pDigAll->At(0);
289 TClonesArray *pCluCh=(TClonesArray*)pCluAll->At(0);
290 pDigO->Fill(pDigCh->GetEntriesFast());
291 pCluO->Fill(pCluCh->GetEntriesFast());
293 for(Int_t iDig=0;iDig<pDigCh->GetEntriesFast();iDig++){//digits loop
294 AliHMPIDDigit *pDig=(AliHMPIDDigit*)pDigCh->UncheckedAt(iDig);
295 pDigQ->Fill(pDig->Q());
296 pDigM->Fill(pDig->LorsX(),pDig->LorsY());
299 for(Int_t iClu=0;iClu<pCluCh->GetEntriesFast();iClu++){//clusters loop
300 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluCh->UncheckedAt(iClu);
301 pCluQ->Fill(pClu->Q());
302 if(pClu->Q()>qmax) qmax=pClu->Q();
303 pCluM->Fill(pClu->X(),pClu->Y());
305 pCluQmax->Fill(qmax);
308 TCanvas *pC=new TCanvas("comic","cosmic"); pC->Divide(2,3);
310 pC->cd(1); pDigM->Draw(); pC->cd(2); pCluM->Draw();
311 pC->cd(3); gPad->SetLogy(); pDigQ->Draw(); pC->cd(4); gPad->SetLogy(); pCluQ->Draw(); pCluQmax->Draw("same");
312 pC->cd(5); gPad->SetLogy(); pDigO->Draw(); pC->cd(6); gPad->SetLogy(); pCluO->Draw();
314 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++