]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/AliHFEextraEventCuts.cxx
Make a Makechain function static
[u/mrichter/AliRoot.git] / PWGHF / hfe / AliHFEextraEventCuts.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 // Cut on the Event at reconstructed level: for the moment 
16 // the requirements on the number of charged tracks and on 
17 // the vertex position and resolution are implemented
18 // The argument of IsSelected member function (passed object) is cast into 
19 // an AliESDEvent. In the future may be modified to use AliVEvent interface
20 // and include more cut variables.
21 // The class derives from AliCFCutBase
22 // Author:S.Arcelli Silvia.Arcelli@cern.ch
23 //
24 //
25 #include "TH1F.h"
26 #include "TList.h"
27 #include "TBits.h"
28 #include "AliLog.h"
29 #include "AliESDEvent.h"
30 #include "AliESDVertex.h"
31 #include "AliAODEvent.h"
32 #include "AliAODVertex.h"
33 #include "AliHFEextraEventCuts.h"
34 ClassImp(AliHFEextraEventCuts) 
35 //____________________________________________________________________
36 AliHFEextraEventCuts::AliHFEextraEventCuts() : 
37   AliCFCutBase(),
38   fRequireVtxCuts(kFALSE),
39   fVtxZMax(1.e99),
40   fVtxZMin(-1.e99),
41   fVtxNCtrbMin(0),
42   fVtxMixed(0),
43   fBitMap(0x0)
44 {
45   //
46   //ctor
47   //
48 }
49
50 //____________________________________________________________________
51 AliHFEextraEventCuts::AliHFEextraEventCuts(Char_t* name, Char_t* title) : 
52   AliCFCutBase(name,title),
53   fRequireVtxCuts(kFALSE),
54   fVtxZMax(1.e99),
55   fVtxZMin(-1.e99),
56   fVtxNCtrbMin(0),
57   fVtxMixed(0),
58   fBitMap(0x0)
59  {
60   //
61   //ctor
62   //
63   fBitMap=new TBits(0);
64   Initialise();
65  }
66
67 //____________________________________________________________________
68 AliHFEextraEventCuts::AliHFEextraEventCuts(const AliHFEextraEventCuts& c) : 
69   AliCFCutBase(c),
70   fRequireVtxCuts(c.fRequireVtxCuts),
71   fVtxZMax(c.fVtxZMax),
72   fVtxZMin(c.fVtxZMin),
73   fVtxNCtrbMin(c.fVtxNCtrbMin),
74   fVtxMixed(c.fVtxMixed),
75   fBitMap(c.fBitMap)
76 {
77   //
78   //copy constructor
79   //
80   for (Int_t i=0; i<c.kNCuts; i++){
81     for (Int_t j=0; j<c.kNStepQA; j++){
82       if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
83     }
84   }
85
86 }
87
88 //____________________________________________________________________
89 AliHFEextraEventCuts::~AliHFEextraEventCuts() {
90   //
91   //dtor
92   //
93
94   for (Int_t i=0; i<kNCuts; i++){
95     for (Int_t j=0; j<kNStepQA; j++){
96       if(fhQA[i][j]) delete fhQA[i][j];
97     }
98   }
99
100   if(fBitMap)delete fBitMap;
101
102 }
103 //_____________________________________________________________________________
104 void AliHFEextraEventCuts::Initialise()
105 {
106
107   //
108   //initialization
109   //
110
111   //
112   // sets pointers to histos to zero
113   //
114
115   for(Int_t i=0; i<kNCuts; i++){
116     for(Int_t j =0; j<kNStepQA; j++){
117       fhQA[i][j]=0x0;
118     }
119   }
120 }
121
122 //____________________________________________________________________
123 AliHFEextraEventCuts& AliHFEextraEventCuts::operator=(const AliHFEextraEventCuts& c)
124 {
125   //
126   // Assignment operator
127   //
128   if (this != &c) {
129     AliCFCutBase::operator=(c) ;
130     fRequireVtxCuts=c.fRequireVtxCuts;
131     fVtxZMax=c.fVtxZMax;
132     fVtxZMin=c.fVtxZMin;
133     fVtxNCtrbMin=c.fVtxNCtrbMin;
134     fVtxMixed=c.fVtxMixed;
135     fBitMap=c.fBitMap;
136   }
137
138   for (Int_t i=0; i<c.kNCuts; i++){
139     for (Int_t j=0; j<c.kNStepQA; j++){
140       if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
141     }
142   }
143
144
145   return *this ;
146 }
147
148 //____________________________________________________________________
149 Bool_t AliHFEextraEventCuts::IsSelected(TObject* obj) {
150   //
151   //Check if the requested cuts are passed
152   //
153
154
155   SelectionBitMap(obj);
156
157   if (fIsQAOn) FillHistograms(obj,0);
158   Bool_t isSelected = kTRUE;
159
160   for (UInt_t icut=0; icut<fBitMap->GetNbits();icut++)
161         if(!fBitMap->TestBitNumber(icut)) isSelected = kFALSE;
162
163   if (!isSelected) return kFALSE ;
164   if (fIsQAOn) FillHistograms(obj,1);
165   return kTRUE;
166
167 }
168 //____________________________________________________________________
169 void AliHFEextraEventCuts::SelectionBitMap(TObject* obj) {
170   //
171   //cut on the number of charged tracks and on the event vertex.
172   //so far specific to AliESDEvents
173   //
174
175   //Check if the requested cuts are passed and return a bitmap
176   for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE);
177
178
179
180   AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
181   AliAODEvent* aod = dynamic_cast<AliAODEvent *>(obj);
182
183
184   // ESD
185
186   if ( esd ) {
187     
188     //now start checking the cuts,
189     //first assume the event will be accepted: 
190     for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kTRUE);
191     
192     
193     
194     if(fRequireVtxCuts){
195       const AliESDVertex* vtxESD = 0x0;
196       if      (fVtxMixed) {
197         vtxESD = esd->GetPrimaryVertexTracks() ;
198         if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
199           vtxESD = esd->GetPrimaryVertexSPD() ;
200           if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
201             for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
202             AliWarning("Cannot get vertex, skipping event");
203             return;
204           }
205         }
206       }
207       else {   
208         vtxESD = esd->GetPrimaryVertexTracks() ;
209       }
210       
211       if(!vtxESD){
212         for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
213         AliWarning("Cannot get vertex, skipping event");
214         return;
215       }
216       
217       // Pick up the position and uncertainties
218       vtxESD = esd->GetPrimaryVertex();
219       if(!vtxESD){
220         for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
221         AliWarning("Cannot get vertex, skipping event");
222         return;
223       }
224       Double_t vtxPos[3];
225       vtxPos[0] = vtxESD->GetXv();
226       vtxPos[1] = vtxESD->GetYv();
227       vtxPos[2] = vtxESD->GetZv();
228       
229       Int_t nCtrb = vtxESD->GetNContributors();
230       
231       // Apply the cut
232       
233       if (vtxPos[2]>fVtxZMax || vtxPos[2]<fVtxZMin)
234         fBitMap->SetBitNumber(0,kFALSE); 
235       if (nCtrb<fVtxNCtrbMin)
236         fBitMap->SetBitNumber(1,kFALSE);
237       
238     }  
239     return;
240     
241     
242   }
243
244   //
245
246   // AOD
247
248   if ( aod ) {
249     
250     //now start checking the cuts,
251     //first assume the event will be accepted: 
252     for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kTRUE);
253     
254     
255     
256     if(fRequireVtxCuts){
257       const AliAODVertex* vtxAOD = 0x0;
258       if      (fVtxMixed) {
259         vtxAOD = aod->GetPrimaryVertex();
260         if((!vtxAOD) || (vtxAOD->GetNContributors() <= 0)) {
261           vtxAOD = aod->GetPrimaryVertexSPD() ;
262           if((!vtxAOD) || (vtxAOD->GetNContributors() <= 0)) {
263             for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
264             AliWarning("Cannot get vertex, skipping event");
265             return;
266           }
267         }
268       }
269       else {   
270         vtxAOD = aod->GetPrimaryVertex() ;
271       }
272       
273       if(!vtxAOD){
274         for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
275         AliWarning("Cannot get vertex, skipping event");
276         return;
277       }
278       
279       // Pick up the position and uncertainties
280       vtxAOD = aod->GetPrimaryVertex();
281       if(!vtxAOD){
282         for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
283         AliWarning("Cannot get vertex, skipping event");
284         return;
285       }
286       Double_t vtxPos[3];
287       vtxPos[0] = vtxAOD->GetX();
288       vtxPos[1] = vtxAOD->GetY();
289       vtxPos[2] = vtxAOD->GetZ();
290       
291       Int_t nCtrb = vtxAOD->GetNContributors();
292       //printf("AliHFEextraCuts:: %d, %f, %f, %f\n",nCtrb,vtxPos[0],vtxPos[1],vtxPos[2]);
293       
294       // Apply the cut
295       
296       if (vtxPos[2]>fVtxZMax || vtxPos[2]<fVtxZMin)
297         fBitMap->SetBitNumber(0,kFALSE); 
298       if (nCtrb<fVtxNCtrbMin)
299         fBitMap->SetBitNumber(1,kFALSE);
300       
301     }  
302     return;
303     
304     
305   }
306   
307   return;
308   
309 }
310
311 //_____________________________________________________________________________
312 void AliHFEextraEventCuts::FillHistograms(TObject* obj, Bool_t b)
313 {
314   //
315   // fill the QA histograms
316   //
317
318   if(!fIsQAOn) return;
319   // cast TObject into VParticle
320   AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
321   if (!esd ) return  ;
322
323   // index = 0: fill histograms before cuts
324   // index = 1: fill histograms after cuts
325   Int_t index = -1;
326   index = ((b) ? 1 : 0);
327
328
329   //look at vertex parameters:
330   const AliESDVertex* vtxESD = 0x0;
331   if      (fVtxMixed) {
332     vtxESD = esd->GetPrimaryVertexTracks() ;
333     if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
334       vtxESD = esd->GetPrimaryVertexSPD() ;
335       if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
336         return;
337       }
338     }
339   }
340   else {   
341     vtxESD = esd->GetPrimaryVertexTracks() ;
342   }
343   if(!vtxESD)return;
344   // vertex position and uncertainties
345   fhQA[kVtxPosZ] [index]->Fill(vtxESD->GetZv());
346   fhQA[kVtxNCtrb][index]->Fill(vtxESD->GetNContributors());
347   
348 }
349
350 //____________________________________________________________________
351 void AliHFEextraEventCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
352 {
353   //
354   //setting x-axis bin limits of QA histogram fhQA[index] 
355   // 
356
357   for(Int_t i=0;i<kNStepQA;i++){
358     if(!fhQA[index][i]){AliWarning("non-existing histogram!");
359     return;
360     }
361     fhQA[index][i]->GetXaxis()->Set(nbins,bins);
362   }
363 }
364 //____________________________________________________________________
365 void AliHFEextraEventCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
366 {
367   //
368   //setting x-axis bins and range of QA histogram fhQA[index] 
369   // 
370
371   for(Int_t i=0;i<kNStepQA;i++){
372     if(!fhQA[index][i]){AliWarning("non-existing histogram!");
373     return;
374     }
375     fhQA[index][i]->GetXaxis()->Set(nbins,xmin,xmax);
376   }
377 }
378
379 //_____________________________________________________________________________
380  void AliHFEextraEventCuts::DefineHistograms() {
381   //
382   // histograms for cut variables
383   //
384   Int_t color = 2;
385
386   if(!fIsQAOn) {
387     AliInfo(Form("No QA histos requested, Please first set the QA flag on!"));
388     return;
389   }  
390   
391   // book QA histograms
392
393   Char_t str[5];
394   for (Int_t i=0; i<kNStepQA; i++) {
395     if (i==0) snprintf(str,5," ");
396     else snprintf(str,5,"_cut");
397
398     fhQA[kVtxPosZ][i]   = new  TH1F(Form("%s_Vtx_Pos_Z%s",GetName(),str),               "",200,-50.,50.);
399     fhQA[kVtxNCtrb][i]  = new  TH1F(Form("%s_Vtx_N_Ctrb%s",GetName(),str),              "",1000,0.,1000.);
400  
401     fhQA[kVtxPosZ][i]   ->SetXTitle("Vertex Position Z (cm)");
402     fhQA[kVtxNCtrb][i]  ->SetXTitle("Number of contributors");
403   }
404
405   for(Int_t i=0; i<kNCuts; i++) fhQA[i][1]->SetLineColor(color);
406
407 }
408
409 //_____________________________________________________________________________
410 void AliHFEextraEventCuts::AddQAHistograms(TList *qaList) {
411   //
412   // saves the histograms in a TList
413   //
414
415   DefineHistograms();
416
417   for (Int_t j=0; j<kNStepQA; j++) {
418     for(Int_t i=0; i<kNCuts; i++)
419         qaList->Add(fhQA[i][j]);
420   }
421 }