]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnValueEvent.cxx
Migration of PWG2/RESONANCES -> PWGLF/RESONANCES
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnValueEvent.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 "AliVVertex.h"
41 #include "AliMultiplicity.h"
42 #include "AliESDtrackCuts.h"
43 #include "AliESDpid.h"
44 #include "AliAODPid.h"
45 #include "AliCentrality.h"
46 #include "AliESDUtils.h"
47
48 #include "AliRsnEvent.h"
49 #include "AliRsnDaughter.h"
50 #include "AliRsnMother.h"
51 #include "AliRsnPairDef.h"
52 #include "AliRsnDaughterDef.h"
53
54 #include "AliRsnValueEvent.h"
55
56 ClassImp(AliRsnValueEvent)
57
58 //_____________________________________________________________________________
59 AliRsnValueEvent::AliRsnValueEvent(const char *name, EType type) :
60    AliRsnValue(name, AliRsnTarget::kEvent),
61    fType(type)
62 {
63 //
64 // Constructor
65 //
66 }
67
68 //_____________________________________________________________________________
69 AliRsnValueEvent::AliRsnValueEvent(const AliRsnValueEvent &copy) :
70    AliRsnValue(copy),
71    fType(copy.fType)
72 {
73 //
74 // Copy constructor
75 //
76 }
77
78 //_____________________________________________________________________________
79 AliRsnValueEvent &AliRsnValueEvent::operator=(const AliRsnValueEvent &copy)
80 {
81 //
82 // Assignment operator.
83 // Works like copy constructor.
84 //
85
86    AliRsnValue::operator=(copy);
87    if (this == &copy)
88       return *this;
89    fType = copy.fType;
90
91    return (*this);
92 }
93
94 //_____________________________________________________________________________
95 const char *AliRsnValueEvent::GetTypeName() const
96 {
97 //
98 // This method returns a string to give a name to each possible
99 // computation value.
100 //
101
102    switch (fType) {
103       case kLeadingPt:       return "EventLeadingPt";
104       case kMult:            return "EventMult";
105       case kMultMC:          return "EventMultMC";
106       case kMultESDCuts:     return "EventMultESDCuts";
107       case kMultSPD:         return "EventMultSPD";
108       case kVz:              return "EventVz";
109       case kCentralityV0:    return "EventCentralityV0";
110       case kCentralityTrack: return "EventCentralityTrack";
111       case kCentralityCL1:   return "EventCentralityCL1";
112       default:               return "Undefined";
113    }
114 }
115
116 //_____________________________________________________________________________
117 Bool_t AliRsnValueEvent::Eval(TObject *object)
118 {
119 //
120 // Evaluation of the required value.
121 // In this implementation, fills the member 4-vectors with data
122 // coming from the object passed as argument, and then returns the value
123 //
124
125    // coherence check, which also casts object
126    // to AliRsnTarget data members and returns kFALSE
127    // in case the object is NULL
128    if (!TargetOK(object)) return kFALSE;
129    if (!fEvent->GetRef()) {
130       AliWarning("NULL ref");
131       return kFALSE;
132    }
133
134    // declare support variables
135    AliCentrality *centrality = fEvent->GetRef()->GetCentrality();
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    switch (fType) {
143       case kMult:
144          fComputedValue = (Double_t)fEvent->GetRef()->GetNumberOfTracks();
145          return (fComputedValue >= 0);
146       case kMultMC:
147          fComputedValue = -999.0;
148          if (fEvent->GetRefMC()) {
149             if (fEvent->IsESD())
150                fComputedValue = (Double_t)fEvent->GetRefMC()->GetNumberOfTracks();
151             else {
152                AliAODEvent *aod = (AliAODEvent *)fEvent->GetRefMC();
153                TClonesArray *mcArray = (TClonesArray *)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
154                if (mcArray) fComputedValue = (Double_t)mcArray->GetEntries();
155             }
156          }
157          return (fComputedValue >= 0);
158       case kMultESDCuts:
159          fComputedValue = -999.0;
160          if (fEvent->IsESD()) {
161             fComputedValue = AliESDtrackCuts::GetReferenceMultiplicity(fEvent->GetRefESD(), kTRUE);
162          } else {
163             AliWarning("Cannot compute ESD cuts multiplicity in AOD");
164             return kFALSE;
165          }
166          return (fComputedValue >= 0);
167       case kMultSPD:
168          fComputedValue = -999.0;
169          if (fEvent->IsESD()) {
170             const AliMultiplicity *mult = fEvent->GetRefESD()->GetMultiplicity();
171             Float_t nClusters[6] = {0.0,0.0,0.0,0.0,0.0,0.0};
172             for(Int_t ilay = 0; ilay < 6; ilay++) nClusters[ilay] = (Float_t)mult->GetNumberOfITSClusters(ilay);
173             fComputedValue = AliESDUtils::GetCorrSPD2(nClusters[1], fEvent->GetRef()->GetPrimaryVertex()->GetZ());
174          } else {
175             AliWarning("Cannot compute SPD multiplicity with AOD");
176             return kFALSE;
177          }
178          return (fComputedValue >= 0);
179       case kLeadingPt:
180          if (fEvent->GetLeadingIndex() >= 0) {
181             AliRsnDaughter leadingPart;
182             fEvent->SetLeadingParticle(leadingPart);
183             fComputedValue = leadingPart.GetRef()->Pt();
184             return kTRUE;
185          } else {
186             AliError("Not found good leading particle");
187             return kFALSE;
188          }
189       case kVz:
190          fComputedValue = fEvent->GetRef()->GetPrimaryVertex()->GetZ();
191          return kTRUE;
192       case kCentralityV0:
193          if (centrality) {
194             fComputedValue = centrality->GetCentralityPercentile("V0M");
195             return kTRUE;
196          } else {
197             AliError("Centrality undefined");
198             return kFALSE;
199          }
200       case kCentralityTrack:
201          if (centrality) {
202             fComputedValue = centrality->GetCentralityPercentile("TRK");
203             return kTRUE;
204          } else {
205             AliError("Centrality undefined");
206             return kFALSE;
207          }
208       case kCentralityCL1:
209          if (centrality) {
210             fComputedValue = centrality->GetCentralityPercentile("CL1");
211             return kTRUE;
212          } else {
213             AliError("Centrality undefined");
214             return kFALSE;
215          }
216       default:
217          AliError(Form("[%s] Invalid value type for this computation", GetName()));
218          return kFALSE;
219    }
220 }