]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnMiniValue.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnMiniValue.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
42 #include "AliLog.h"
43
44 #include "AliRsnMiniPair.h"
45 #include "AliRsnMiniEvent.h"
46 #include "AliRsnMiniParticle.h"
47
48 #include "AliRsnMiniValue.h"
49
50 ClassImp(AliRsnMiniValue)
51
52 //_____________________________________________________________________________
53 AliRsnMiniValue::AliRsnMiniValue(EType type, Bool_t useMC) :
54    TNamed(ValueName(type, useMC), ""),
55    fType(type),
56    fUseMCInfo(useMC)
57 {
58 //
59 // Constructor
60 //
61 }
62
63 //_____________________________________________________________________________
64 AliRsnMiniValue::AliRsnMiniValue(const AliRsnMiniValue &copy) :
65    TNamed(copy),
66    fType(copy.fType),
67    fUseMCInfo(copy.fUseMCInfo)
68 {
69 //
70 // Copy constructor
71 //
72 }
73
74 //_____________________________________________________________________________
75 AliRsnMiniValue &AliRsnMiniValue::operator=(const AliRsnMiniValue &copy)
76 {
77 //
78 // Assignment operator.
79 // Works like copy constructor.
80 //
81    TNamed::operator=(copy);
82    if (this == &copy)
83       return *this;
84    fType = copy.fType;
85    fUseMCInfo = copy.fUseMCInfo;
86
87    return (*this);
88 }
89
90 //_____________________________________________________________________________
91 const char *AliRsnMiniValue::TypeName(EType type)
92 {
93 //
94 // This method returns a string to give a name to each possible
95 // computation value.
96 //
97
98    switch (type) {
99       case kVz:           return "EventVz";
100       case kMult:         return "EventMult";
101       case kPlaneAngle:   return "EventPlane";
102       case kLeadingPt:    return "EventLeadingPt";
103       case kPt:           return "Pt";
104       case kPz:           return "Pz";
105       case kInvMass:      return "InvMass";
106       case kInvMassRes:   return "InvMassResolution";
107       case kInvMassDiff:  return "InvMassDifference";
108       case kEta:          return "Eta";
109       case kMt:           return "Mt";
110       case kY:            return "Y";
111       case kPtRatio:      return "PtRatio";
112       case kDipAngle:     return "DipAngle";
113       case kCosThetaStar: return "CosThetaStar";
114       case kAngleLeading: return "AngleToLeading";
115       case kFirstDaughterPt: return "FirstDaughterPt";
116       case kSecondDaughterPt: return "SecondDaughterPt";
117       case kFirstDaughterP: return "FirstDaughterP";
118       case kSecondDaughterP: return "SecondDaughterP";
119       default:            return "Undefined";
120    }
121 }
122
123 //_____________________________________________________________________________
124 Float_t AliRsnMiniValue::Eval(AliRsnMiniPair *pair, AliRsnMiniEvent *event)
125 {
126 //
127 // Evaluation of the required value.
128 // In this implementation, fills the member 4-vectors with data
129 // coming from the object passed as argument, and then returns the value
130 //
131
132    if (!pair && fType > kEventCuts) {
133       AliError("Null pair passed!");
134       return 1E20;
135    }
136
137    // compute value depending on types in the enumeration
138    // if the type does not match any available choice, or if
139    // the computation is not doable due to any problem
140    // (not initialized support object, wrong values, risk of floating point errors)
141    // the method returng kFALSE and sets the computed value to a meaningless number
142    Double_t p3[3]= {0.,0.,0.};
143    AliRsnMiniParticle *l;
144    TLorentzVector v;
145    switch (fType) {
146          // ---- event values -------------------------------------------------------------------------
147       case kVz:
148          return event->Vz();
149       case kMult:
150          return event->Mult();
151       case kPlaneAngle:
152          return event->Angle();
153       case kLeadingPt:
154          l = event->LeadingParticle();
155          if (l) {
156             l->Set4Vector(v,-1.0,fUseMCInfo);
157             return v.Pt();
158          }
159          return 0.0;
160       case kPt:
161          return pair->Pt(fUseMCInfo);
162       case kInvMass:
163          return pair->InvMass(fUseMCInfo);
164       case kEta:
165          return pair->Eta(fUseMCInfo);
166       case kInvMassRes:
167          return pair->InvMassRes();
168       case kInvMassDiff:
169          return pair->InvMassDiff();
170       case kMt:
171          return pair->Mt(fUseMCInfo);
172       case kY:
173          return pair->Y(fUseMCInfo);
174       case kPtRatio:
175          return pair->PtRatio(fUseMCInfo);
176       case kDipAngle:
177          return pair->DipAngle(fUseMCInfo);
178       case kCosThetaStar:
179          return pair->CosThetaStar(fUseMCInfo);
180       case kAngleLeading:
181          l = event->LeadingParticle();
182          if (l) {
183             l->Set4Vector(v,-1.0,fUseMCInfo);
184             Double_t angle = v.Phi() - pair->Sum(fUseMCInfo).Phi();
185
186             //return angle w.r.t. leading particle in the range -pi/2, 3/2pi
187             while (angle >= 1.5 * TMath::Pi()) angle -= 2 * TMath::Pi();
188             while (angle < -0.5 * TMath::Pi()) angle += 2 * TMath::Pi();
189             return angle;
190          }
191 //         AliWarning("This method is not yet implemented");
192          return 1E20;
193       case kFirstDaughterPt:
194          return pair->DaughterPt(0,fUseMCInfo);
195       case kSecondDaughterPt:
196          return pair->DaughterPt(1,fUseMCInfo);
197       case kFirstDaughterP:
198          pair->DaughterPxPyPz(0,fUseMCInfo, p3);
199          return TMath::Sqrt(p3[0]*p3[0]+p3[1]*p3[1]+p3[2]*p3[2]);
200       case kSecondDaughterP:
201          pair->DaughterPxPyPz(1,fUseMCInfo, p3);
202          return TMath::Sqrt(p3[0]*p3[0]+p3[1]*p3[1]+p3[2]*p3[2]);
203       default:
204          AliError("Invalid value type");
205          return 1E20;
206    }
207 }