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