]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFAcceptanceCuts.cxx
Fix for the loophole in the magnets currents check: the L3Off/DipON was passing the...
[u/mrichter/AliRoot.git] / CORRFW / AliCFAcceptanceCuts.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 ///////////////////////////////////////////////////////////////////////////
17 //          ----   CORRECTION FRAMEWORK   ----
18 // AliCFAcceptanceCuts implementation
19 // Class to cut on the number of AliTrackReference's 
20 // for each detector
21 ///////////////////////////////////////////////////////////////////////////
22 // author : R. Vernet (renaud.vernet@cern.ch)
23 ///////////////////////////////////////////////////////////////////////////
24
25 #include "AliLog.h"
26 #include "AliMCParticle.h"
27 #include "AliCFAcceptanceCuts.h"
28 #include "TH1F.h"
29 #include "TH2F.h"
30 #include "TList.h"
31 #include "TBits.h"
32
33 ClassImp(AliCFAcceptanceCuts)
34
35 //______________________________
36 AliCFAcceptanceCuts::AliCFAcceptanceCuts() : 
37   AliCFCutBase(),
38   fMCInfo(0x0),
39   fMinNHitITS(0),
40   fMinNHitTPC(0),
41   fMinNHitTRD(0),
42   fMinNHitTOF(0),
43   fMinNHitMUON(0),
44   fhCutStatistics(0x0),
45   fhCutCorrelation(0x0),
46   fBitmap(new TBits(0))
47 {
48   //
49   //ctor
50   //
51   for (int i=0; i<kNCuts; i++) for (int j=0; j<kNStepQA; j++) fhQA[i][j]=0x0;
52 }
53
54 //______________________________
55 AliCFAcceptanceCuts::AliCFAcceptanceCuts(const Char_t* name, const Char_t* title) : 
56   AliCFCutBase(name,title),
57   fMCInfo(0x0),
58   fMinNHitITS(0),
59   fMinNHitTPC(0),
60   fMinNHitTRD(0),
61   fMinNHitTOF(0),
62   fMinNHitMUON(0),
63   fhCutStatistics(0x0),
64   fhCutCorrelation(0x0),
65   fBitmap(new TBits(0))
66 {
67   //
68   //ctor
69   //
70   for (int i=0; i<kNCuts; i++) for (int j=0; j<kNStepQA; j++) fhQA[i][j]=0x0;
71 }
72
73 //______________________________
74 AliCFAcceptanceCuts::AliCFAcceptanceCuts(const AliCFAcceptanceCuts& c) : 
75   AliCFCutBase(c),
76   fMCInfo(c.fMCInfo),
77   fMinNHitITS(c.fMinNHitITS),
78   fMinNHitTPC(c.fMinNHitTPC),
79   fMinNHitTRD(c.fMinNHitTRD),
80   fMinNHitTOF(c.fMinNHitTOF),
81   fMinNHitMUON(c.fMinNHitMUON),
82   fhCutStatistics(c.fhCutStatistics),
83   fhCutCorrelation(c.fhCutCorrelation),
84   fBitmap(c.fBitmap)
85 {
86   //
87   //copy ctor
88   //
89   for (int i=0; i<kNCuts; i++) for (int j=0; j<kNStepQA; j++) fhQA[i][j]=c.fhQA[i][j];
90 }
91
92 //______________________________
93 AliCFAcceptanceCuts& AliCFAcceptanceCuts::operator=(const AliCFAcceptanceCuts& c)
94 {
95   //
96   // Assignment operator
97   //
98   if (this != &c) {
99     AliCFCutBase::operator=(c) ;
100     fMCInfo=c.fMCInfo;
101     fMinNHitITS=c.fMinNHitITS;
102     fMinNHitTPC=c.fMinNHitTPC;
103     fMinNHitTRD=c.fMinNHitTRD;
104     fMinNHitTOF=c.fMinNHitTOF;
105     fMinNHitMUON=c.fMinNHitMUON;
106     fhCutStatistics  = c.fhCutStatistics ;
107     fhCutCorrelation = c.fhCutCorrelation ;
108     fBitmap          = c.fBitmap ;
109     for (int i=0; i<kNCuts; i++) for (int j=0; j<kNStepQA; j++) fhQA[i][j]=c.fhQA[i][j];
110   }
111   return *this ;
112 }
113
114 //______________________________________________________________
115 Bool_t AliCFAcceptanceCuts::IsSelected(TObject *obj) {
116   //
117   // check if selections on 'obj' are passed
118   // 'obj' must be an AliMCParticle
119   //
120   
121   SelectionBitMap(obj);
122
123   if (fIsQAOn) FillHistograms(obj,kFALSE);
124
125   for (UInt_t icut=0; icut<fBitmap->GetNbits(); icut++)
126     if (!fBitmap->TestBitNumber(icut)) return kFALSE ;
127   
128   if (fIsQAOn) FillHistograms(obj,kTRUE);
129   return kTRUE;
130 }
131
132 //______________________________
133 void AliCFAcceptanceCuts::SelectionBitMap(TObject* obj) {
134   //
135   // checks the number of track references associated to 'obj'
136   // 'obj' must be an AliMCParticle
137   //
138
139   for (UInt_t i=0; i<kNCuts; i++) fBitmap->SetBitNumber(i,kFALSE);
140
141   if (!obj) return;
142   TString className(obj->ClassName());
143   if (className.CompareTo("AliMCParticle") != 0) {
144     AliError("obj must point to an AliMCParticle !");
145     return ;
146   }
147
148   AliMCParticle * part = dynamic_cast<AliMCParticle*>(obj);
149   if(!part) return ;
150
151   Int_t nHitsITS=0, nHitsTPC=0, nHitsTRD=0, nHitsTOF=0, nHitsMUON=0 ;
152   for (Int_t iTrackRef=0; iTrackRef<part->GetNumberOfTrackReferences(); iTrackRef++) {
153     AliTrackReference * trackRef = part->GetTrackReference(iTrackRef);
154     if(trackRef){
155       Int_t detectorId = trackRef->DetectorId();
156       switch(detectorId) {
157       case AliTrackReference::kITS  : nHitsITS++  ; break ;
158       case AliTrackReference::kTPC  : nHitsTPC++  ; break ;
159       case AliTrackReference::kTRD  : nHitsTRD++  ; break ;
160       case AliTrackReference::kTOF  : nHitsTOF++  ; break ;
161       case AliTrackReference::kMUON : nHitsMUON++ ; break ;
162       default : break ;
163       }
164     }
165   }
166   
167   Int_t iCutBit = 0;
168
169   if (nHitsITS  >= fMinNHitITS  ) fBitmap->SetBitNumber(iCutBit,kTRUE);
170   iCutBit++;
171   if (nHitsTPC  >= fMinNHitTPC  ) fBitmap->SetBitNumber(iCutBit,kTRUE);
172   iCutBit++;
173   if (nHitsTRD  >= fMinNHitTRD  ) fBitmap->SetBitNumber(iCutBit,kTRUE);
174   iCutBit++;
175   if (nHitsTOF  >= fMinNHitTOF  ) fBitmap->SetBitNumber(iCutBit,kTRUE);
176   iCutBit++;
177   if (nHitsMUON  >= fMinNHitMUON) fBitmap->SetBitNumber(iCutBit,kTRUE);
178 }
179
180
181 void AliCFAcceptanceCuts::SetMCEventInfo(const TObject* mcInfo) {
182   //
183   // Sets pointer to MC event information (AliMCEvent)
184   //
185
186   if (!mcInfo) {
187     AliError("Pointer to AliMCEvent !");
188     return;
189   }
190   
191   TString className(mcInfo->ClassName());
192   if (className.CompareTo("AliMCEvent") != 0) {
193     AliError("argument must point to an AliMCEvent !");
194     return ;
195   }
196   
197   fMCInfo = (AliMCEvent*) mcInfo ;
198 }
199
200
201 //__________________________________________________________________________________
202 void AliCFAcceptanceCuts::AddQAHistograms(TList *qaList) {
203   //
204   // saves the histograms in a TList
205   //
206
207   DefineHistograms();
208
209   qaList->Add(fhCutStatistics);
210   qaList->Add(fhCutCorrelation);
211
212   for (Int_t j=0; j<kNStepQA; j++) {
213     for(Int_t i=0; i<kNCuts; i++)
214       qaList->Add(fhQA[i][j]);
215   }
216 }
217
218 //__________________________________________________________________________________
219 void AliCFAcceptanceCuts::DefineHistograms() {
220   //
221   // histograms for cut variables, cut statistics and cut correlations
222   //
223   Int_t color = 2;
224
225   // book cut statistics and cut correlation histograms
226   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()),"",kNCuts,0.5,kNCuts+0.5);
227   fhCutStatistics->SetLineWidth(2);
228   int k = 1;
229   fhCutStatistics->GetXaxis()->SetBinLabel(k,"hits ITS") ; k++;
230   fhCutStatistics->GetXaxis()->SetBinLabel(k,"hits TPC") ; k++;
231   fhCutStatistics->GetXaxis()->SetBinLabel(k,"hits TRD") ; k++;
232   fhCutStatistics->GetXaxis()->SetBinLabel(k,"hits TOF") ; k++;
233   fhCutStatistics->GetXaxis()->SetBinLabel(k,"hits MUON"); k++;
234
235
236   fhCutCorrelation = new TH2F(Form("%s_cut_correlation",GetName()),"",kNCuts,0.5,kNCuts+0.5,kNCuts,0.5,kNCuts+0.5);
237   fhCutCorrelation->SetLineWidth(2);
238   for (k=1; k<=kNCuts; k++) {
239     fhCutCorrelation->GetXaxis()->SetBinLabel(k,fhCutStatistics->GetXaxis()->GetBinLabel(k));
240     fhCutCorrelation->GetYaxis()->SetBinLabel(k,fhCutStatistics->GetXaxis()->GetBinLabel(k));
241   }
242
243   Char_t str[256];
244   for (int i=0; i<kNStepQA; i++) {
245     if (i==0) sprintf(str," ");
246     else sprintf(str,"_cut");
247     fhQA[kCutHitsITS] [i] = new TH1F(Form("%s_HitsITS%s"  ,GetName(),str),"",10,0,10);
248     fhQA[kCutHitsTPC] [i] = new TH1F(Form("%s_HitsTPC%s"  ,GetName(),str),"",5,0,5);
249     fhQA[kCutHitsTRD] [i] = new TH1F(Form("%s_HitsTRD%s"  ,GetName(),str),"",20,0,20);
250     fhQA[kCutHitsTOF] [i] = new TH1F(Form("%s_HitsTOF%s"  ,GetName(),str),"",5,0,5);
251     fhQA[kCutHitsMUON][i] = new TH1F(Form("%s_HitsMUON%s" ,GetName(),str),"",5,0,5);
252   }
253   for(Int_t i=0; i<kNCuts; i++) fhQA[i][1]->SetLineColor(color);
254 }
255
256 //__________________________________________________________________________________
257 void AliCFAcceptanceCuts::FillHistograms(TObject* obj, Bool_t afterCuts)
258 {
259   //
260   // fill the QA histograms
261   //
262   AliMCParticle* part = dynamic_cast<AliMCParticle *>(obj);
263
264   Int_t nHitsITS=0, nHitsTPC=0, nHitsTRD=0, nHitsTOF=0, nHitsMUON=0 ;
265   for (Int_t iTrackRef=0; iTrackRef<part->GetNumberOfTrackReferences(); iTrackRef++) {
266     AliTrackReference * trackRef = part->GetTrackReference(iTrackRef);
267     if(trackRef){
268       Int_t detectorId = trackRef->DetectorId();
269       switch(detectorId) {
270       case AliTrackReference::kITS  : nHitsITS++  ; break ;
271       case AliTrackReference::kTPC  : nHitsTPC++  ; break ;
272       case AliTrackReference::kTRD  : nHitsTRD++  ; break ;
273       case AliTrackReference::kTOF  : nHitsTOF++  ; break ;
274       case AliTrackReference::kMUON : nHitsMUON++ ; break ;
275       default : break ;
276       }
277     }
278   }
279   
280   fhQA[kCutHitsITS ][afterCuts]->Fill(nHitsITS );
281   fhQA[kCutHitsTPC ][afterCuts]->Fill(nHitsTPC );
282   fhQA[kCutHitsTRD ][afterCuts]->Fill(nHitsTRD );
283   fhQA[kCutHitsTOF ][afterCuts]->Fill(nHitsTOF );
284   fhQA[kCutHitsMUON][afterCuts]->Fill(nHitsMUON);
285
286   // fill cut statistics and cut correlation histograms with information from the bitmap
287   if (afterCuts) return;
288
289   // Number of single cuts in this class
290   UInt_t ncuts = fBitmap->GetNbits();
291   for(UInt_t bit=0; bit<ncuts;bit++) {
292     if (!fBitmap->TestBitNumber(bit)) {
293       fhCutStatistics->Fill(bit+1);
294       for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
295         if (!fBitmap->TestBitNumber(bit2)) 
296           fhCutCorrelation->Fill(bit+1,bit2+1);
297       }
298     }
299   }
300 }