]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnValuePair.cxx
Updated macros for Sigma* analysis (M.Venaruzzo)
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnValuePair.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 "AliVVertex.h"
42 #include "AliMultiplicity.h"
43 #include "AliESDtrackCuts.h"
44 #include "AliESDpid.h"
45 #include "AliAODPid.h"
46 #include "AliCentrality.h"
47 #include "AliESDUtils.h"
48
49 #include "AliRsnEvent.h"
50 #include "AliRsnDaughter.h"
51 #include "AliRsnMother.h"
52 #include "AliRsnPairDef.h"
53 #include "AliRsnDaughterDef.h"
54
55 #include "AliRsnValuePair.h"
56
57 ClassImp(AliRsnValuePair)
58
59 //_____________________________________________________________________________
60 AliRsnValuePair::AliRsnValuePair(const char *name, EType type) :
61    AliRsnValue(name, AliRsnTarget::kMother),
62    fType(type)
63 {
64 //
65 // Constructor
66 //
67 }
68
69 //_____________________________________________________________________________
70 AliRsnValuePair::AliRsnValuePair(const AliRsnValuePair &copy) :
71    AliRsnValue(copy),
72    fType(copy.fType)
73 {
74 //
75 // Copy constructor
76 //
77 }
78
79 //_____________________________________________________________________________
80 AliRsnValuePair &AliRsnValuePair::operator=(const AliRsnValuePair &copy)
81 {
82 //
83 // Assignment operator.
84 // Works like copy constructor.
85 //
86
87    AliRsnValue::operator=(copy);
88    if (this == &copy)
89       return *this;
90    fType = copy.fType;
91
92    return (*this);
93 }
94
95 //_____________________________________________________________________________
96 const char *AliRsnValuePair::GetTypeName() const
97 {
98 //
99 // This method returns a string to give a name to each possible
100 // computation value.
101 //
102
103    switch (fType) {
104       case kPt:           return "PairPt";
105       case kPz:           return "PairPz";
106       case kInvMass:      return "PairInvMass";
107       case kInvMassRes:   return "PairInvMassResolution";
108       case kEta:          return "PairEta";
109       case kMt:           return "PairMt";
110       case kY:            return "PairY";
111       case kPtRatio:      return "PairPtRatio";
112       case kDipAngle:     return "PairDipAngle";
113       case kCosThetaStar: return "PairCosThetaStar";
114       case kAngleLeading: return "PairAngleToLeading";
115       default:            return "Undefined";
116    }
117 }
118
119 //_____________________________________________________________________________
120 Bool_t AliRsnValuePair::Eval(TObject *object)
121 {
122 //
123 // Evaluation of the required value.
124 // In this implementation, fills the member 4-vectors with data
125 // coming from the object passed as argument, and then returns the value
126 //
127
128    // coherence check, which also casts object
129    // to AliRsnTarget data members and returns kFALSE
130    // in case the object is NULL
131    if (!TargetOK(object)) return kFALSE;
132
133    // set a reference to the mother momentum
134    TLorentzVector &sum   = fMother->Sum(fUseMCInfo);
135    TLorentzVector &ref   = fMother->Ref(fUseMCInfo);
136    TLorentzVector &p1    = fMother->GetDaughter(0)->P(fUseMCInfo);
137    TLorentzVector &p2    = fMother->GetDaughter(1)->P(fUseMCInfo);
138
139    // utility variables
140    Bool_t success;
141
142    // compute value depending on types in the enumeration
143    // if the type does not match any available choice, or if
144    // the computation is not doable due to any problem
145    // (not initialized support object, wrong values, risk of floating point errors)
146    // the method returng kFALSE and sets the computed value to a meaningless number
147    switch (fType) {
148       case kPt:
149          fComputedValue = sum.Perp();
150          return kTRUE;
151       case kInvMass:
152          fComputedValue = sum.M();
153          return kTRUE;
154       case kEta:
155          fComputedValue = sum.Eta();
156          return kTRUE;
157       case kInvMassRes:
158          fComputedValue  = fMother->Sum(kFALSE).M() - fMother->Sum(kTRUE).M();
159          fComputedValue /= fMother->Sum(kTRUE).M();
160          return kTRUE;
161       case kMt:
162          fComputedValue = ref.Mt();
163          return kTRUE;
164       case kY:
165          fComputedValue = ref.Rapidity();
166          return kTRUE;
167       case kPtRatio:
168          fComputedValue  = TMath::Abs(p1.Perp() - p2.Perp());
169          fComputedValue /= TMath::Abs(p1.Perp() + p2.Perp());
170          return kTRUE;
171       case kDipAngle:
172          fComputedValue  = p1.Perp() * p2.Perp() + p1.Z() * p2.Z();
173          fComputedValue /= p1.Mag() * p2.Mag();
174          return kTRUE;
175       case kCosThetaStar:
176          fComputedValue = fMother->CosThetaStar();
177          return kTRUE;
178       case kAngleLeading:
179          fComputedValue = fMother->AngleToLeading(success);
180          return success;
181       default:
182          AliError(Form("[%s] Invalid value type for this computation", GetName()));
183          return kFALSE;
184    }
185 }