1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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
29 #include "AliESDEvent.h"
30 #include "AliESDVertex.h"
31 #include "AliCFEventRecCuts.h"
32 ClassImp(AliCFEventRecCuts)
33 //____________________________________________________________________
34 AliCFEventRecCuts::AliCFEventRecCuts() :
38 fRequireVtxCuts(kFALSE),
57 //____________________________________________________________________
58 AliCFEventRecCuts::AliCFEventRecCuts(Char_t* name, Char_t* title) :
59 AliCFCutBase(name,title),
62 fRequireVtxCuts(kFALSE),
81 //____________________________________________________________________
82 AliCFEventRecCuts::AliCFEventRecCuts(const AliCFEventRecCuts& c) :
84 fNTracksMin(c.fNTracksMin),
85 fNTracksMax(c.fNTracksMax),
86 fRequireVtxCuts(c.fRequireVtxCuts),
93 fVtxXResMax(c.fVtxXResMax),
94 fVtxYResMax(c.fVtxYResMax),
95 fVtxZResMax(c.fVtxZResMax),
101 for (Int_t i=0; i<c.kNCuts; i++){
102 for (Int_t j=0; j<c.kNStepQA; j++){
103 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
109 //____________________________________________________________________
110 AliCFEventRecCuts::~AliCFEventRecCuts() {
115 for (Int_t i=0; i<kNCuts; i++){
116 for (Int_t j=0; j<kNStepQA; j++){
117 if(fhQA[i][j]) delete fhQA[i][j];
121 if(fBitMap)delete fBitMap;
124 //_____________________________________________________________________________
125 void AliCFEventRecCuts::Initialise()
133 // sets pointers to histos to zero
136 for(Int_t i=0; i<kNCuts; i++){
137 for(Int_t j =0; j<kNStepQA; j++){
143 //____________________________________________________________________
144 AliCFEventRecCuts& AliCFEventRecCuts::operator=(const AliCFEventRecCuts& c)
147 // Assignment operator
150 AliCFCutBase::operator=(c) ;
151 fNTracksMin=c.fNTracksMin;
152 fNTracksMax=c.fNTracksMax;
153 fRequireVtxCuts=c.fRequireVtxCuts;
160 fVtxXResMax=c.fVtxXResMax;
161 fVtxYResMax=c.fVtxYResMax;
162 fVtxZResMax=c.fVtxZResMax;
166 for (Int_t i=0; i<c.kNCuts; i++){
167 for (Int_t j=0; j<c.kNStepQA; j++){
168 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
176 //____________________________________________________________________
177 Bool_t AliCFEventRecCuts::IsSelected(TObject* obj) {
179 //Check if the requested cuts are passed
183 SelectionBitMap(obj);
185 if (fIsQAOn) FillHistograms(obj,0);
186 Bool_t isSelected = kTRUE;
188 for (UInt_t icut=0; icut<fBitMap->GetNbits();icut++)
189 if(!fBitMap->TestBitNumber(icut)) isSelected = kFALSE;
191 if (!isSelected) return kFALSE ;
192 if (fIsQAOn) FillHistograms(obj,1);
196 //____________________________________________________________________
197 void AliCFEventRecCuts::SelectionBitMap(TObject* obj) {
199 //cut on the number of charged tracks and on the event vertex.
200 //so far specific to AliESDEvents
203 //Check if the requested cuts are passed and return a bitmap
204 for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE);
205 AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
208 //now start checking the cuts,
209 //first assume the event will be accepted:
210 for(Int_t j=0;j<kNCuts;j++)fBitMap->SetBitNumber(j,kTRUE);
212 //Number of charged tracks:
213 Int_t nTracks = esd->GetNumberOfTracks();
214 if(nTracks<fNTracksMin || nTracks>fNTracksMax)
215 fBitMap->SetBitNumber(0,kFALSE);
218 const AliESDVertex* vtxESD = esd->GetVertex();
220 for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE);
223 // Require the vertex to have been reconstructed successfully
224 if (strcmp(vtxESD->GetName(), "default")==0){
225 AliWarning(Form(" No reconstructed vertex found, skip event"));
226 for(Int_t j=1;j<kNCuts;j++)fBitMap->SetBitNumber(j,kFALSE);
229 // Pick up the position and uncertainties
232 vtxPos[0] = vtxESD->GetXv();
233 vtxPos[1] = vtxESD->GetYv();
234 vtxPos[2] = vtxESD->GetZv();
237 vtxRes[0] = vtxESD->GetXRes();
238 vtxRes[1] = vtxESD->GetYRes();
239 vtxRes[2] = vtxESD->GetZRes();
243 if (vtxPos[0]>fVtxXMax || vtxPos[0]<fVtxXMin)
244 fBitMap->SetBitNumber(1,kFALSE);
245 if (vtxPos[1]>fVtxYMax || vtxPos[1]<fVtxYMin)
246 fBitMap->SetBitNumber(2,kFALSE);
247 if (vtxPos[2]>fVtxZMax || vtxPos[2]<fVtxZMin)
248 fBitMap->SetBitNumber(3,kFALSE);
249 if (vtxRes[0]==0 || vtxRes[0]>fVtxXResMax)
250 fBitMap->SetBitNumber(4,kFALSE);
251 if (vtxRes[1]==0 || vtxRes[1]>fVtxYResMax)
252 fBitMap->SetBitNumber(5,kFALSE);
253 if (vtxRes[2]==0 || vtxRes[2]>fVtxZResMax)
254 fBitMap->SetBitNumber(6,kFALSE);
259 //_____________________________________________________________________________
260 void AliCFEventRecCuts::FillHistograms(TObject* obj, Bool_t b)
263 // fill the QA histograms
267 // cast TObject into VParticle
268 AliESDEvent* esd = dynamic_cast<AliESDEvent *>(obj);
271 // index = 0: fill histograms before cuts
272 // index = 1: fill histograms after cuts
274 index = ((b) ? 1 : 0);
277 //number of charged tracks:
278 Int_t nTracks = esd->GetNumberOfTracks();
279 fhQA[kNTracks][index]->Fill(nTracks);
281 //look at vertex parameters:
282 const AliESDVertex* vtxESD = esd->GetVertex();
284 // Require the vertex to have been reconstructed successfully
285 if (strcmp(vtxESD->GetName(), "default")==0)return;
286 // vertex position and uncertainties
287 fhQA[kVtxPosX][index]->Fill(vtxESD->GetXv());
288 fhQA[kVtxPosY][index]->Fill(vtxESD->GetYv());
289 fhQA[kVtxPosZ][index]->Fill(vtxESD->GetZv());
290 fhQA[kVtxResX][index]->Fill(vtxESD->GetXRes());
291 fhQA[kVtxResY][index]->Fill(vtxESD->GetYRes());
292 fhQA[kVtxResZ][index]->Fill(vtxESD->GetZRes());
296 //____________________________________________________________________
297 void AliCFEventRecCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
300 //setting x-axis bin limits of QA histogram fhQA[index]
303 for(Int_t i=0;i<kNStepQA;i++){
304 if(!fhQA[index][i]){AliWarning("non-existing histogram!");
307 fhQA[index][i]->GetXaxis()->Set(nbins,bins);
310 //____________________________________________________________________
311 void AliCFEventRecCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
314 //setting x-axis bins and range of QA histogram fhQA[index]
317 for(Int_t i=0;i<kNStepQA;i++){
318 if(!fhQA[index][i]){AliWarning("non-existing histogram!");
321 fhQA[index][i]->GetXaxis()->Set(nbins,xmin,xmax);
325 //_____________________________________________________________________________
326 void AliCFEventRecCuts::DefineHistograms() {
328 // histograms for cut variables
333 AliInfo(Form("No QA histos requested, Please first set the QA flag on!"));
337 // book QA histograms
340 for (Int_t i=0; i<kNStepQA; i++) {
341 if (i==0) sprintf(str," ");
342 else sprintf(str,"_cut");
344 fhQA[kNTracks][i] = new TH1F(Form("%s_NTracks%s",GetName(),str), "",501,-0.5,500.5);
345 fhQA[kVtxPosX][i] = new TH1F(Form("%s_Vtx_Pos_X%s",GetName(),str), "",100,-5.,5.);
346 fhQA[kVtxPosY][i] = new TH1F(Form("%s_Vtx_Pos_Y%s",GetName(),str), "",100,-5.,5.);
347 fhQA[kVtxPosZ][i] = new TH1F(Form("%s_Vtx_Pos_Z%s",GetName(),str), "",200,-50.,50.);
349 fhQA[kVtxResX][i] = new TH1F(Form("%s_Vtx_Res_X%s",GetName(),str), "",100,-1.,1.);
350 fhQA[kVtxResY][i] = new TH1F(Form("%s_Vtx_Res_Y%s",GetName(),str), "",100,-1.,1.);
351 fhQA[kVtxResZ][i] = new TH1F(Form("%s_Vtx_Res_Z%s",GetName(),str), "",100,-1.,1.);
353 fhQA[kNTracks][i] ->SetXTitle("Number of ESD tracks");
354 fhQA[kVtxPosX][i] ->SetXTitle("Vertex Position X (cm)");
355 fhQA[kVtxPosY][i] ->SetXTitle("Vertex Position Y (cm)");
356 fhQA[kVtxPosZ][i] ->SetXTitle("Vertex Position Z (cm)");
357 fhQA[kVtxResX][i] ->SetXTitle("Vertex Resolution X (cm)");
358 fhQA[kVtxResY][i] ->SetXTitle("Vertex Resolution Y (cm)");
359 fhQA[kVtxResZ][i] ->SetXTitle("Vertex Resolution Z (cm)");
363 for(Int_t i=0; i<kNCuts; i++) fhQA[i][1]->SetLineColor(color);
367 //_____________________________________________________________________________
368 void AliCFEventRecCuts::AddQAHistograms(TList *qaList) {
370 // saves the histograms in a TList
375 for (Int_t j=0; j<kNStepQA; j++) {
376 for(Int_t i=0; i<kNCuts; i++)
377 qaList->Add(fhQA[i][j]);