]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliFlowEventCuts.cxx
from Carlos Eugenio Perez Lara:
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliFlowEventCuts.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 /* $Id$ */ 
17
18 // AliFlowEventCuts:
19 // An event cut class for the flow framework
20 //
21 // origin: Mikolaj Krzewicki (mikolaj.krzewicki@cern.ch)
22
23 #include <limits.h>
24 #include <float.h>
25 #include <TList.h>
26 #include <TH1F.h>
27 #include <TH2F.h>
28 #include <TBrowser.h>
29 #include "TMath.h"
30 #include "TNamed.h"
31 #include "AliVVertex.h"
32 #include "AliVEvent.h"
33 #include "AliESDEvent.h"
34 #include "AliAODEvent.h"
35 #include "AliAODHeader.h"
36 #include "AliCentrality.h"
37 #include "AliESDVZERO.h"
38 #include "AliMultiplicity.h"
39 #include "AliMCEvent.h"
40 #include "AliFlowEventCuts.h"
41 #include "AliFlowTrackCuts.h"
42 #include "AliTriggerAnalysis.h"
43
44 ClassImp(AliFlowEventCuts)
45
46 //-----------------------------------------------------------------------
47 AliFlowEventCuts::AliFlowEventCuts():
48   TNamed(),
49   fQA(NULL),
50   fCutNumberOfTracks(kFALSE),
51   fNumberOfTracksMax(INT_MAX),
52   fNumberOfTracksMin(INT_MIN),
53   fCutRefMult(kFALSE),
54   fRefMultMethod(kTPConly),
55   fUseAliESDtrackCutsRefMult(kFALSE),
56   fRefMultMethodAliESDtrackCuts(AliESDtrackCuts::kTrackletsITSTPC),
57   fRefMultMax(INT_MAX),
58   fRefMultMin(INT_MIN),
59   fRefMultCuts(NULL),
60   fMeanPtCuts(NULL),
61   fStandardTPCcuts(NULL),
62   fStandardGlobalCuts(NULL),
63   fCutPrimaryVertexX(kFALSE),
64   fPrimaryVertexXmax(INT_MAX),
65   fPrimaryVertexXmin(INT_MIN),
66   fCutPrimaryVertexY(kFALSE),
67   fPrimaryVertexYmax(INT_MAX),
68   fPrimaryVertexYmin(INT_MIN),
69   fCutPrimaryVertexZ(kFALSE),
70   fPrimaryVertexZmax(INT_MAX),
71   fPrimaryVertexZmin(INT_MIN),
72   fCutNContributors(kFALSE),
73   fNContributorsMax(INT_MAX),
74   fNContributorsMin(INT_MIN),
75   fCutMeanPt(kFALSE),
76   fMeanPtMax(-DBL_MAX),
77   fMeanPtMin(DBL_MAX),
78   fCutSPDvertexerAnomaly(kFALSE),
79   fCutSPDTRKVtxZ(kFALSE),
80   fCutTPCmultiplicityOutliers(kFALSE),
81   fCutCentralityPercentile(kFALSE),
82   fUseCentralityUnchecked(kFALSE),
83   fCentralityPercentileMethod(kTPConly),
84   fCentralityPercentileMax(100.),
85   fCentralityPercentileMin(0.),
86   fCutZDCtiming(kFALSE),
87   fTrigAna()
88 {
89   //constructor 
90 }
91
92 //-----------------------------------------------------------------------
93 AliFlowEventCuts::AliFlowEventCuts(const char* name, const char* title):
94   TNamed(name, title),
95   fQA(NULL),
96   fCutNumberOfTracks(kFALSE),
97   fNumberOfTracksMax(INT_MAX),
98   fNumberOfTracksMin(INT_MIN),
99   fCutRefMult(kFALSE),
100   fRefMultMethod(kTPConly),
101   fUseAliESDtrackCutsRefMult(kFALSE),
102   fRefMultMethodAliESDtrackCuts(AliESDtrackCuts::kTrackletsITSTPC),
103   fRefMultMax(INT_MAX),
104   fRefMultMin(INT_MIN),
105   fRefMultCuts(NULL),
106   fMeanPtCuts(NULL),
107   fStandardTPCcuts(AliFlowTrackCuts::GetStandardTPCStandaloneTrackCuts2010()),
108   fStandardGlobalCuts(AliFlowTrackCuts::GetStandardGlobalTrackCuts2010()),
109   fCutPrimaryVertexX(kFALSE),
110   fPrimaryVertexXmax(INT_MAX),
111   fPrimaryVertexXmin(INT_MIN),
112   fCutPrimaryVertexY(kFALSE),
113   fPrimaryVertexYmax(INT_MAX),
114   fPrimaryVertexYmin(INT_MIN),
115   fCutPrimaryVertexZ(kFALSE),
116   fPrimaryVertexZmax(INT_MAX),
117   fPrimaryVertexZmin(INT_MIN),
118   fCutNContributors(kFALSE),
119   fNContributorsMax(INT_MAX),
120   fNContributorsMin(INT_MIN),
121   fCutMeanPt(kFALSE),
122   fMeanPtMax(-DBL_MAX),
123   fMeanPtMin(DBL_MAX),
124   fCutSPDvertexerAnomaly(kFALSE),
125   fCutSPDTRKVtxZ(kFALSE),
126   fCutTPCmultiplicityOutliers(kFALSE),
127   fCutCentralityPercentile(kFALSE),
128   fUseCentralityUnchecked(kFALSE),
129   fCentralityPercentileMethod(kTPConly),
130   fCentralityPercentileMax(100.),
131   fCentralityPercentileMin(0.),
132   fCutZDCtiming(kFALSE),
133   fTrigAna()
134 {
135   //constructor 
136 }
137
138 ////-----------------------------------------------------------------------
139 AliFlowEventCuts::AliFlowEventCuts(const AliFlowEventCuts& that):
140   TNamed(that),
141   fQA(NULL),
142   fCutNumberOfTracks(that.fCutNumberOfTracks),
143   fNumberOfTracksMax(that.fNumberOfTracksMax),
144   fNumberOfTracksMin(that.fNumberOfTracksMin),
145   fCutRefMult(that.fCutRefMult),
146   fRefMultMethod(that.fRefMultMethod),
147   fUseAliESDtrackCutsRefMult(that.fUseAliESDtrackCutsRefMult),
148   fRefMultMethodAliESDtrackCuts(that.fRefMultMethodAliESDtrackCuts),
149   fRefMultMax(that.fRefMultMax),
150   fRefMultMin(that.fRefMultMin),
151   fRefMultCuts(NULL),
152   fMeanPtCuts(NULL),
153   fStandardTPCcuts(NULL),
154   fStandardGlobalCuts(NULL),
155   fCutPrimaryVertexX(that.fCutPrimaryVertexX),
156   fPrimaryVertexXmax(that.fPrimaryVertexXmax),
157   fPrimaryVertexXmin(that.fPrimaryVertexXmin),
158   fCutPrimaryVertexY(that.fCutPrimaryVertexX),
159   fPrimaryVertexYmax(that.fPrimaryVertexYmax),
160   fPrimaryVertexYmin(that.fPrimaryVertexYmin),
161   fCutPrimaryVertexZ(that.fCutPrimaryVertexX),
162   fPrimaryVertexZmax(that.fPrimaryVertexZmax),
163   fPrimaryVertexZmin(that.fPrimaryVertexZmin),
164   fCutNContributors(that.fCutNContributors),
165   fNContributorsMax(that.fNContributorsMax),
166   fNContributorsMin(that.fNContributorsMin),
167   fCutMeanPt(that.fCutMeanPt),
168   fMeanPtMax(that.fMeanPtMax),
169   fMeanPtMin(that.fMeanPtMin),
170   fCutSPDvertexerAnomaly(that.fCutSPDvertexerAnomaly),
171   fCutSPDTRKVtxZ(that.fCutSPDTRKVtxZ),
172   fCutTPCmultiplicityOutliers(that.fCutTPCmultiplicityOutliers),
173   fCutCentralityPercentile(that.fCutCentralityPercentile),
174   fUseCentralityUnchecked(that.fUseCentralityUnchecked),
175   fCentralityPercentileMethod(that.fCentralityPercentileMethod),
176   fCentralityPercentileMax(that.fCentralityPercentileMax),
177   fCentralityPercentileMin(that.fCentralityPercentileMin),
178   fCutZDCtiming(that.fCutZDCtiming),
179   fTrigAna()
180 {
181   if (that.fQA) DefineHistograms();
182   //copy constructor 
183   if (that.fRefMultCuts)
184     fRefMultCuts = new AliFlowTrackCuts(*(that.fRefMultCuts));
185   if (that.fMeanPtCuts)
186     fMeanPtCuts = new AliFlowTrackCuts(*(that.fMeanPtCuts));
187   fStandardTPCcuts = AliFlowTrackCuts::GetStandardTPCStandaloneTrackCuts2010();
188   fStandardGlobalCuts = AliFlowTrackCuts::GetStandardGlobalTrackCuts2010();
189 }
190
191 ////-----------------------------------------------------------------------
192 AliFlowEventCuts::~AliFlowEventCuts()
193 {
194   //dtor
195   delete fMeanPtCuts;
196   delete fRefMultCuts;
197   delete fStandardGlobalCuts;
198   delete fStandardTPCcuts;
199   if (fQA) { fQA->SetOwner(); fQA->Delete(); delete fQA; }
200 }
201
202 ////-----------------------------------------------------------------------
203 AliFlowEventCuts& AliFlowEventCuts::operator=(const AliFlowEventCuts& that)
204 {
205   //assignment
206   if (this==&that) return *this;
207
208   if (that.fQA)
209   {
210     if (fQA)
211     {
212       fQA->Delete();
213       delete fQA;
214     }
215     fQA = static_cast<TList*>(that.fQA->Clone());
216   }
217   else
218   {
219     fQA->Delete();
220     delete fQA;
221     fQA=NULL;
222   }
223
224   fCutNumberOfTracks=that.fCutNumberOfTracks;
225   fNumberOfTracksMax=that.fNumberOfTracksMax;
226   fNumberOfTracksMin=that.fNumberOfTracksMin;
227   fCutRefMult=that.fCutRefMult;
228   fRefMultMethod=that.fRefMultMethod;
229   fUseAliESDtrackCutsRefMult=that.fUseAliESDtrackCutsRefMult;
230   fRefMultMethodAliESDtrackCuts=that.fRefMultMethodAliESDtrackCuts;
231   fRefMultMax=that.fRefMultMax;
232   fRefMultMin=that.fRefMultMin;
233   if (that.fRefMultCuts) *fRefMultCuts=*(that.fRefMultCuts);
234   if (that.fMeanPtCuts) *fMeanPtCuts=*(that.fMeanPtCuts);
235   fStandardTPCcuts = AliFlowTrackCuts::GetStandardTPCStandaloneTrackCuts2010();
236   fStandardGlobalCuts = AliFlowTrackCuts::GetStandardGlobalTrackCuts2010();
237   fCutPrimaryVertexX=that.fCutPrimaryVertexX;
238   fPrimaryVertexXmax=that.fPrimaryVertexXmax;
239   fPrimaryVertexXmin=that.fPrimaryVertexXmin;
240   fCutPrimaryVertexY=that.fCutPrimaryVertexY;
241   fPrimaryVertexYmax=that.fPrimaryVertexYmax;
242   fPrimaryVertexYmin=that.fPrimaryVertexYmin;
243   fCutPrimaryVertexZ=that.fCutPrimaryVertexZ;
244   fPrimaryVertexZmax=that.fPrimaryVertexZmax;
245   fPrimaryVertexZmin=that.fPrimaryVertexZmin;
246   fCutNContributors=that.fCutNContributors;
247   fNContributorsMax=that.fNContributorsMax;
248   fNContributorsMin=that.fNContributorsMin;
249   fCutMeanPt=that.fCutMeanPt;
250   fMeanPtMax=that.fMeanPtMax;
251   fMeanPtMin=that.fMeanPtMin;
252   fCutSPDvertexerAnomaly=that.fCutSPDvertexerAnomaly;
253   fCutSPDTRKVtxZ=that.fCutSPDTRKVtxZ;
254   fCutTPCmultiplicityOutliers=that.fCutTPCmultiplicityOutliers;
255   fCutCentralityPercentile=that.fCutCentralityPercentile;
256   fUseCentralityUnchecked=that.fUseCentralityUnchecked;
257   fCentralityPercentileMethod=that.fCentralityPercentileMethod;
258   fCentralityPercentileMax=that.fCentralityPercentileMax;
259   fCentralityPercentileMin=that.fCentralityPercentileMin;
260   fCutZDCtiming=that.fCutZDCtiming;
261   return *this;
262 }
263
264 //----------------------------------------------------------------------- 
265 Bool_t AliFlowEventCuts::IsSelected(TObject* obj)
266 {
267   //check cuts
268   AliVEvent* vevent = dynamic_cast<AliVEvent*>(obj);
269   if (vevent) return PassesCuts(vevent);
270   return kFALSE;  //when passed wrong type of object
271 }
272 //----------------------------------------------------------------------- 
273 Bool_t AliFlowEventCuts::PassesCuts(AliVEvent *event)
274 {
275   ///check if event passes cuts
276   const AliVVertex* pvtx=event->GetPrimaryVertex();
277   Double_t pvtxx = pvtx->GetX();
278   Double_t pvtxy = pvtx->GetY();
279   Double_t pvtxz = pvtx->GetZ();
280   Int_t ncontrib = pvtx->GetNContributors();
281   Bool_t pass=kTRUE;
282   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*>(event);
283   AliAODEvent* aodevent = dynamic_cast<AliAODEvent*>(event);
284   Int_t multTPC = 0;
285   Int_t multGlobal = 0; 
286   if (fQA)
287   {
288     multTPC = fStandardTPCcuts->Count(event);
289     multGlobal = fStandardGlobalCuts->Count(event);
290     QAbefore(0)->Fill(pvtxz);
291     QAbefore(1)->Fill(multGlobal,multTPC);
292   }
293   if (fCutTPCmultiplicityOutliers && esdevent)
294   {
295     //this is pretty slow as we check the event track by track twice
296     //this cut will work for 2010 PbPb data and is dependent on
297     //TPC and ITS reco efficiency (e.g. geometry, calibration etc)
298     if (!fQA)
299     {
300       multTPC = fStandardTPCcuts->Count(event);
301       multGlobal = fStandardGlobalCuts->Count(event);
302     }
303     if (multTPC > ( 23+1.216*multGlobal)) {pass=kFALSE;}
304     if (multTPC < (-20+1.087*multGlobal)) {pass=kFALSE;}
305   }
306   if (fCutNContributors)
307   {
308     if (ncontrib < fNContributorsMin || ncontrib >= fNContributorsMax) pass=kFALSE;
309   }
310   if (fCutPrimaryVertexX)
311   {
312     if (pvtxx < fPrimaryVertexXmin || pvtxx >= fPrimaryVertexXmax) pass=kFALSE;
313   }
314   if (fCutPrimaryVertexY)
315   {
316     if (pvtxy < fPrimaryVertexYmin || pvtxy >= fPrimaryVertexYmax) pass=kFALSE;
317   }
318   if (fCutPrimaryVertexZ)
319   {
320     if (pvtxz < fPrimaryVertexZmin || pvtxz >= fPrimaryVertexZmax)
321       pass=kFALSE;
322   }
323   if (fCutCentralityPercentile&&esdevent)
324   {
325     AliCentrality* centr = esdevent->GetCentrality();
326     if (fUseCentralityUnchecked)
327     {
328       if (!centr->IsEventInCentralityClassUnchecked( fCentralityPercentileMin,
329                                                      fCentralityPercentileMax,
330                                                      CentrMethName(fCentralityPercentileMethod) ))
331       {
332         pass=kFALSE;
333       }
334     }
335     else
336     {
337       if (!centr->IsEventInCentralityClass( fCentralityPercentileMin,
338                                             fCentralityPercentileMax,
339                                             CentrMethName(fCentralityPercentileMethod) ))
340       {
341         pass=kFALSE;
342       }
343     }
344   }
345   if (fCutSPDvertexerAnomaly&&esdevent)
346   {
347     const AliESDVertex* sdpvertex = esdevent->GetPrimaryVertexSPD();
348     if (sdpvertex->GetNContributors()<1) pass=kFALSE;
349     if (sdpvertex->GetDispersion()>0.04) pass=kFALSE;
350     if (sdpvertex->GetZRes()>0.25) pass=kFALSE;
351     const AliESDVertex* tpcvertex = esdevent->GetPrimaryVertexTPC();
352     if (tpcvertex->GetNContributors()<1) pass=kFALSE;
353     const AliMultiplicity* tracklets = esdevent->GetMultiplicity();
354     if (tpcvertex->GetNContributors()<(-10.0+0.25*tracklets->GetNumberOfITSClusters(0)))
355     {
356       pass=kFALSE;
357     }
358   }
359   if (fCutZDCtiming&&esdevent)
360   {
361     if (!fTrigAna.ZDCTimeTrigger(esdevent))
362     {
363       pass=kFALSE;
364     }
365   }
366   if(fCutNumberOfTracks) {if ( event->GetNumberOfTracks() < fNumberOfTracksMin ||
367                                event->GetNumberOfTracks() >= fNumberOfTracksMax ) pass=kFALSE;}
368   if(fCutRefMult&&esdevent)
369   {
370     //reference multiplicity still to be defined
371     Double_t refMult = RefMult(event);
372     if (refMult < fRefMultMin || refMult >= fRefMultMax )
373     {
374       pass=kFALSE;
375     }
376   }
377
378   // Handles AOD event
379   if(aodevent) {
380     if(fCutSPDTRKVtxZ) {
381       Double_t tVtxZ = aodevent->GetPrimaryVertex()->GetZ();
382       Double_t tSPDVtxZ = aodevent->GetPrimaryVertexSPD()->GetZ();
383       if( TMath::Abs(tVtxZ-tSPDVtxZ) > 0.5 ) pass = kFALSE;
384     }
385     AliCentrality* centr = aodevent->GetHeader()->GetCentralityP();
386     if(fCutTPCmultiplicityOutliers){
387       Double_t v0Centr  = centr->GetCentralityPercentile("V0M");
388       Double_t trkCentr = centr->GetCentralityPercentile("TRK"); 
389       if(TMath::Abs(v0Centr-trkCentr) > 5) pass = kFALSE;
390     }
391     if (fCutCentralityPercentile) {
392       if (fUseCentralityUnchecked) {
393         if (!centr->IsEventInCentralityClassUnchecked( fCentralityPercentileMin,
394                                                        fCentralityPercentileMax,
395                                                        CentrMethName(fCentralityPercentileMethod) )) {
396           pass = kFALSE;
397         }
398       } else {
399         if (!centr->IsEventInCentralityClass( fCentralityPercentileMin,
400                                               fCentralityPercentileMax,
401                                               CentrMethName(fCentralityPercentileMethod) )) {
402           pass = kFALSE;
403         }
404       }
405     }
406   }
407
408   if (fCutMeanPt)
409   {
410     Float_t meanpt=0.0;
411     Int_t ntracks=event->GetNumberOfTracks();
412     Int_t nselected=0;
413     for (Int_t i=0; i<ntracks; i++)
414     {
415       AliVParticle* track = event->GetTrack(i);
416       if (!track) continue;
417       Bool_t localpass=kTRUE;
418       if (fMeanPtCuts) localpass=fMeanPtCuts->IsSelected(track);
419       if (localpass) 
420       {
421         meanpt += track->Pt();
422         nselected++;
423       }
424     }
425     if (nselected) meanpt=meanpt/nselected;
426     if (meanpt<fMeanPtMin || meanpt >= fMeanPtMax) pass=kFALSE;
427   }
428   if (fQA&&pass) 
429   {
430     QAafter(1)->Fill(multGlobal,multTPC);
431     QAafter(0)->Fill(pvtxz);
432   }
433   return pass;
434 }
435
436 //----------------------------------------------------------------------- 
437 const char* AliFlowEventCuts::CentrMethName(refMultMethod method) const
438 {
439   //get the string for refmultmethod, for use with AliCentrality in
440   //the cut on centrality percentile
441   switch (method)
442   {
443     case kSPDtracklets:
444       return "TKL";
445     case kSPD1clusters:
446       return "CL1";
447     case kTPConly:
448       return "TRK";
449     case kV0:
450       return "V0M";
451     default:
452       return "";
453   }
454 }
455 //----------------------------------------------------------------------- 
456 AliFlowEventCuts* AliFlowEventCuts::StandardCuts()
457 {
458   //make a set of standard event cuts, caller becomes owner
459   AliFlowEventCuts* cuts = new AliFlowEventCuts();
460   return cuts;
461 }
462
463 //----------------------------------------------------------------------- 
464 Int_t AliFlowEventCuts::RefMult(AliVEvent* event)
465 {
466   //calculate the reference multiplicity, if all fails return 0
467   AliESDVZERO* vzero = NULL;
468   AliESDEvent* esdevent = dynamic_cast<AliESDEvent*>(event);
469
470   if (fUseAliESDtrackCutsRefMult && esdevent)
471   {
472     //use the standard ALICE reference multiplicity with the default eta range
473     return AliESDtrackCuts::GetReferenceMultiplicity(esdevent, fRefMultMethodAliESDtrackCuts);
474   }
475
476   if (fRefMultMethod==kTPConly && !fRefMultCuts)
477   {
478     fRefMultCuts = AliFlowTrackCuts::GetStandardTPCStandaloneTrackCuts();
479     fRefMultCuts->SetEtaRange(-0.8,0.8);
480     fRefMultCuts->SetPtMin(0.15);
481   }
482   else if (fRefMultMethod==kSPDtracklets && !fRefMultCuts)
483   {
484     fRefMultCuts = new AliFlowTrackCuts("tracklet refmult cuts");
485     fRefMultCuts->SetParamType(AliFlowTrackCuts::kSPDtracklet);
486     fRefMultCuts->SetEtaRange(-0.8,0.8);
487   }
488   else if (fRefMultMethod==kV0)
489   {
490     if (!esdevent) return 0;
491     vzero=esdevent->GetVZEROData();
492     if (!vzero) return 0;
493     return TMath::Nint(vzero->GetMTotV0A()+vzero->GetMTotV0C());
494   }
495   else if (fRefMultMethod==kSPD1clusters)
496   {
497     if (!esdevent) return 0;
498     const AliMultiplicity* mult = esdevent->GetMultiplicity();
499     if (!mult) return 0;
500     return mult->GetNumberOfITSClusters(1);
501   }
502
503   Int_t refmult=0;
504   fRefMultCuts->SetEvent(event);
505   for (Int_t i=0; i<fRefMultCuts->GetNumberOfInputObjects(); i++)
506   {
507     if (fRefMultCuts->IsSelected(fRefMultCuts->GetInputObject(i),i))
508       refmult++;
509   }
510   return refmult;
511 }
512 //_____________________________________________________________________________
513 void AliFlowEventCuts::DefineHistograms()
514 {
515   //define QA histos
516   if (fQA) return;
517
518   Bool_t adddirstatus = TH1::AddDirectoryStatus();
519   TH1::AddDirectory(kFALSE);
520   fQA = new TList(); fQA->SetOwner();
521   fQA->SetName(Form("%s QA",GetName()));
522   TList* before = new TList(); before->SetOwner();
523   before->SetName("before");
524   TList* after = new TList(); after->SetOwner();
525   after->SetName("after");
526   fQA->Add(before);
527   fQA->Add(after);
528   before->Add(new TH1F("zvertex",";z;event cout",500,-15.,15.)); //0
529   after->Add(new TH1F("zvertex",";z;event cout",500,-15.,15.)); //0
530   before->Add(new TH2F("fTPCvsGlobalMult","TPC only vs Global track multiplicity;global;TPC only",500,0,2500,500,0,3500));//1
531   after->Add(new TH2F("fTPCvsGlobalMult","TPC only vs Global track multiplicity;global;TPC only",500,0,2500,500,0,3500));//1
532   TH1::AddDirectory(adddirstatus);
533 }
534
535 //---------------------------------------------------------------//
536 void AliFlowEventCuts::Browse(TBrowser* b)
537 {
538   //some browsing capabilities
539   if (fQA) b->Add(fQA);
540 }
541
542 //---------------------------------------------------------------//
543 Long64_t AliFlowEventCuts::Merge(TCollection* list)
544 {
545   //merge
546   Int_t number=0;
547   AliFlowEventCuts* obj;
548   if (!list) return 0;
549   if (list->GetEntries()<1) return 0;
550   TIter next(list);
551   while ( (obj = dynamic_cast<AliFlowEventCuts*>(next())) )
552   {
553     if (obj==this) continue;
554     TList listwrapper;
555     listwrapper.Add(obj->GetQA());
556     fQA->Merge(&listwrapper);
557     number++;
558   }
559   return number;
560 }
561