]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnValue.cxx
fix for Coverity (B.Hippolyte)
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnValue.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 ////////////////////////////////////////////////////////////////////////////////
17 //
18 //  This class contains all code which is used to compute any of the values
19 //  which can be of interest within a resonance analysis. Besides the obvious
20 //  invariant mass, it allows to compute other utility values on all possible
21 //  targets, in order to allow a wide spectrum of binning and checks.
22 //  When needed, this object can also define a binning in the variable which
23 //  it is required to compute, which is used for initializing axes of output
24 //  histograms (see AliRsnFunction).
25 //  The value computation requires this object to be passed the object whose
26 //  informations will be used. This object can be of any allowed input type
27 //  (track, pair, event), then this class must inherit from AliRsnTarget.
28 //  Then, when value computation is attempted, a check on target type is done
29 //  and computation is successful only if expected target matches that of the
30 //  passed object.
31 //  In some cases, the value computation can require a support external object,
32 //  which must then be passed to this class. It can be of any type inheriting
33 //  from TObject.
34 //
35 //  authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
36 //           M. Vala (martin.vala@cern.ch)
37 //
38 ////////////////////////////////////////////////////////////////////////////////
39
40 #include <Riostream.h>
41 #include "AliESDtrackCuts.h"
42 #include "AliESDpid.h"
43 #include "AliAODPid.h"
44 #include "AliCentrality.h"
45
46 #include "AliRsnEvent.h"
47 #include "AliRsnDaughter.h"
48 #include "AliRsnMother.h"
49 #include "AliRsnPairDef.h"
50 #include "AliRsnDaughterDef.h"
51
52 #include "AliRsnValue.h"
53
54 ClassImp(AliRsnValue)
55
56 //_____________________________________________________________________________
57 AliRsnValue::AliRsnValue() :
58    AliRsnTarget(),
59    fComputedValue(0),
60    fValueType(kValueTypes),
61    fBinArray(0),
62    fSupportObject(0x0)
63 {
64 //
65 // Default constructor without arguments.
66 // Initialize data members to meaningless values.
67 // This method is provided for ROOT streaming,
68 // but should never be used directly by a user.
69 //
70 }
71
72 //_____________________________________________________________________________
73 AliRsnValue::AliRsnValue
74 (const char *name, EValueType type, Int_t nbins, Double_t min, Double_t max) :
75    AliRsnTarget(name, TargetType(type)),
76    fComputedValue(0.0),
77    fValueType(type),
78    fBinArray(0),
79    fSupportObject(0x0)
80 {
81 //
82 // Main constructor (version 1).
83 // This constructor defines in meaningful way all data members,
84 // and defined a fixed binnings, subdividing the specified interval
85 // into that many bins as specified in the integer argument.
86 // ---
87 // This method is also the entry point for all instances
88 // of this class which don't need to do binning (e.g.: TNtuple inputs),
89 // since arguments 3 to 5 have default values which don't create any
90 // binning array, in order not to allocate memory when this is useless.
91 //
92
93    SetBins(nbins, min, max);
94 }
95
96 //_____________________________________________________________________________
97 AliRsnValue::AliRsnValue
98 (const char *name, EValueType type, Double_t min, Double_t max, Double_t step) :
99    AliRsnTarget(name, TargetType(type)),
100    fComputedValue(0.0),
101    fValueType(type),
102    fBinArray(0),
103    fSupportObject(0x0)
104 {
105 //
106 // Main constructor (version 2).
107 // This constructor defines in meaningful way all data members
108 // and creates enough equal bins of the specified size to cover
109 // the required interval.
110 //
111
112    SetBins(min, max, step);
113 }
114
115 //_____________________________________________________________________________
116 AliRsnValue::AliRsnValue
117 (const char *name, EValueType type, Int_t nbins, Double_t *array) :
118    AliRsnTarget(name, TargetType(type)),
119    fComputedValue(0.0),
120    fValueType(type),
121    fBinArray(0),
122    fSupportObject(0x0)
123 {
124 //
125 // Main constructor (version 3).
126 // This constructor defines in meaningful way all data members
127 // and creates a set of variable bins delimited by the passed array.
128 //
129
130    SetBins(nbins, array);
131 }
132
133 //_____________________________________________________________________________
134 AliRsnValue::AliRsnValue(const AliRsnValue& copy) :
135    AliRsnTarget(copy),
136    fComputedValue(copy.fComputedValue),
137    fValueType(copy.fValueType),
138    fBinArray(copy.fBinArray),
139    fSupportObject(copy.fSupportObject)
140 {
141 //
142 // Copy constructor.
143 // Duplicates the binning array and copies all settings.
144 //
145 }
146
147 //_____________________________________________________________________________
148 AliRsnValue& AliRsnValue::operator=(const AliRsnValue& copy)
149 {
150 //
151 // Assignment operator.
152 // Works like copy constructor.
153 //
154
155    AliRsnTarget::operator=(copy);
156
157    fComputedValue = copy.fComputedValue;
158    fBinArray = copy.fBinArray;
159    fSupportObject = copy.fSupportObject;
160    fValueType = copy.fValueType;
161
162    return (*this);
163 }
164
165 //_____________________________________________________________________________
166 void AliRsnValue::SetBins(Int_t nbins, Double_t min, Double_t max)
167 {
168 //
169 // Set binning for the axis in equally spaced bins
170 // where the number of bins, minimum and maximum are given.
171 //
172
173    if (!nbins) {
174       fBinArray.Set(0);
175       return;
176    }
177
178    fBinArray.Set(nbins + 1);
179
180    Double_t mymax = TMath::Max(min, max);
181    Double_t mymin = TMath::Min(min, max);
182
183    Int_t    k = 0;
184    Double_t binSize = (mymax - mymin) / ((Double_t)nbins);
185
186    fBinArray[0] = mymin;
187    for (k = 1; k <= nbins; k++) fBinArray[k] = fBinArray[k - 1] + binSize;
188 }
189
190 //_____________________________________________________________________________
191 void AliRsnValue::SetBins(Double_t min, Double_t max, Double_t step)
192 {
193 //
194 // Set binning for the axis in equally spaced bins
195 // where the bin size, minimum and maximum are given.
196 //
197
198    Double_t dblNbins = TMath::Abs(max - min) / step;
199    Int_t    intNbins = ((Int_t)dblNbins) + 1;
200
201    SetBins(intNbins, min, max);
202 }
203
204 //_____________________________________________________________________________
205 void AliRsnValue::SetBins(Int_t nbins, Double_t *array)
206 {
207 //
208 // Set binning for the axis in unequally spaced bins
209 // using the same way it is done in TAxis
210 //
211
212    if (!nbins) {
213       fBinArray.Set(0);
214       return;
215    }
216
217    fBinArray.Adopt(nbins, array);
218 }
219
220 //_____________________________________________________________________________
221 const char* AliRsnValue::GetValueTypeName() const
222 {
223 //
224 // This method returns a string to give a name to each possible
225 // computation value.
226 //
227
228    switch (fValueType) {
229       case kTrackP:               return "SingleTrackPtot";
230       case kTrackPt:              return "SingleTrackPt";
231       case kTrackEta:             return "SingleTrackEta";
232       case kTrackY:               return "SingleTrackRapidity";
233       case kTrackITSsignal:       return "SingleTrackITSsignal";
234       case kTrackTPCsignal:       return "SingleTrackTPCsignal";
235       case kTrackTOFsignal:       return "SingleTrackTOFsignal";
236       case kTrackLength:          return "SingleTrackLength";
237       case kPairP1:               return "PairPtotDaughter1";
238       case kPairP2:               return "PairPtotDaughter2";
239       case kPairP1t:              return "PairPtDaughter1";
240       case kPairP2t:              return "PairPtDaughter2";
241       case kPairP1z:              return "PairPzDaughter1";
242       case kPairP2z:              return "PairPzDaughter2";
243       case kPairInvMass:          return "PairInvMass";
244       case kPairInvMassMC:        return "PairInvMassMC";
245       case kPairInvMassRes:       return "PairInvMassResolution";
246       case kPairPt:               return "PairPt";
247       case kPairPz:               return "PairPz";
248       case kPairEta:              return "PairEta";
249       case kPairMt:               return "PairMt";
250       case kPairY:                return "PairY";
251       case kPairPhi:              return "PairPhi";
252       case kPairPhiMC:            return "PairPhiMC";
253       case kPairPtRatio:          return "PairPtRatio";
254       case kPairDipAngle:         return "PairDipAngle";
255       case kPairCosThetaStar:     return "PairCosThetaStar";
256       case kPairQInv:             return "PairQInv";
257       case kPairAngleToLeading:   return "PairAngleToLeading";
258       case kEventLeadingPt:       return "EventLeadingPt";
259       case kEventMult:            return "EventMult";
260       case kEventMultMC:          return "EventMultMC";
261       case kEventMultESDCuts:     return "EventMultESDCuts";
262       case kEventVz:              return "EventVz";
263       case kEventCentralityV0:    return "EventCentralityV0";
264       case kEventCentralityTrack: return "EventCentralityTrack";
265       case kEventCentralityCL1:   return "EventCentralityCL1";
266       default:                    return "Undefined";
267    }
268 }
269
270 //_____________________________________________________________________________
271 Bool_t AliRsnValue::Eval(TObject *object, Bool_t useMC)
272 {
273 //
274 // Evaluation of the required value.
275 // Checks that the passed object is of the right type
276 // and if this check is successful, computes the required value.
277 // The output of the function tells if computing was successful,
278 // and the values must be taken with GetValue().
279 //
280
281    // coherence check 
282    // (which also casts object to AliRsnTarget data members)
283    if (!TargetOK(object)) return kFALSE;
284    if (IsAllNull()) {
285       AliError("TargetOK passed but all pointers are NULL");
286       return kFALSE;
287    }
288
289    // cast the input to the allowed types
290    AliESDEvent *esdev  = fgCurrentEvent->GetRefESD();
291    AliESDtrack *esdt   = 0x0;
292    AliAODTrack *aodt   = 0x0;
293    AliAODPid   *pidObj = 0x0;
294    
295    // conditional initializations
296    if (fDaughter) {
297       esdt = fDaughter->GetRefESDtrack();
298       aodt = fDaughter->GetRefAODtrack();
299       if (aodt) pidObj = aodt->GetDetPid();
300    }
301
302    // common variables
303    TLorentzVector pRec;   // 4-momentum for single track or pair sum (reco)
304    TLorentzVector pSim;   // 4-momentum for single track or pair sum (MC)
305    TLorentzVector pRec0;  // 4-momentum of first daughter (reco)
306    TLorentzVector pSim0;  // 4-momentum of first daughter (MC)
307    TLorentzVector pRec1;  // 4-momentum of second daughter (reco)
308    TLorentzVector pSim1;  // 4-momentum of second daughter (MC)
309
310    // initialize the above 4-vectors according to the 
311    // expected target type (which has been checked above)
312    switch (fTargetType) {
313       case AliRsnTarget::kDaughter:
314          pRec = fDaughter->Psim();
315          pSim = fDaughter->Prec();
316          break;
317       case AliRsnTarget::kMother:
318          pRec  = fMother->Sum();
319          pSim  = fMother->SumMC();
320          pRec0 = fMother->GetDaughter(0)->Prec();
321          pRec1 = fMother->GetDaughter(1)->Prec();
322          pSim0 = fMother->GetDaughter(0)->Psim();
323          pSim1 = fMother->GetDaughter(1)->Psim();
324          break;
325       case AliRsnTarget::kEvent:
326          if (!fgCurrentEvent) {
327             AliError(Form("[%s] current event not initialized", GetName()));
328             return kFALSE;
329          }
330          break;
331       default:
332          AliError(Form("[%s] Wrong type", GetName()));
333          return kFALSE;
334    }
335
336    // cast the support object to the types which could be needed
337    AliRsnPairDef     *pairDef     = 0x0;
338    AliRsnDaughterDef *daughterDef = 0x0;
339    AliESDpid         *esdPID      = 0x0;
340    if (fSupportObject) {
341       if (fSupportObject->InheritsFrom(AliRsnPairDef    ::Class())) pairDef     = static_cast<AliRsnPairDef*>(fSupportObject);
342       if (fSupportObject->InheritsFrom(AliRsnDaughterDef::Class())) daughterDef = static_cast<AliRsnDaughterDef*>(fSupportObject);
343       if (fSupportObject->InheritsFrom(AliESDpid        ::Class())) esdPID      = static_cast<AliESDpid*>(fSupportObject);
344    }
345    
346    // compute value depending on types in the enumeration
347    // if the type does not match any available choice, or if
348    // the computation is not doable due to any problem
349    // (not initialized support object, wrong values, risk of floating point errors)
350    // the method returng kFALSE and sets the computed value to a very large number
351    switch (fValueType) {
352       case kTrackP:
353          // single track:
354          // total momentum 
355          fComputedValue = useMC ? pSim.Mag() : pRec.Mag();
356          return kTRUE;
357       case kTrackPt:
358          // single track:
359          // transverse momentum
360          fComputedValue = useMC ? pSim.Perp() : pRec.Perp();
361          return kTRUE;
362       case kTrackEta:
363          // single track:
364          // pseudo-rapidity
365          fComputedValue = useMC ? pSim.Eta() : pRec.Eta();
366          return kTRUE;
367       case kTrackY:
368          // single track:
369          // rapidity (requires an AliRsnDaughterDef support)
370          if (daughterDef) {
371             pRec.SetXYZM(pRec.X(), pRec.Y(), pRec.Z(), daughterDef->GetMass());
372             pSim.SetXYZM(pSim.X(), pSim.Y(), pSim.Z(), daughterDef->GetMass());
373             fComputedValue = useMC ? pSim.Rapidity() : pRec.Rapidity();
374             return kTRUE;
375          }
376          else {
377             AliError(Form("[%s] Required a correctly initialized AliRsnDaughterDef support object to compute this value", GetName()));
378             fComputedValue = fgkVeryBig;
379             return kFALSE;
380          } 
381       case kTrackITSsignal:
382          // single track:
383          // ITS signal (successful only for tracks)
384          if (esdt) {
385             fComputedValue = esdt->GetITSsignal();
386             return kTRUE;
387          }
388          else if (aodt && pidObj) {
389             fComputedValue = pidObj->GetITSsignal();
390             return kTRUE;
391          }
392          else {
393             AliError(Form("[%s] Detector signals can be computed only on reconstructed tracks", GetName()));
394             fComputedValue = fgkVeryBig;
395             return kFALSE;
396          }
397       case kTrackTPCsignal:
398          // single track:
399          // TPC signal (successful only for tracks)
400          if (esdt) {
401             fComputedValue = esdt->GetTPCsignal();
402             return kTRUE;
403          }
404          else if (aodt && pidObj) {
405             fComputedValue = pidObj->GetTPCsignal();
406             return kTRUE;
407          }
408          else {
409             AliError(Form("[%s] Detector signals can be computed only on reconstructed tracks", GetName()));
410             fComputedValue = fgkVeryBig;
411             return kFALSE;
412          }
413       case kTrackTOFsignal:
414          // single track:
415          // TOF signal (successful only for tracks, for ESD requires an AliESDpid support)
416          if (esdt) {
417             if (!esdPID) {
418                AliError(Form("[%s] Required a correctly initialized AliESDpid support object to compute this value with ESDs", GetName()));
419                fComputedValue = fgkVeryBig;
420                return kFALSE;
421             }
422             else {
423                fComputedValue = (esdt->GetTOFsignal() - esdPID->GetTOFResponse().GetStartTime(esdt->P()));
424                return kTRUE;
425             }
426          }
427          else if (aodt && pidObj) {
428             fComputedValue = pidObj->GetTOFsignal();
429             return kTRUE;
430          }
431          else {
432             AliError(Form("[%s] Detector signals can be computed only on reconstructed tracks", GetName()));
433             fComputedValue = fgkVeryBig;
434             return kFALSE;
435          }
436       case kTrackLength:
437          // single tracks:
438          // integrated length (computed only on ESDs)
439          if (esdt) {
440             fComputedValue = esdt->GetIntegratedLength();
441             return kTRUE;
442          }
443          else {
444             AliError(Form("[%s] Length information not available in AODs", GetName()));
445             fComputedValue = fgkVeryBig;
446             return kFALSE;
447          }
448       //---------------------------------------------------------------------------------------------------------------------
449       case kPairP1:
450          // pair:
451          // momentum of first daughter (which matches definition #1 in pairDef)
452          fComputedValue = useMC ? pSim0.Mag() : pRec0.Mag();
453          return kTRUE;
454       case kPairP2:
455          // pair:
456          // momentum of second daughter (which matches definition #2 in pairDef)
457          fComputedValue = useMC ? pSim1.Mag() : pRec1.Mag();
458          return kTRUE;
459       case kPairP1t:
460          // pair:
461          // transverse momentum of first daughter 
462          fComputedValue = useMC ? pSim0.Perp() : pRec0.Perp();
463          return kTRUE;
464       case kPairP2t:
465          // pair:
466          // transverse momentum of second daughter 
467          fComputedValue = useMC ? pSim1.Perp() : pRec1.Perp();
468          return kTRUE;
469       case kPairP1z:
470          // pair:
471          // longitudinal momentum of first daughter 
472          fComputedValue = useMC ? pSim0.Z() : pRec0.Z();
473          return kTRUE;
474       case kPairP2z:
475          // pair:
476          // longitudinal momentum of second daughter 
477          fComputedValue = useMC ? pSim1.Z() : pRec1.Z();
478          return kTRUE;
479       case kPairInvMass:
480          // pair:
481          // invariant mass
482          fComputedValue = useMC ? pSim.M() : pRec.M();
483          return kTRUE;
484       case kPairInvMassRes:
485          // pair:
486          // invariant mass resolution (requires MC)
487          if (TMath::Abs(pSim.M()) > 0.0) {
488             fComputedValue = (pSim.M() - pRec.M()) / pSim.M();
489             return kTRUE;
490          }
491          else {
492             AliError(Form("[%s] Caught a null MC mass", GetName()));
493             return kFALSE;
494          }
495       case kPairPt:
496          // pair:
497          // total transverse momentum
498          fComputedValue = useMC ? pSim.Perp() : pRec.Perp();
499          return kTRUE;
500       case kPairEta:
501          // pair:
502          // pseudo-rapidiry
503          fComputedValue = useMC ? pSim.Eta() : pRec.Eta();
504          return kTRUE;
505       case kPairMt:
506          // pair:
507          // transverse mass (requires an AliRsnPairDef to get mass hypothesis)
508          if (pairDef) {
509             pRec.SetXYZM(pRec.X(), pRec.Y(), pRec.Z(), pairDef->GetMotherMass());
510             pSim.SetXYZM(pSim.X(), pSim.Y(), pSim.Z(), pairDef->GetMotherMass());
511             fComputedValue = useMC ? pSim.Mt() : pRec.Mt();
512             return kTRUE;
513          }
514          else {
515             AliError(Form("[%s] Required a correctly initialized AliRsnPairDef support object to compute this value", GetName()));
516             fComputedValue = fgkVeryBig;
517             return kFALSE;
518          }
519       case kPairY:
520          // pair:
521          // rapidity (requires an AliRsnPairDef to get mass hypothesis)
522          if (pairDef) {
523             pRec.SetXYZM(pRec.X(), pRec.Y(), pRec.Z(), pairDef->GetMotherMass());
524             pSim.SetXYZM(pSim.X(), pSim.Y(), pSim.Z(), pairDef->GetMotherMass());
525             fComputedValue = useMC ? pSim.Rapidity() : pRec.Rapidity();
526             return kTRUE;
527          }
528          else {
529             AliError(Form("[%s] Required a correctly initialized AliRsnPairDef support object to compute this value", GetName()));
530             fComputedValue = fgkVeryBig;
531             return kFALSE;
532          }
533       case kPairPhi:
534          // pair:
535          // phi angle of total momentum
536          fComputedValue = useMC ? pSim.Phi() : pRec.Phi();
537          return kTRUE;
538       case kPairPtRatio:
539          // pair:
540          // ratio of relative sum and difference of daughter transverse momenta
541          if (useMC) {
542             fComputedValue  = TMath::Abs(pSim0.Perp() - pSim1.Perp());
543             fComputedValue /= TMath::Abs(pSim0.Perp() + pSim1.Perp());
544          } else {
545             fComputedValue  = TMath::Abs(pRec0.Perp() - pRec1.Perp());
546             fComputedValue /= TMath::Abs(pRec0.Perp() + pRec1.Perp());
547          }
548          return kTRUE;
549       case kPairDipAngle:
550          // pair:
551          // dip-angle in the transverse-Z plane
552          // (used to check conversion electrons)
553          if (useMC) {
554             fComputedValue  = pSim0.Perp() * pSim1.Perp() + pSim0.Z() * pSim1.Z();
555             fComputedValue /= pSim0.Mag() * pSim1.Mag();
556          } else {
557             fComputedValue  = pRec0.Perp() * pRec1.Perp() + pRec0.Z() * pRec1.Z();
558             fComputedValue /= pRec0.Mag() * pRec1.Mag();
559          }
560          fComputedValue = TMath::Abs(TMath::ACos(fComputedValue));
561          return kTRUE;
562       case kPairCosThetaStar:
563          // pair:
564          // cosine of theta star angle
565          // (= angle of first daughter to the total momentum, in resonance rest frame)
566          fComputedValue = fMother->CosThetaStar(useMC);
567          return kTRUE;
568       case kPairQInv:
569          // pair:
570          // Q-invariant
571          pSim0 -= pSim1;
572          pRec0 -= pRec1;
573          fComputedValue = useMC ? pSim0.M() : pRec0.M();
574          return kTRUE;
575       case kPairAngleToLeading:
576          // pair:
577          // angle w.r. to leading particle (if any)
578          {
579             int ID1 = (fMother->GetDaughter(0))->GetID();
580             int ID2 = (fMother->GetDaughter(1))->GetID();
581             Int_t leadingID = fgCurrentEvent->GetLeadingParticleID();
582             if (leadingID == ID1 || leadingID == ID2) {
583                fComputedValue = -99.;
584                return kFALSE;
585             }
586             AliRsnDaughter leadingPart = fgCurrentEvent->GetDaughter(leadingID);
587             AliVParticle  *ref = leadingPart.GetRef();
588             fComputedValue = ref->Phi() - fMother->Sum().Phi();
589             //return angle w.r.t. leading particle in the range -pi/2, 3/2pi
590             while (fComputedValue >= 1.5 * TMath::Pi()) fComputedValue -= 2 * TMath::Pi();
591             while (fComputedValue < -0.5 * TMath::Pi()) fComputedValue += 2 * TMath::Pi();
592          }
593          return kTRUE;
594       //---------------------------------------------------------------------------------------------------------------------
595       case kEventMult:
596          // event:
597          // multiplicity of tracks
598          fComputedValue = (Double_t)fgCurrentEvent->GetMultiplicity(0x0);
599          return kTRUE;
600       case kEventMultMC:
601          // event:
602          // multiplicity of MC tracks
603          fComputedValue = (Double_t)fgCurrentEvent->GetMultiplicityMC();
604       case kEventMultESDCuts:
605          // event:
606          // multiplicity of good quality tracks 
607          if (esdev) {
608             fComputedValue = AliESDtrackCuts::GetReferenceMultiplicity(esdev, kTRUE);
609             return kTRUE;
610          }
611          else {
612             AliError(Form("[%s] Can be computed only on ESDs", GetName()));
613             fComputedValue = fgkVeryBig;
614             return kFALSE;
615          }
616       case kEventLeadingPt: 
617          // event:
618          // transverse momentum of leading particle
619          {
620             int leadingID = fgCurrentEvent->GetLeadingParticleID(); //fEvent->SelectLeadingParticle(0);
621             if (leadingID >= 0) {
622                AliRsnDaughter leadingPart = fgCurrentEvent->GetDaughter(leadingID);
623                AliVParticle *ref = leadingPart.GetRef();
624                fComputedValue = ref->Pt();
625             } else fComputedValue = 0;
626          }
627          return kTRUE;
628       case kEventVz:
629          // event:
630          // Z position of primary vertex
631          fComputedValue = fgCurrentEvent->GetRef()->GetPrimaryVertex()->GetZ();
632          return kTRUE;
633       case kEventCentralityV0:
634          // event:
635          // centrality using V0 method
636          if (esdev) {
637             AliCentrality *centrality = esdev->GetCentrality();
638             fComputedValue = centrality->GetCentralityPercentile("V0M");
639             return kTRUE;
640          } else {
641             AliError(Form("[%s] Centrality computation is implemented for ESDs only up to now", GetName()));
642             return kFALSE;
643          }
644       case kEventCentralityTrack:
645          // event:
646          // centrality using tracks method
647          if (esdev) {
648             AliCentrality *centrality = esdev->GetCentrality();
649             fComputedValue = centrality->GetCentralityPercentile("TRK");
650             return kTRUE;
651          } else {
652             AliError(Form("[%s] Centrality computation is implemented for ESDs only up to now", GetName()));
653             return kFALSE;
654          }
655       case kEventCentralityCL1:
656          // event:
657          // centrality using CL1 method
658          if (esdev) {
659             AliCentrality *centrality = esdev->GetCentrality();
660             fComputedValue = centrality->GetCentralityPercentile("CL1");
661             return kTRUE;
662          } else {
663             AliError(Form("[%s] Centrality computation is implemented for ESDs only up to now", GetName()));
664             return kFALSE;
665          }
666       default:
667          AliError(Form("[%s] Invalid value type for this computation", GetName()));
668          return kFALSE;
669    }
670 }
671
672 //_____________________________________________________________________________
673 void AliRsnValue::Print(Option_t * /*option */) const
674 {
675 //
676 // Print informations about this object
677 //
678
679    AliInfo("=== VALUE INFO =================================================");
680    AliInfo(Form(" Name                  : %s", GetName()));
681    AliInfo(Form(" Type                  : %s", GetValueTypeName()));
682    AliInfo(Form(" Current computed value: %f", fComputedValue));
683    Int_t i;
684    for (i = 0; i < fBinArray.GetSize(); i++) {
685       AliInfo(Form(" Bin limit #%d         = %f", i, fBinArray[i]));
686    }
687    AliInfo(Form(" Support object        : %s", (fSupportObject ? fSupportObject->ClassName() : " NO SUPPORT")));
688    AliInfo("=== END VALUE INFO =============================================");
689 }
690
691 //_____________________________________________________________________________
692 RSNTARGET AliRsnValue::TargetType(EValueType type)
693 {
694 //
695 // This method assigns the target to be expected by this object
696 // in the computation, depending on its type chosen in the enum.
697 //
698
699    if (type < kTrackValues)
700       return AliRsnTarget::kDaughter;
701    else if (type < kPairValues)
702       return AliRsnTarget::kMother;
703    else
704       return AliRsnTarget::kEvent;
705 }