]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliHFEextraEventCuts.cxx
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWG3 / 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 "AliHFEextraEventCuts.h"
32 ClassImp(AliHFEextraEventCuts) 
33 //____________________________________________________________________
34 AliHFEextraEventCuts::AliHFEextraEventCuts() : 
35   AliCFCutBase(),
36   fRequireVtxCuts(kFALSE),
37   fVtxZMax(1.e99),
38   fVtxZMin(-1.e99),
39   fVtxNCtrbMin(0),
40   fVtxMixed(0),
41   fBitMap(0x0)
42 {
43   //
44   //ctor
45   //
46 }
47
48 //____________________________________________________________________
49 AliHFEextraEventCuts::AliHFEextraEventCuts(Char_t* name, Char_t* title) : 
50   AliCFCutBase(name,title),
51   fRequireVtxCuts(kFALSE),
52   fVtxZMax(1.e99),
53   fVtxZMin(-1.e99),
54   fVtxNCtrbMin(0),
55   fVtxMixed(0),
56   fBitMap(0x0)
57  {
58   //
59   //ctor
60   //
61   fBitMap=new TBits(0);
62   Initialise();
63  }
64
65 //____________________________________________________________________
66 AliHFEextraEventCuts::AliHFEextraEventCuts(const AliHFEextraEventCuts& c) : 
67   AliCFCutBase(c),
68   fRequireVtxCuts(c.fRequireVtxCuts),
69   fVtxZMax(c.fVtxZMax),
70   fVtxZMin(c.fVtxZMin),
71   fVtxNCtrbMin(c.fVtxNCtrbMin),
72   fVtxMixed(c.fVtxMixed),
73   fBitMap(c.fBitMap)
74 {
75   //
76   //copy constructor
77   //
78   for (Int_t i=0; i<c.kNCuts; i++){
79     for (Int_t j=0; j<c.kNStepQA; j++){
80       if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
81     }
82   }
83
84 }
85
86 //____________________________________________________________________
87 AliHFEextraEventCuts::~AliHFEextraEventCuts() {
88   //
89   //dtor
90   //
91
92   for (Int_t i=0; i<kNCuts; i++){
93     for (Int_t j=0; j<kNStepQA; j++){
94       if(fhQA[i][j]) delete fhQA[i][j];
95     }
96   }
97
98   if(fBitMap)delete fBitMap;
99
100 }
101 //_____________________________________________________________________________
102 void AliHFEextraEventCuts::Initialise()
103 {
104
105   //
106   //initialization
107   //
108
109   //
110   // sets pointers to histos to zero
111   //
112
113   for(Int_t i=0; i<kNCuts; i++){
114     for(Int_t j =0; j<kNStepQA; j++){
115       fhQA[i][j]=0x0;
116     }
117   }
118 }
119
120 //____________________________________________________________________
121 AliHFEextraEventCuts& AliHFEextraEventCuts::operator=(const AliHFEextraEventCuts& c)
122 {
123   //
124   // Assignment operator
125   //
126   if (this != &c) {
127     AliCFCutBase::operator=(c) ;
128     fRequireVtxCuts=c.fRequireVtxCuts;
129     fVtxZMax=c.fVtxZMax;
130     fVtxZMin=c.fVtxZMin;
131     fVtxNCtrbMin=c.fVtxNCtrbMin;
132     fVtxMixed=c.fVtxMixed;
133     fBitMap=c.fBitMap;
134   }
135
136   for (Int_t i=0; i<c.kNCuts; i++){
137     for (Int_t j=0; j<c.kNStepQA; j++){
138       if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
139     }
140   }
141
142
143   return *this ;
144 }
145
146 //____________________________________________________________________
147 Bool_t AliHFEextraEventCuts::IsSelected(TObject* obj) {
148   //
149   //Check if the requested cuts are passed
150   //
151
152
153   SelectionBitMap(obj);
154
155   if (fIsQAOn) FillHistograms(obj,0);
156   Bool_t isSelected = kTRUE;
157
158   for (UInt_t icut=0; icut<fBitMap->GetNbits();icut++)
159         if(!fBitMap->TestBitNumber(icut)) isSelected = kFALSE;
160
161   if (!isSelected) return kFALSE ;
162   if (fIsQAOn) FillHistograms(obj,1);
163   return kTRUE;
164
165 }
166 //____________________________________________________________________
167 void AliHFEextraEventCuts::SelectionBitMap(TObject* obj) {
168   //
169   //cut on the number of charged tracks and on the event vertex.
170   //so far specific to AliESDEvents
171   //
172
173   //Check if the requested cuts are passed and return a bitmap
174   for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE);
175   AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
176   if ( !esd ) return;
177
178   //now start checking the cuts,
179   //first assume the event will be accepted: 
180   for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kTRUE);
181
182
183   
184   if(fRequireVtxCuts){
185     const AliESDVertex* vtxESD = 0x0;
186     if      (fVtxMixed) {
187       vtxESD = esd->GetPrimaryVertexTracks() ;
188       if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
189         vtxESD = esd->GetPrimaryVertexSPD() ;
190         if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
191           for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
192           AliWarning("Cannot get vertex, skipping event");
193           return;
194         }
195       }
196     }
197     else {   
198       vtxESD = esd->GetPrimaryVertexTracks() ;
199     }
200
201     if(!vtxESD){
202       for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
203       AliWarning("Cannot get vertex, skipping event");
204       return;
205     }
206     
207     // Pick up the position and uncertainties
208     vtxESD = esd->GetPrimaryVertex();
209     if(!vtxESD){
210       for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE); 
211       AliWarning("Cannot get vertex, skipping event");
212       return;
213     }
214     Double_t vtxPos[3];
215     vtxPos[0] = vtxESD->GetXv();
216     vtxPos[1] = vtxESD->GetYv();
217     vtxPos[2] = vtxESD->GetZv();
218     
219     Int_t nCtrb = vtxESD->GetNContributors();
220
221     // Apply the cut
222     
223     if (vtxPos[2]>fVtxZMax || vtxPos[2]<fVtxZMin)
224       fBitMap->SetBitNumber(0,kFALSE); 
225     if (nCtrb<fVtxNCtrbMin)
226       fBitMap->SetBitNumber(1,kFALSE);
227
228   }  
229   return;
230 }
231
232 //_____________________________________________________________________________
233 void AliHFEextraEventCuts::FillHistograms(TObject* obj, Bool_t b)
234 {
235   //
236   // fill the QA histograms
237   //
238
239   if(!fIsQAOn) return;
240   // cast TObject into VParticle
241   AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
242   if (!esd ) return  ;
243
244   // index = 0: fill histograms before cuts
245   // index = 1: fill histograms after cuts
246   Int_t index = -1;
247   index = ((b) ? 1 : 0);
248
249
250   //look at vertex parameters:
251   const AliESDVertex* vtxESD = 0x0;
252   if      (fVtxMixed) {
253     vtxESD = esd->GetPrimaryVertexTracks() ;
254     if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
255       vtxESD = esd->GetPrimaryVertexSPD() ;
256       if((!vtxESD) || (vtxESD->GetNContributors() <= 0)) {
257         return;
258       }
259     }
260   }
261   else {   
262     vtxESD = esd->GetPrimaryVertexTracks() ;
263   }
264   if(!vtxESD)return;
265   // vertex position and uncertainties
266   fhQA[kVtxPosZ] [index]->Fill(vtxESD->GetZv());
267   fhQA[kVtxNCtrb][index]->Fill(vtxESD->GetNContributors());
268   
269 }
270
271 //____________________________________________________________________
272 void AliHFEextraEventCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
273 {
274   //
275   //setting x-axis bin limits of QA histogram fhQA[index] 
276   // 
277
278   for(Int_t i=0;i<kNStepQA;i++){
279     if(!fhQA[index][i]){AliWarning("non-existing histogram!");
280     return;
281     }
282     fhQA[index][i]->GetXaxis()->Set(nbins,bins);
283   }
284 }
285 //____________________________________________________________________
286 void AliHFEextraEventCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
287 {
288   //
289   //setting x-axis bins and range of QA histogram fhQA[index] 
290   // 
291
292   for(Int_t i=0;i<kNStepQA;i++){
293     if(!fhQA[index][i]){AliWarning("non-existing histogram!");
294     return;
295     }
296     fhQA[index][i]->GetXaxis()->Set(nbins,xmin,xmax);
297   }
298 }
299
300 //_____________________________________________________________________________
301  void AliHFEextraEventCuts::DefineHistograms() {
302   //
303   // histograms for cut variables
304   //
305   Int_t color = 2;
306
307   if(!fIsQAOn) {
308     AliInfo(Form("No QA histos requested, Please first set the QA flag on!"));
309     return;
310   }  
311   
312   // book QA histograms
313
314   Char_t str[5];
315   for (Int_t i=0; i<kNStepQA; i++) {
316     if (i==0) snprintf(str,5," ");
317     else snprintf(str,5,"_cut");
318
319     fhQA[kVtxPosZ][i]   = new  TH1F(Form("%s_Vtx_Pos_Z%s",GetName(),str),               "",200,-50.,50.);
320     fhQA[kVtxNCtrb][i]  = new  TH1F(Form("%s_Vtx_N_Ctrb%s",GetName(),str),              "",1000,0.,1000.);
321  
322     fhQA[kVtxPosZ][i]   ->SetXTitle("Vertex Position Z (cm)");
323     fhQA[kVtxNCtrb][i]  ->SetXTitle("Number of contributors");
324   }
325
326   for(Int_t i=0; i<kNCuts; i++) fhQA[i][1]->SetLineColor(color);
327
328 }
329
330 //_____________________________________________________________________________
331 void AliHFEextraEventCuts::AddQAHistograms(TList *qaList) {
332   //
333   // saves the histograms in a TList
334   //
335
336   DefineHistograms();
337
338   for (Int_t j=0; j<kNStepQA; j++) {
339     for(Int_t i=0; i<kNCuts; i++)
340         qaList->Add(fhQA[i][j]);
341   }
342 }