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