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