]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliEventTagCuts.cxx
Adding the EMCAL and PHOS cuts as well as some minor modifications in the code.
[u/mrichter/AliRoot.git] / STEER / AliEventTagCuts.cxx
1 /**************************************************************************
2  * Author: Panos Christakoglou.                                           *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 /* $Id$ */
15
16 //-----------------------------------------------------------------
17 //           AliEventTagCuts class
18 //   This is the class to deal with the event tag level cuts
19 //   Origin: Panos Christakoglou, UOA-CERN, Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
21
22 class AliLog;
23 class AliESD;
24
25 #include "AliEventTag.h"
26 #include "AliEventTagCuts.h"
27
28 ClassImp(AliEventTagCuts)
29
30
31 //___________________________________________________________________________
32 AliEventTagCuts::AliEventTagCuts() {
33   //Default constructor which calls the Reset method.
34   Reset();
35 }
36
37 //___________________________________________________________________________
38 AliEventTagCuts::~AliEventTagCuts() {  
39   //Defaut destructor.
40 }
41
42 //___________________________________________________________________________
43 void AliEventTagCuts::Reset() {
44   //Sets dummy values to every private member.
45   fVxFlag = kFALSE;
46   fVyFlag = kFALSE;
47   fVzFlag = kFALSE;
48   fParticipantsFlag = kFALSE;
49   fImpactParamFlag = kFALSE;
50   fPVFlag = kFALSE;
51
52   fPVzErrorFlag = kFALSE;
53   fTriggerMaskFlag = kFALSE;
54   fTriggerClusterFlag = kFALSE;
55
56   fZDCNeutron1EnergyFlag = kFALSE;
57   fZDCProton1EnergyFlag = kFALSE;
58   fZDCNeutron2EnergyFlag = kFALSE;
59   fZDCProton2EnergyFlag = kFALSE;
60   fZDCEMEnergyFlag = kFALSE;
61   fT0VertexZFlag = kFALSE;
62   fMultFlag = kFALSE;
63   fMultPosFlag = kFALSE;
64   fMultNegFlag = kFALSE;
65   fMultNeutrFlag = kFALSE;
66   fV0sFlag = kFALSE;
67   fCascadesFlag = kFALSE;
68   fkinksFlag = kFALSE;
69
70   fPMDTracksFlag = kFALSE;
71   fFMDTracksFlag = kFALSE;
72   fPHOSClustersFlag = kFALSE;
73   fEMCALClustersFlag = kFALSE;
74   fJetCandidatesFlag = kFALSE;
75
76   fMaxJetEnergyFlag = kFALSE;
77   fNHardPhotonsCandidatesFlag = kFALSE;
78   fMaxNeutralFlag = kFALSE;
79   fChargedAbove1GeVFlag = kFALSE;
80   fChargedAbove3GeVFlag = kFALSE;
81   fChargedAbove10GeVFlag = kFALSE;
82   fMuonsAbove1GeVFlag = kFALSE;
83   fMuonsAbove3GeVFlag = kFALSE;
84   fMuonsAbove10GeVFlag = kFALSE;
85   fElectronsAbove1GeVFlag = kFALSE;
86   fElectronsAbove3GeVFlag = kFALSE;
87   fElectronsAbove10GeVFlag = kFALSE;
88   fElectronsFlag = kFALSE;
89   fMuonsFlag = kFALSE;
90   fPionsFlag = kFALSE;
91   fKaonsFlag = kFALSE;
92   fProtonsFlag = kFALSE;
93   fLambdasFlag = kFALSE;
94   fPhotonFlag = kFALSE;
95   fPi0sFlag = kFALSE;
96   fNeutronsFlag = kFALSE;
97   fKaon0sFlag = kFALSE;
98   fTotalPFlag = kFALSE;
99   fMeanPtFlag = kFALSE;
100   fMaxPtFlag = kFALSE;
101   fTotalNeutralPFlag = kFALSE;
102   fMeanNeutralPtFlag = kFALSE;
103   fMaxNeutralPtFlag = kFALSE;
104   fEventPlaneAngleFlag = kFALSE;
105   fHBTRadiiFlag = kFALSE;
106   
107   fVxMin = -1000.0; fVxMax = 1000.0; 
108   fVyMin = -1000.0; fVyMax = 1000.0;  
109   fVzMin = -1000.0; fVzMax = 1000.0;
110   fParticipantsMin = -1; fParticipantMax = 10000;
111   fImpactParamMin = -1.0; fImpactParamMax = 1000.0;
112   fPrimaryVertexFlag = 1;
113
114   fPrimaryVertexZErrorMin = -10000.; fPrimaryVertexZErrorMax = 10000.;
115   fTriggerMask = 0;
116   fTriggerCluster = 0;
117  
118   fZDCNeutron1EnergyMin = -1.0; fZDCNeutron1EnergyMax = 100000.0;
119   fZDCProton1EnergyMin = -1.0; fZDCProton1EnergyMax = 100000.0;
120   fZDCNeutron2EnergyMin = -1.0; fZDCNeutron2EnergyMax = 100000.0;
121   fZDCProton2EnergyMin = -1.0; fZDCProton2EnergyMax = 100000.0;
122   fZDCEMEnergyMin = -1.0; fZDCEMEnergyMax = 100000.0;
123   fT0VertexZMin = -10000.0; fT0VertexZMax = 10000.0;  
124   fMultMin = 0; fMultMax = 100000;  
125   fMultPosMin = -1; fMultPosMax = 100000;
126   fMultNegMin = -1; fMultNegMax = 100000;
127   fMultNeutrMin = -1; fMultNeutrMax = 100000;
128   fV0sMin = -1; fV0sMax = 1000000;
129   fCascadesMin = -1; fCascadesMax = 100000;
130   fkinksMin = -1; fkinksMax = 1000000;
131
132   fPMDTracksMin = -1, fPMDTracksMax = 100000;
133   fFMDTracksMin = -1, fFMDTracksMax = 100000;
134   fPHOSClustersMin = -1, fPHOSClustersMax = 100000;
135   fEMCALClustersMin = -1, fEMCALClustersMax = 100000;
136   fJetCandidatesMin = -1, fJetCandidatesMax = 100000;
137
138   fMaxJetEnergy = -1.0; 
139   fNHardPhotonsCandidatesMin = -1; fNHardPhotonsCandidatesMax = 100000;
140   fMaxNeutralEnergy = -1.0; 
141   fChargedAbove1GeVMin = -1; fChargedAbove1GeVMax = 100000;
142   fChargedAbove3GeVMin = -1; fChargedAbove3GeVMax = 100000;
143   fChargedAbove10GeVMin = -1; fChargedAbove10GeVMax = 100000;
144   fMuonsAbove1GeVMin = -1; fMuonsAbove1GeVMax = 100000;
145   fMuonsAbove3GeVMin = -1; fMuonsAbove3GeVMax = 100000;
146   fMuonsAbove10GeVMin = -1; fMuonsAbove10GeVMax = 100000; 
147   fElectronsAbove1GeVMin = -1; fElectronsAbove1GeVMax = 100000;
148   fElectronsAbove3GeVMin = -1; fElectronsAbove3GeVMax = 100000;
149   fElectronsAbove10GeVMin = -1; fElectronsAbove10GeVMax = 100000;
150   fElectronsMin = -1; fElectronsMax = 100000;
151   fMuonsMin = -1; fMuonsMax = 100000;
152   fPionsMin = -1; fPionsMax = 100000;
153   fKaonsMin = -1; fKaonsMax = 100000;
154   fProtonsMin = -1; fProtonsMax = 100000;
155   fLambdasMin = -1; fLambdasMax = 100000;
156   fPhotonsMin = -1; fPhotonsMax = 100000;
157   fPi0sMin = -1; fPi0sMax = 100000; 
158   fNeutronsMin = -1; fNeutronsMax = 100000; 
159   fKaon0sMin = -1; fKaon0sMax = 100000; 
160   fTotalPMin = -1.0; fTotalPMax = 1000000.0;
161   fMeanPtMin = -1.0; fMeanPtMax = 100000.0;
162   fMaxPt = -1.0; fTotalNeutralPMin = -1.0;
163   fTotalNeutralPMax = 1000000.0;   
164   fMeanNeutralPtMin = -1.0; fMeanNeutralPtMax = 1000000.0; 
165   fMaxNeutralPt = -1.0; 
166   fEventPlaneAngleMin = -10000000.0; fEventPlaneAngleMax = 10000000.0; 
167   fHBTRadiiMin = -1.0; fHBTRadiiMax = 100000.0; 
168 }
169
170 //___________________________________________________________________________
171 void AliEventTagCuts::SetPrimaryVertexXRange(Float_t r1, Float_t r2) {
172   //Sets the primary vertex x range 
173   //and the corresponding flag to kTRUE if the cut is used.
174   fVxMin = r1;
175   fVxMax = r2; 
176   fVxFlag = kTRUE;
177 }
178
179 //___________________________________________________________________________
180 void AliEventTagCuts::SetPrimaryVertexYRange(Float_t r1, Float_t r2) {
181   //Sets the primary vertex y range 
182   //and the corresponding flag to kTRUE if the cut is used.
183   fVyMin = r1;
184   fVyMax = r2; 
185   fVyFlag = kTRUE;
186 }
187
188 //___________________________________________________________________________
189 void AliEventTagCuts::SetPrimaryVertexZRange(Float_t r1, Float_t r2) {
190   //Sets the primary vertex z range 
191   //and the corresponding flag to kTRUE if the cut is used.
192   fVzMin = r1;
193   fVzMax = r2; 
194   fVzFlag = kTRUE;
195 }
196
197 //___________________________________________________________________________
198 void AliEventTagCuts::SetMultiplicityRange(Int_t n1, Int_t n2) {
199   //Sets the primary multiplicity range 
200   //and the corresponding flag to kTRUE if the cut is used.
201   fMultMin = n1;
202   fMultMax = n2;
203   fMultFlag = kTRUE;
204 }
205
206 //___________________________________________________________________________
207 void AliEventTagCuts::SetParticipantsRange(Int_t i1, Int_t i2) {
208   //Sets the number of participants range 
209   //and the corresponding flag to kTRUE if the cut is used.
210   fParticipantsMin = i1;
211   fParticipantMax = i2;
212   fParticipantsFlag = kTRUE;
213 }
214
215 //___________________________________________________________________________
216 void AliEventTagCuts::SetImpactParamRange(Float_t r1, Float_t r2) {
217   //Sets the impact parameter range 
218   //and the corresponding flag to kTRUE if the cut is used.
219   fImpactParamMin = r1;
220   fImpactParamMax = r2;
221   fImpactParamFlag = kTRUE;
222 }
223  
224
225 //___________________________________________________________________________
226 void AliEventTagCuts::SetPrimaryVertexFlag(Int_t i) {
227   //Sets the primary vertex flag cut 
228   //and the corresponding flag to kTRUE if the cut is used.
229   fPrimaryVertexFlag = i;
230   fPVFlag = kTRUE;
231 }
232
233 //___________________________________________________________________________
234 void AliEventTagCuts::SetZDCNeutr1Range(Float_t r1, Float_t r2) {
235   //Sets the ZDC's neutron energy range 
236   //and the corresponding flag to kTRUE if the cut is used.
237   fZDCNeutron1EnergyMin = r1;
238   fZDCNeutron1EnergyMax = r2;
239   fZDCNeutron1EnergyFlag = kTRUE;
240 }
241
242 //___________________________________________________________________________
243 void AliEventTagCuts::SetZDCProt1Range(Float_t r1, Float_t r2) {
244   //Sets the ZDC's proton energy range 
245   //and the corresponding flag to kTRUE if the cut is used.
246   fZDCProton1EnergyMin = r1;
247   fZDCProton1EnergyMax = r2;
248   fZDCProton1EnergyFlag = kTRUE;
249 }
250 //___________________________________________________________________________
251 void AliEventTagCuts::SetZDCNeutr2Range(Float_t r1, Float_t r2) {
252   //Sets the ZDC's neutron energy range 
253   //and the corresponding flag to kTRUE if the cut is used.
254   fZDCNeutron2EnergyMin = r1;
255   fZDCNeutron2EnergyMax = r2;
256   fZDCNeutron2EnergyFlag = kTRUE;
257 }
258 //___________________________________________________________________________
259 void AliEventTagCuts::SetZDCProt2Range(Float_t r1, Float_t r2) {
260   //Sets the ZDC's proton energy range 
261   //and the corresponding flag to kTRUE if the cut is used.
262   fZDCProton2EnergyMin = r1;
263   fZDCProton2EnergyMax = r2;
264   fZDCProton2EnergyFlag = kTRUE;
265 }
266 //___________________________________________________________________________
267 void AliEventTagCuts::SetZDCEMRange(Float_t r1, Float_t r2) {
268   //Sets the ZDC's e/m energy range 
269   //and the corresponding flag to kTRUE if the cut is used.
270   fZDCEMEnergyMin = r1;
271   fZDCEMEnergyMax = r2;
272   fZDCEMEnergyFlag = kTRUE;
273 }
274
275 //___________________________________________________________________________
276 void AliEventTagCuts::SetT0VertexZRange(Float_t r1, Float_t r2) {
277   //Sets the T0's Vz range 
278   //and the corresponding flag to kTRUE if the cut is used.
279   fT0VertexZMin = r1;
280   fT0VertexZMax = r2;
281   fT0VertexZFlag = kTRUE;
282 }
283
284 //___________________________________________________________________________
285 void AliEventTagCuts::SetPosMultiplicityRange(Int_t n1, Int_t n2) {
286   //Sets the positive multiplicity range 
287   //and the corresponding flag to kTRUE if the cut is used.
288   fMultPosMin = n1;
289   fMultPosMax = n2;
290   fMultPosFlag = kTRUE;
291 }
292
293
294 //___________________________________________________________________________
295 void AliEventTagCuts::SetNegMultiplicityRange(Int_t n1, Int_t n2) {
296   //Sets the negative multiplicity range 
297   //and the corresponding flag to kTRUE if the cut is used.
298   fMultNegMin = n1;
299   fMultNegMax = n2;
300   fMultNegFlag = kTRUE;
301 }
302
303
304 //___________________________________________________________________________
305 void AliEventTagCuts::SetNeutrMultiplicityRange(Int_t n1, Int_t n2) {
306   //Sets the neutral particle multiplicity range 
307   //and the corresponding flag to kTRUE if the cut is used.
308   fMultNeutrMin = n1;
309   fMultNeutrMax = n2;
310   fMultNeutrFlag = kTRUE;
311 }
312
313 //___________________________________________________________________________
314 void AliEventTagCuts::SetV0sRange(Int_t n1, Int_t n2) {
315   //Sets the v0s multiplicity range 
316   //and the corresponding flag to kTRUE if the cut is used.
317   fV0sMin = n1;
318   fV0sMax = n2;
319   fV0sFlag = kTRUE;
320 }
321
322 //___________________________________________________________________________
323 void AliEventTagCuts::SetCascadesRange(Int_t n1, Int_t n2) {
324   //Sets the cascades multiplicity range 
325   //and the corresponding flag to kTRUE if the cut is used.
326   fCascadesMin = n1;
327   fCascadesMax = n2;
328   fCascadesFlag = kTRUE;
329 }
330
331 //___________________________________________________________________________
332 void AliEventTagCuts::SetKinksRange(Int_t n1, Int_t n2) {
333   //Sets the kinks multipliicity range 
334   //and the corresponding flag to kTRUE if the cut is used.
335   fkinksMin = n1;
336   fkinksMax = n2;
337   fkinksFlag = kTRUE;
338 }
339
340 //___________________________________________________________________________
341 void AliEventTagCuts::SetMaxJetEnergy(Float_t r1) {
342   //Sets the lower limit of the maximum jet energy
343   //and the corresponding flag to kTRUE if the cut is used.
344   fMaxJetEnergy = r1; 
345   fMaxJetEnergyFlag = kTRUE;
346 }
347 //___________________________________________________________________________
348 void AliEventTagCuts::SetMaxNeutralEnergy(Float_t r1) {
349   //Sets the lower limit of the maximum neutral jet energy
350   //and the corresponding flag to kTRUE if the cut is used.
351   fMaxNeutralEnergy = r1; 
352   fMaxNeutralFlag = kTRUE;
353 }
354 //___________________________________________________________________________
355 void AliEventTagCuts::SetHardPhotonsRange(Int_t i1, Int_t i2) {
356   //Sets the hard photons multiplicity range
357   //and the corresponding flag to kTRUE if the cut is used.
358   fNHardPhotonsCandidatesMin = i1;
359   fNHardPhotonsCandidatesMax = i2;
360   fNHardPhotonsCandidatesFlag = kTRUE;
361
362
363 //___________________________________________________________________________
364 void AliEventTagCuts::SetNChargedAbove1GeVRange(Int_t i1, Int_t i2) {
365   //Sets the number of charged above 1GeV range
366   //and the corresponding flag to kTRUE if the cut is used.
367   fChargedAbove1GeVMin = i1;
368   fChargedAbove1GeVMax = i2;
369   fChargedAbove1GeVFlag = kTRUE;
370 }
371
372 //___________________________________________________________________________
373  void AliEventTagCuts::SetNChargedAbove3GeVRange(Int_t i1, Int_t i2) {
374   //Sets the number of charged above 3GeV range
375   //and the corresponding flag to kTRUE if the cut is used.
376   fChargedAbove3GeVMin = i1;
377   fChargedAbove3GeVMax = i2;
378   fChargedAbove3GeVFlag = kTRUE;
379 }
380
381
382 //___________________________________________________________________________
383 void AliEventTagCuts::SetNChargedAbove10GeVRange(Int_t i1, Int_t i2) {
384   //Sets the number of charged above 10GeV range
385   //and the corresponding flag to kTRUE if the cut is used.
386   fChargedAbove10GeVMin = i1;
387   fChargedAbove10GeVMax = i2;
388   fChargedAbove10GeVFlag = kTRUE;
389 }
390
391
392 //___________________________________________________________________________
393 void AliEventTagCuts::SetNMuonsAbove1GeVRange(Int_t i1, Int_t i2) {
394   //Sets the number of muons above 1GeV range
395   //and the corresponding flag to kTRUE if the cut is used.
396   fMuonsAbove1GeVMin = i1;
397   fMuonsAbove1GeVMax = i2;
398   fMuonsAbove1GeVFlag = kTRUE;
399 }
400
401
402 //___________________________________________________________________________
403 void AliEventTagCuts::SetNMuonsAbove3GeVRange(Int_t i1, Int_t i2) {
404   //Sets the number of muons above 3GeV range
405   //and the corresponding flag to kTRUE if the cut is used.
406   fMuonsAbove3GeVMin = i1;
407   fMuonsAbove3GeVMax = i2;
408   fMuonsAbove3GeVFlag = kTRUE;
409
410
411 //___________________________________________________________________________
412 void AliEventTagCuts::SetNMuonsAbove10GeVRange(Int_t i1, Int_t i2) {
413   //Sets the number of muons above 10GeV range
414   //and the corresponding flag to kTRUE if the cut is used.
415   fMuonsAbove10GeVMin = i1;
416   fMuonsAbove10GeVMax = i2; 
417   fMuonsAbove10GeVFlag = kTRUE;
418 }
419
420
421 //___________________________________________________________________________
422 void AliEventTagCuts::SetNElectronsAbove1GeVRange(Int_t i1, Int_t i2) {
423   //Sets the number of electrons above 1GeV range
424   //and the corresponding flag to kTRUE if the cut is used.
425   fElectronsAbove1GeVMin = i1;
426   fElectronsAbove1GeVMax = i2;
427   fElectronsAbove1GeVFlag = kTRUE;
428 }
429
430 //___________________________________________________________________________
431 void AliEventTagCuts::SetNElectronsAbove3GeVRange(Int_t i1, Int_t i2) {
432   //Sets the number of electrons above 3GeV range
433   //and the corresponding flag to kTRUE if the cut is used.
434   fElectronsAbove3GeVMin = i1;
435   fElectronsAbove3GeVMax = i2;
436   fElectronsAbove3GeVFlag = kTRUE;
437 }
438
439 //___________________________________________________________________________
440 void AliEventTagCuts::SetNElectronsAbove10GeVRange(Int_t i1, Int_t i2) {  
441   //Sets the number of electrons above 10GeV range
442   //and the corresponding flag to kTRUE if the cut is used.
443   fElectronsAbove10GeVMin = i1;
444   fElectronsAbove10GeVMax = i2;
445   fElectronsAbove10GeVFlag = kTRUE;
446 }
447 //___________________________________________________________________________
448 void AliEventTagCuts::SetNElectronRange(Int_t n1, Int_t n2) {
449   //Sets the electron multiplicity range
450   //and the corresponding flag to kTRUE if the cut is used.
451   fElectronsMin = n1;
452   fElectronsMax = n2;
453   fElectronsFlag = kTRUE;
454 }
455 //___________________________________________________________________________
456 void AliEventTagCuts::SetNMuonRange(Int_t n1, Int_t n2) {
457   //Sets the muon multiplicity range
458   //and the corresponding flag to kTRUE if the cut is used.
459   fMuonsMin = n1;
460   fMuonsMax = n2;
461   fMuonsFlag = kTRUE;
462 }
463
464 //___________________________________________________________________________
465 void AliEventTagCuts::SetNPionRange(Int_t n1, Int_t n2) {
466   //Sets the pion multiplicity range
467   //and the corresponding flag to kTRUE if the cut is used.
468   fPionsMin = n1;
469   fPionsMax = n2;
470   fPionsFlag = kTRUE;
471
472
473 //___________________________________________________________________________
474 void AliEventTagCuts::SetNKaonRange(Int_t n1, Int_t n2) {
475   //Sets the kaon multiplicity range
476   //and the corresponding flag to kTRUE if the cut is used.
477   fKaonsMin = n1;
478   fKaonsMax = n2;
479   fKaonsFlag = kTRUE;
480 }
481
482 //___________________________________________________________________________
483 void AliEventTagCuts::SetNProtonRange(Int_t n1, Int_t n2) {
484   //Sets the proton multiplicity range
485   //and the corresponding flag to kTRUE if the cut is used.
486   fProtonsMin = n1;
487   fProtonsMax = n2;
488   fProtonsFlag = kTRUE;
489
490
491 //___________________________________________________________________________
492 void AliEventTagCuts::SetNLambdaRange(Int_t n1, Int_t n2) {
493   //Sets the lambda multiplicity range
494   //and the corresponding flag to kTRUE if the cut is used.
495   fLambdasMin = n1;
496   fLambdasMax = n2;
497   fLambdasFlag = kTRUE;
498
499 //___________________________________________________________________________
500 void AliEventTagCuts::SetNPhotonRange(Int_t n1, Int_t n2) {
501   //Sets the photon multiplicity range
502   //and the corresponding flag to kTRUE if the cut is used.
503   fPhotonsMin = n1;
504   fPhotonsMax = n2;
505   fPhotonFlag = kTRUE;
506
507 //___________________________________________________________________________
508 void AliEventTagCuts::SetNPi0Range(Int_t n1, Int_t n2) {
509   //Sets the pi0 multiplicity range
510   //and the corresponding flag to kTRUE if the cut is used.
511   fPi0sMin = n1;
512   fPi0sMax = n2; 
513   fPi0sFlag = kTRUE;
514 }  
515
516 //___________________________________________________________________________
517 void AliEventTagCuts::SetNNeutronRange(Int_t n1, Int_t n2) {
518   //Sets the neutron multiplicity range
519   //and the corresponding flag to kTRUE if the cut is used.
520   fNeutronsMin = n1;
521   fNeutronsMax = n2; 
522   fNeutronsFlag = kTRUE;
523 }
524
525 //___________________________________________________________________________
526 void AliEventTagCuts::SetNKaon0Range(Int_t n1, Int_t n2) {  
527   //Sets the K0s multiplicity range
528   //and the corresponding flag to kTRUE if the cut is used.
529   fKaon0sMin = n1;
530   fKaon0sMax = n2; 
531   fKaon0sFlag = kTRUE;
532 }
533
534 //___________________________________________________________________________
535 void AliEventTagCuts::SetTotalPRange(Float_t r1, Float_t r2) {
536   //Sets the total momentum range
537   //and the corresponding flag to kTRUE if the cut is used.
538   fTotalPMin = r1;
539   fTotalPMax = r2;
540   fTotalPFlag = kTRUE;
541 }
542
543 //___________________________________________________________________________
544 void AliEventTagCuts::SetMeanPtRange(Float_t r1, Float_t r2) {
545   //Sets the mean Pt range
546   //and the corresponding flag to kTRUE if the cut is used.
547   fMeanPtMin = r1;
548   fMeanPtMax = r2;
549   fMeanPtFlag = kTRUE;
550 }  
551
552 //___________________________________________________________________________
553 void AliEventTagCuts::SetMaxPt(Float_t r1) {
554   //Sets the lower limit of the max Pt value
555   //and the corresponding flag to kTRUE if the cut is used.
556   fMaxPt = r1; 
557   fMaxPtFlag = kTRUE;
558 }
559
560 //___________________________________________________________________________
561 void AliEventTagCuts::SetTotalNeutralPRange(Float_t r1, Float_t r2) {  
562   //Sets the total momentum of neutral particles range
563   //and the corresponding flag to kTRUE if the cut is used.
564   fTotalNeutralPMin =r1 ;
565   fTotalNeutralPMax = r2;  
566   fTotalNeutralPFlag = kTRUE;
567 }
568 //___________________________________________________________________________
569 void AliEventTagCuts::SetMeanNeutralPtPRange(Float_t r1, Float_t r2) {  
570   //Sets the mean Pt of neutral particles range
571   //and the corresponding flag to kTRUE if the cut is used.
572   fMeanNeutralPtMin = r1;
573   fMeanNeutralPtMax = r2; 
574   fMeanNeutralPtFlag = kTRUE;
575
576 //___________________________________________________________________________
577 void AliEventTagCuts::SetMaxNeutralPt(Float_t r1) {  
578   //Sets the lower limit of the maximum Pt of neutral particles
579   //and the corresponding flag to kTRUE if the cut is used.
580   fMaxNeutralPt = r1; 
581   fMaxNeutralPtFlag = kTRUE;
582 }
583
584 //___________________________________________________________________________
585 void AliEventTagCuts::SetEvPlaneAngleRange(Float_t r1, Float_t r2) {
586   //Sets the event plane range
587   //and the corresponding flag to kTRUE if the cut is used.
588   fEventPlaneAngleMin = r1;
589   fEventPlaneAngleMax = r2; 
590   fEventPlaneAngleFlag = kTRUE;
591 }
592
593 //___________________________________________________________________________
594 void AliEventTagCuts::SetHBTRadiiRange(Float_t r1, Float_t r2) {
595   //Sets the HBT radii range
596   //and the corresponding flag to kTRUE if the cut is used.
597   fHBTRadiiMin = r1;
598   fHBTRadiiMax = r2; 
599   fHBTRadiiFlag = kTRUE;
600 }
601
602 //___________________________________________________________________________
603 Bool_t AliEventTagCuts::IsAccepted(AliEventTag *EvTag) const {
604   //Returns true if the event is accepted otherwise false.
605   if(fVzFlag)
606     if((EvTag->GetVertexZ() < fVzMin) || (EvTag->GetVertexZ() > fVzMax))
607       return kFALSE;
608   
609   if(fVyFlag)
610     if((EvTag->GetVertexY() < fVyMin) || (EvTag->GetVertexY() > fVyMax))
611       return kFALSE;
612   
613   if(fVxFlag)
614     if((EvTag->GetVertexX() < fVxMin) || (EvTag->GetVertexX() > fVxMax))
615       return kFALSE;
616   
617   if(fParticipantsFlag)
618     if((EvTag->GetNumOfParticipants() < fParticipantsMin) || (EvTag->GetNumOfParticipants() > fParticipantMax))
619       return kFALSE; 
620   
621   if(fImpactParamFlag)
622     if((EvTag->GetImpactParameter() < fImpactParamMin) || (EvTag->GetImpactParameter() > fImpactParamMax))
623       return kFALSE; 
624   
625   if(fPVFlag)
626     if((EvTag->GetVertexFlag() != fPrimaryVertexFlag))
627       return kFALSE; 
628   
629   if(fPVzErrorFlag)
630     if((EvTag->GetVertexZError() < fPrimaryVertexZErrorMin) || (EvTag->GetVertexZError() > fPrimaryVertexZErrorMax))
631       return kFALSE; 
632   if(fTriggerMaskFlag)
633     if((EvTag->GetTriggerMask() != fTriggerMask))
634       return kFALSE; 
635   if(fTriggerClusterFlag)
636     if((EvTag->GetTriggerMask() != fTriggerMask))
637       return kFALSE; 
638
639   if(fZDCNeutron1EnergyFlag)
640     if((EvTag->GetZDCNeutron1Energy() < fZDCNeutron1EnergyMin) || (EvTag->GetZDCNeutron1Energy() > fZDCNeutron1EnergyMax))
641       return kFALSE; 
642   
643   if(fZDCProton1EnergyFlag)
644     if((EvTag->GetZDCProton1Energy() < fZDCProton1EnergyMin) || (EvTag->GetZDCProton1Energy() > fZDCProton1EnergyMax))
645       return kFALSE; 
646   
647   if(fZDCNeutron2EnergyFlag)
648     if((EvTag->GetZDCNeutron2Energy() < fZDCNeutron2EnergyMin) || (EvTag->GetZDCNeutron2Energy() > fZDCNeutron2EnergyMax))
649       return kFALSE; 
650   
651   if(fZDCProton2EnergyFlag)
652     if((EvTag->GetZDCProton2Energy() < fZDCProton2EnergyMin) || (EvTag->GetZDCProton2Energy() > fZDCProton2EnergyMax))
653       return kFALSE; 
654   
655   if(fZDCEMEnergyFlag)
656     if((EvTag->GetZDCEMEnergy() < fZDCEMEnergyMin) || (EvTag->GetZDCEMEnergy() > fZDCEMEnergyMax))
657       return kFALSE; 
658   
659   if(fT0VertexZFlag)
660     if((EvTag->GetT0VertexZ() < fT0VertexZMin) || (EvTag->GetT0VertexZ() > fT0VertexZMax))
661       return kFALSE; 
662   
663   if(fMultFlag)
664     if((EvTag->GetNumOfTracks() < fMultMin) || (EvTag->GetNumOfTracks() > fMultMax))
665       return kFALSE; 
666   if(fMultPosFlag)
667     if((EvTag->GetNumOfPosTracks() < fMultPosMin) || (EvTag->GetNumOfPosTracks() > fMultPosMax))
668       return kFALSE; 
669   
670   if(fMultNegFlag)
671     if((EvTag->GetNumOfNegTracks() < fMultNegMin) || (EvTag->GetNumOfNegTracks() > fMultNegMax))
672       return kFALSE; 
673   
674   if(fMultNeutrFlag)
675     if((EvTag->GetNumOfNeutrTracks() < fMultNeutrMin) || (EvTag->GetNumOfNeutrTracks() > fMultNeutrMax))
676       return kFALSE; 
677   
678   if(fV0sFlag)
679     if((EvTag->GetNumOfV0s() < fV0sMin) || (EvTag->GetNumOfV0s() > fV0sMax))
680       return kFALSE; 
681   
682   if(fCascadesFlag)
683     if((EvTag->GetNumOfCascades() < fCascadesMin) || (EvTag->GetNumOfCascades() > fCascadesMax))
684       return kFALSE; 
685   
686   if(fkinksFlag)
687     if((EvTag->GetNumOfKinks() < fkinksMin) || (EvTag->GetNumOfKinks() > fkinksMax))
688       return kFALSE; 
689
690
691   if(fPMDTracksFlag)
692     if((EvTag->GetNumOfPMDTracks() < fPMDTracksMin) || (EvTag->GetNumOfPMDTracks() > fPMDTracksMax))
693       return kFALSE; 
694   if(fFMDTracksFlag)
695     if((EvTag->GetNumOfFMDTracks() < fFMDTracksMin) || (EvTag->GetNumOfFMDTracks() > fFMDTracksMax))
696       return kFALSE; 
697   if(fPHOSClustersFlag)
698     if((EvTag->GetNumOfPHOSClusters() < fPHOSClustersMin) || (EvTag->GetNumOfPHOSClusters() > fPHOSClustersMax))
699       return kFALSE; 
700   if(fEMCALClustersFlag)
701     if((EvTag->GetNumOfEMCALClusters() < fEMCALClustersMin) || (EvTag->GetNumOfEMCALClusters() > fEMCALClustersMax))
702       return kFALSE; 
703   if(fJetCandidatesFlag)
704     if((EvTag->GetNumOfJetCandidates() < fJetCandidatesMin) || (EvTag->GetNumOfJetCandidates() > fJetCandidatesMax))
705       return kFALSE; 
706
707
708   if(fMaxJetEnergyFlag)
709     if((EvTag->GetMaxJetEnergy() < fMaxJetEnergy))
710       return kFALSE; 
711   
712   if(fNHardPhotonsCandidatesFlag)
713     if((EvTag->GetNumOfHardPhotonsCandidates() < fNHardPhotonsCandidatesMin) || (EvTag->GetNumOfHardPhotonsCandidates() > fNHardPhotonsCandidatesMax))
714       return kFALSE; 
715   
716   if(fMaxNeutralFlag)
717     if((EvTag->GetMaxNeutralEnergy() < fMaxNeutralEnergy))
718       return kFALSE; 
719   
720   if(fChargedAbove1GeVFlag)
721     if((EvTag->GetNumOfChargedAbove1GeV() < fChargedAbove1GeVMin) || (EvTag->GetNumOfChargedAbove1GeV() > fChargedAbove1GeVMax))
722       return kFALSE; 
723   
724   if(fChargedAbove3GeVFlag)
725     if((EvTag->GetNumOfChargedAbove3GeV() < fChargedAbove3GeVMin) || (EvTag->GetNumOfChargedAbove3GeV() > fChargedAbove3GeVMax))
726       return kFALSE; 
727   
728   if(fChargedAbove10GeVFlag)
729     if((EvTag->GetNumOfChargedAbove10GeV() < fChargedAbove10GeVMin) || (EvTag->GetNumOfChargedAbove10GeV() > fChargedAbove10GeVMax))
730       return kFALSE; 
731   
732   if(fMuonsAbove1GeVFlag)
733     if((EvTag->GetNumOfMuonsAbove1GeV() < fMuonsAbove1GeVMin) || (EvTag->GetNumOfMuonsAbove1GeV() > fMuonsAbove1GeVMax))
734       return kFALSE; 
735   
736   if(fMuonsAbove3GeVFlag)
737     if((EvTag->GetNumOfMuonsAbove3GeV() < fMuonsAbove3GeVMin) || (EvTag->GetNumOfMuonsAbove3GeV() > fMuonsAbove3GeVMax))
738       return kFALSE; 
739   
740   if(fMuonsAbove10GeVFlag)
741     if((EvTag->GetNumOfMuonsAbove10GeV() < fMuonsAbove10GeVMin) || (EvTag->GetNumOfMuonsAbove10GeV() > fMuonsAbove10GeVMax))
742       return kFALSE; 
743   
744   if(fElectronsAbove1GeVFlag)
745     if((EvTag->GetNumOfElectronsAbove1GeV()  < fElectronsAbove1GeVMin) || (EvTag->GetNumOfElectronsAbove1GeV()  > fElectronsAbove1GeVMax))
746       return kFALSE; 
747   
748   if(fElectronsAbove3GeVFlag)
749     if((EvTag->GetNumOfElectronsAbove3GeV() < fElectronsAbove3GeVMin) || (EvTag->GetNumOfElectronsAbove3GeV() > fElectronsAbove3GeVMax))
750       return kFALSE; 
751   
752   if(fElectronsAbove10GeVFlag)
753     if((EvTag->GetNumOfElectronsAbove10GeV() < fElectronsAbove10GeVMin) || (EvTag->GetNumOfElectronsAbove10GeV() > fElectronsAbove10GeVMax))
754       return kFALSE; 
755   
756   if(fElectronsFlag)
757     if((EvTag->GetNumOfElectrons() < fElectronsMin) || (EvTag->GetNumOfElectrons() > fElectronsMax))
758       return kFALSE; 
759   
760   if(fMuonsFlag)
761     if((EvTag->GetNumOfMuons() < fMuonsMin) || (EvTag->GetNumOfMuons() > fMuonsMax))
762       return kFALSE; 
763   
764   if(fPionsFlag)
765     if((EvTag->GetNumOfPions() < fPionsMin) || (EvTag->GetNumOfPions() > fPionsMax))
766       return kFALSE; 
767   
768   if(fKaonsFlag)
769     if((EvTag->GetNumOfKaons() < fKaonsMin) || (EvTag->GetNumOfKaons() > fKaonsMax))
770       return kFALSE; 
771   
772   if(fProtonsFlag)
773     if((EvTag->GetNumOfProtons() < fProtonsMin) || (EvTag->GetNumOfProtons() > fProtonsMax))
774       return kFALSE; 
775   
776   if(fLambdasFlag)
777     if((EvTag->GetNumOfLambdas() < fLambdasMin) || (EvTag->GetNumOfLambdas() > fLambdasMax))
778       return kFALSE; 
779   
780   if(fPhotonFlag)
781     if((EvTag->GetNumOfPhotons() < fPhotonsMin) || (EvTag->GetNumOfPhotons() > fPhotonsMax))
782       return kFALSE; 
783   
784   if(fPi0sFlag)
785     if((EvTag->GetNumOfPi0s() < fPi0sMin) || (EvTag->GetNumOfPi0s() > fPi0sMax))
786       return kFALSE; 
787   
788   if(fNeutronsFlag)
789     if((EvTag->GetNumOfNeutrons() < fNeutronsMin) || (EvTag->GetNumOfNeutrons() > fNeutronsMax))
790       return kFALSE; 
791   
792   if(fKaon0sFlag)
793     if((EvTag->GetNumOfKaon0s() < fKaon0sMin) || (EvTag->GetNumOfKaon0s() > fKaon0sMax))
794       return kFALSE; 
795   
796   if(fTotalPFlag)
797     if((EvTag->GetTotalMomentum() < fTotalPMin) || (EvTag->GetTotalMomentum() > fTotalPMax))
798       return kFALSE; 
799   
800   if(fMeanPtFlag)
801     if((EvTag->GetMeanPt() < fMeanPtMin) || (EvTag->GetMeanPt() > fMeanPtMax))
802       return kFALSE; 
803   
804   if(fMaxPtFlag)
805     if((EvTag->GetMaxPt() < fMaxPt))
806       return kFALSE; 
807   
808   if(fTotalNeutralPFlag)
809     if((EvTag->GetNeutralTotalMomentum() < fTotalNeutralPMin) || (EvTag->GetNeutralTotalMomentum() > fTotalNeutralPMax))
810       return kFALSE; 
811   
812   if(fMeanNeutralPtFlag)
813     if((EvTag->GetNeutralMeanPt() < fMeanNeutralPtMin) || (EvTag->GetNeutralMeanPt() >fMeanNeutralPtMax ))
814       return kFALSE; 
815   
816   if(fMaxNeutralPtFlag)
817     if((EvTag->GetNeutralMaxPt() < fMaxNeutralPt))
818       return kFALSE; 
819   
820   if(fEventPlaneAngleFlag)
821     if((EvTag->GetEventPlaneAngle() < fEventPlaneAngleMin) || (EvTag->GetEventPlaneAngle() > fEventPlaneAngleMax))
822       return kFALSE; 
823   
824   if(fHBTRadiiFlag)
825     if((EvTag->GetHBTRadii() < fHBTRadiiMin) || (EvTag->GetHBTRadii() > fHBTRadiiMax))
826       return kFALSE; 
827   
828   return kTRUE;
829 }