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