]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliAnalysisTaskTrigChEff.cxx
- Clean methods related to matching and hit pattern in AODTrack.
[u/mrichter/AliRoot.git] / PWG3 / muon / AliAnalysisTaskTrigChEff.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 //-----------------------------------------------------------------------------
17 /// \class AliAnalysisTaskSingleMu
18 /// Analysis task for single muons in the spectrometer.
19 /// The output is a list of histograms.
20 /// The macro class can run on AOD or ESDs.
21 /// If Monte Carlo information is present, some basics checks are performed.
22 ///
23 /// \author Diego Stocco
24 //-----------------------------------------------------------------------------
25
26 //----------------------------------------------------------------------------
27 //    Implementation of the class for trigger chamber efficiency determinaltion
28 //----------------------------------------------------------------------------
29
30
31 #define AliAnalysisTaskTrigChEff_cxx
32
33 // ROOT includes
34 #include "TH1.h"
35 #include "TCanvas.h"
36 #include "TROOT.h"
37 #include "TString.h"
38 #include "TList.h"
39
40 // STEER includes
41 #include "AliLog.h"
42 #include "AliESDEvent.h"
43 #include "AliESDMuonTrack.h"
44
45 // ANALYSIS includes
46 #include "AliAnalysisTaskSE.h"
47
48 #include "AliAnalysisTaskTrigChEff.h"
49
50 ClassImp(AliAnalysisTaskTrigChEff)
51
52 //________________________________________________________________________
53 AliAnalysisTaskTrigChEff::AliAnalysisTaskTrigChEff(const char *name) :
54   AliAnalysisTaskSE(name), 
55   fUseGhosts(kFALSE),
56   fList(0)
57 {
58   //
59   /// Constructor.
60   //
61   // Output slot #1 writes into a TObjArray container
62   DefineOutput(1,  TList::Class());
63 }
64
65 //___________________________________________________________________________
66 void AliAnalysisTaskTrigChEff::UserCreateOutputObjects() {
67   //
68   /// Create histograms
69   /// Called once
70   //
71
72   TString cathCode[2] = {"bendPlane", "nonBendPlane"};
73   TString countTypeName[2] = {"CountInCh", "NonCountInCh"};
74
75   const Char_t* yAxisTitle = "counts";
76
77   const Int_t kNboards = 234; //AliMpConstants::NofLocalBoards();
78   const Int_t kFirstTrigCh = 11;//AliMpConstants::NofTrackingChambers()+1;
79
80   Int_t chamberBins = kNchambers;
81   Float_t chamberLow = kFirstTrigCh-0.5, chamberHigh = kFirstTrigCh+kNchambers-0.5;
82   const Char_t* chamberName = "chamber";
83
84   Int_t slatBins = kNslats;
85   Float_t slatLow = 0-0.5, slatHigh = kNslats-0.5;
86   const Char_t* slatName = "slat";
87
88   Int_t boardBins = kNboards;
89   Float_t boardLow = 1-0.5, boardHigh = kNboards+1.-0.5;
90   const Char_t* boardName = "board";
91
92   Int_t angleBins = 280;
93   Float_t angleLow = -70., angleHigh = 70.;
94   const Char_t* angleNameX = "#theta_{x} (deg)";
95   const Char_t* angleNameY = "#theta_{y} (deg)";
96
97   TString baseName, histoName;
98   fList = new TList();
99
100   TH1F* histo;
101
102   histo = new TH1F("nTracksInSlat", "Num. of tracks used for efficiency calculation", 
103                    slatBins, slatLow, slatHigh);
104   histo->GetXaxis()->SetTitle(slatName);
105   histo->GetYaxis()->SetTitle("num of used tracks");
106
107   fList->AddAt(histo, kHtracksInSlat);
108
109   histo = new TH1F("nTracksInBoard", "Num. of tracks used for efficiency calculation", 
110                    boardBins, boardLow, boardHigh);
111   histo->GetXaxis()->SetTitle(boardName);
112   histo->GetYaxis()->SetTitle("num of used tracks");
113
114   fList->AddAt(histo, kHtracksInBoard);
115
116   for(Int_t hType=0; hType<kNcounts; hType++){
117     Int_t hindex = (hType==0) ? kHchamberEff : kHchamberNonEff;
118     for(Int_t cath=0; cath<kNcathodes; cath++){
119       histoName = Form("%sChamber%s", cathCode[cath].Data(), countTypeName[hType].Data());
120       histo = new TH1F(histoName, histoName,
121                        chamberBins, chamberLow, chamberHigh);
122       histo->GetXaxis()->SetTitle(chamberName);
123       histo->GetYaxis()->SetTitle(yAxisTitle);
124         
125       fList->AddAt(histo, hindex + cath);
126     } // loop on cath
127   } // loop on counts
128
129   for(Int_t hType=0; hType<kNcounts; hType++){
130     Int_t hindex = (hType==0) ? kHslatEff : kHslatNonEff;
131     for(Int_t cath=0; cath<kNcathodes; cath++){
132       for(Int_t ch=0; ch<kNchambers; ch++){
133         Int_t chCath = GetPlane(cath, ch);
134         histoName = Form("%sSlat%s%i", cathCode[cath].Data(), countTypeName[hType].Data(), kFirstTrigCh+ch);
135         histo = new TH1F(histoName, histoName,
136                          slatBins, slatLow, slatHigh);
137         histo->GetXaxis()->SetTitle(slatName);
138         histo->GetYaxis()->SetTitle(yAxisTitle);
139
140         fList->AddAt(histo, hindex + chCath);
141       } // loop on chamber
142     } // loop on cath
143   } // loop on counts
144
145   for(Int_t hType=0; hType<kNcounts; hType++){
146     Int_t hindex = (hType==0) ? kHboardEff : kHboardNonEff;
147     for(Int_t cath=0; cath<kNcathodes; cath++){
148       for(Int_t ch=0; ch<kNchambers; ch++){
149         Int_t chCath = GetPlane(cath, ch);
150         histoName = Form("%sBoard%s%i", cathCode[cath].Data(), countTypeName[hType].Data(), kFirstTrigCh+ch);
151         histo = new TH1F(histoName, histoName,
152                          boardBins, boardLow, boardHigh);
153         histo->GetXaxis()->SetTitle(boardName);
154         histo->GetYaxis()->SetTitle(yAxisTitle);
155
156         fList->AddAt(histo, hindex + chCath);
157       } // loop on chamber
158     } // loop on cath
159   } // loop on counts
160
161   histo = new TH1F("thetaX", "Angular distribution",
162                    angleBins, angleLow, angleHigh);
163   histo->GetXaxis()->SetTitle(angleNameX);
164   histo->GetYaxis()->SetTitle("entries");
165   fList->AddAt(histo, kHthetaX);
166
167   histo = new TH1F("thetaY", "Angular distribution",
168                    angleBins, angleLow, angleHigh);
169   histo->GetXaxis()->SetTitle(angleNameY);
170   histo->GetYaxis()->SetTitle("entries");
171   fList->AddAt(histo, kHthetaY);
172
173 }
174
175 //________________________________________________________________________
176 void AliAnalysisTaskTrigChEff::UserExec(Option_t *) {
177   //
178   /// Main loop
179   /// Called for each event
180   //
181   AliESDEvent* esdEvent = dynamic_cast<AliESDEvent*> (InputEvent());
182
183   if (!esdEvent) {
184     Printf("ERROR: esdEvent not available\n");
185     return;
186   }
187
188   Int_t slat = 0, board = 0;
189   UShort_t pattern = 0;
190   AliESDMuonTrack *esdTrack = 0x0;
191
192   const Float_t kRadToDeg = 180./TMath::Pi();
193   Int_t nTracks = esdEvent->GetNumberOfMuonTracks();
194
195   const Int_t kFirstTrigCh = 11; //AliMpConstants::NofTrackingChambers()+1;
196
197   TArrayI othersEfficient(kNchambers);
198
199   for (Int_t itrack = 0; itrack < nTracks; itrack++) {
200     esdTrack = esdEvent->GetMuonTrack(itrack);
201
202     if ( ! esdTrack->ContainTrackerData() && ! fUseGhosts ) continue;
203
204     pattern =  esdTrack->GetHitsPatternInTrigCh();
205     Int_t effFlag = AliESDMuonTrack::GetEffFlag(pattern);
206
207     if(effFlag < AliESDMuonTrack::kChEff) continue; // Track not good for efficiency calculation
208
209     ((TH1F*)fList->At(kHthetaX))->Fill(esdTrack->GetThetaX() * kRadToDeg);
210     ((TH1F*)fList->At(kHthetaY))->Fill(esdTrack->GetThetaY() * kRadToDeg);
211
212     othersEfficient.Reset(1);
213     for(Int_t cath=0; cath<kNcathodes; cath++){
214       for(Int_t ich=0; ich<kNchambers; ich++){
215         if( ! AliESDMuonTrack::IsChamberHit(pattern, cath, ich)){
216           for(Int_t jch=0; jch<kNchambers; jch++){
217             if ( jch != ich) {
218               othersEfficient[jch] = 0;
219               //AliInfo(Form("%s ch %i by New", baseOutString.Data(), jch));
220             }
221           } // loop on other chambers
222           break;
223         } // if chamber not efficient
224       } // loop on chambers
225     } // loop on cathodes
226
227     Bool_t rejectTrack = kTRUE;
228     for (Int_t ich=0; ich<kNchambers; ich++){
229       if ( othersEfficient[ich] > 0 ){
230         rejectTrack = kFALSE;
231         break;
232       }
233     }
234
235     if ( rejectTrack ) continue;
236
237     slat = AliESDMuonTrack::GetSlatOrInfo(pattern);
238     board = esdTrack->LoCircuit();
239
240     if(effFlag >= AliESDMuonTrack::kSlatEff) ((TH1F*)fList->At(kHtracksInSlat))->Fill(slat);
241     if(effFlag >= AliESDMuonTrack::kBoardEff) ((TH1F*)fList->At(kHtracksInBoard))->Fill(board);
242
243     for(Int_t cath=0; cath<kNcathodes; cath++){
244       for(Int_t ch=0; ch<kNchambers; ch++){
245         if ( ! othersEfficient[ch] )
246           continue; // Reject track if the info of the chamber under study 
247                     // is necessary to create the track itself
248
249         Int_t whichType = AliESDMuonTrack::IsChamberHit(pattern, cath, ch) ? kChHit : kChNonHit;
250
251         Int_t iChamber = kFirstTrigCh + ch;
252         Int_t hindex = ( whichType == kChHit ) ? kHchamberEff : kHchamberNonEff;
253         ((TH1F*)fList->At(hindex + cath))->Fill(iChamber);
254
255         if(effFlag < AliESDMuonTrack::kSlatEff) continue; // Track crossed different slats
256         Int_t chCath = GetPlane(cath, ch);
257         hindex = ( whichType == kChHit ) ? kHslatEff : kHslatNonEff;
258         ((TH1F*)fList->At(hindex + chCath))->Fill(slat);
259
260         if(effFlag < AliESDMuonTrack::kBoardEff) continue; // Track crossed different boards
261         hindex = ( whichType == kChHit ) ? kHboardEff : kHboardNonEff;
262         ((TH1F*)fList->At(hindex + chCath))->Fill(board);
263       } // loop on chambers
264     } // loop on cathodes
265   } // loop on tracks
266
267   // Post final data. It will be written to a file with option "RECREATE"
268   PostData(1, fList);
269 }
270
271 //________________________________________________________________________
272 void AliAnalysisTaskTrigChEff::Terminate(Option_t *) {
273   //
274   /// Draw result to the screen
275   /// Called once at the end of the query.
276   //
277   if (!gROOT->IsBatch()) {
278     TCanvas *can[kNcathodes];
279     TH1F *num = 0x0;
280     TH1F *den = 0x0;
281     for(Int_t cath=0; cath<kNcathodes; cath++){
282       TString canName = Form("can%i",cath);
283       can[cath] = new TCanvas(canName.Data(),canName.Data(),10*(1+cath),10*(1+cath),310,310);
284       can[cath]->SetFillColor(10); can[cath]->SetHighLightColor(10);
285       can[cath]->SetLeftMargin(0.15); can[cath]->SetBottomMargin(0.15);  
286       can[cath]->Divide(2,2);
287       for(Int_t ch=0; ch<kNchambers; ch++){
288         Int_t chCath = GetPlane(cath, ch);
289         num = (TH1F*)(fList->At(kHboardEff + chCath)->Clone());
290         den = (TH1F*)(fList->At(kHboardNonEff + chCath)->Clone());
291         den->Add(num);
292         num->Divide(den);
293         can[cath]->cd(ch+1);
294         num->DrawCopy("E");
295       }
296     }
297   }
298 }