]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnValueDaughter.cxx
Fixed all fixable coding conventions violations
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnValueDaughter.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, Daughter, 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 #include "AliPIDResponse.h"
49
50 #include "AliRsnEvent.h"
51 #include "AliRsnDaughter.h"
52 #include "AliRsnMother.h"
53 #include "AliRsnDaughterDef.h"
54 #include "AliRsnDaughterDef.h"
55
56 #include "AliRsnValueDaughter.h"
57
58 ClassImp(AliRsnValueDaughter)
59
60 //_____________________________________________________________________________
61 AliRsnValueDaughter::AliRsnValueDaughter(const char *name, EType type) :
62    AliRsnValue(name, AliRsnTarget::kDaughter),
63    fType(type)
64 {
65 //
66 // Constructor
67 //
68 }
69
70 //_____________________________________________________________________________
71 AliRsnValueDaughter::AliRsnValueDaughter(const AliRsnValueDaughter& copy) :
72    AliRsnValue(copy),
73    fType(copy.fType)
74 {
75 //
76 // Copy constructor
77 //
78 }
79
80 //_____________________________________________________________________________
81 AliRsnValueDaughter& AliRsnValueDaughter::operator=(const AliRsnValueDaughter& copy)
82 {
83 //
84 // Assignment operator.
85 // Works like copy constructor.
86 //
87
88    AliRsnValue::operator=(copy);
89    fType = copy.fType;
90
91    return (*this);
92 }
93
94 //_____________________________________________________________________________
95 const char* AliRsnValueDaughter::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 kP:           return "SingleTrackPtot";
104       case kPt:          return "SingleTrackPt";
105       case kPtpc:        return "SingleTrackPtpc";
106       case kEta:         return "SingleTrackEta";
107       case kITSsignal:   return "SingleTrackITSsignal";
108       case kTPCsignal:   return "SingleTrackTPCsignal";
109       case kTOFsignal:   return "SingleTrackTOFsignal";
110       case kTOFnsigmaPi: return "SingleTrackTOFnsigmaPion";
111       case kTOFnsigmaK:  return "SingleTrackTOFnsigmaKaon";
112       case kTOFnsigmaP:  return "SingleTrackTOFnsigmaProton";
113       default:           return "Undefined";
114    }
115 }
116
117 //_____________________________________________________________________________
118 Bool_t AliRsnValueDaughter::Eval(TObject *object)
119 {
120 //
121 // Evaluation of the required value.
122 // Checks that the passed object is of the right type
123 // and if this check is successful, computes the required value.
124 // The output of the function tells if computing was successful,
125 // and the values must be taken with GetValue().
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    AliVParticle   *ref   = fDaughter->GetRef();
135    AliVParticle   *refMC = fDaughter->GetRefMC();
136    AliVTrack      *track = fDaughter->Ref2Vtrack();
137    if (fUseMCInfo && !refMC) {
138       AliError("No MC");
139       return kFALSE;
140    }
141    if (!fUseMCInfo && !ref) {
142       AliError("No DATA");
143       return kFALSE;
144    }
145    
146    // compute value depending on types in the enumeration
147    // if the type does not match any available choice, or if
148    // the computation is not doable due to any problem
149    // (not initialized support object, wrong values, risk of floating point errors)
150    // the method returng kFALSE and sets the computed value to a meaningless number
151    switch (fType) {
152       case kP:
153          fComputedValue = (fUseMCInfo ? refMC->P() : ref->P());
154          return kTRUE;
155       case kPt:
156          fComputedValue = (fUseMCInfo ? refMC->Pt() : ref->Pt());
157          return kTRUE;
158       case kEta:
159          fComputedValue = (fUseMCInfo ? refMC->Eta() : ref->Eta());
160          return kTRUE;
161       case kPtpc:
162          if (track) {
163             fComputedValue = track->GetTPCmomentum();
164             return kTRUE;
165          } else {
166             AliWarning("Cannot get TPC momentum for non-track object");
167             fComputedValue = 0.0;
168             return kFALSE;
169          }
170       case kITSsignal:
171          if (track) {
172             fComputedValue = track->GetITSsignal();
173             return kTRUE;
174          } else {
175             AliWarning("Cannot get ITS signal for non-track object");
176             fComputedValue = 0.0;
177             return kFALSE;
178          }
179       case kTPCsignal:
180          if (track) {
181             fComputedValue = track->GetTPCsignal();
182             return kTRUE;
183          } else {
184             AliWarning("Cannot get TPC signal for non-track object");
185             fComputedValue = 0.0;
186             return kFALSE;
187          }
188       case kTOFsignal:
189          if (track) {
190             fComputedValue = track->GetTOFsignal();
191             return kTRUE;
192          } else {
193             AliWarning("Cannot get TOF signal for non-track object");
194             fComputedValue = 0.0;
195             return kFALSE;
196          }
197       case kTOFnsigmaPi:
198          if (track) {
199             AliPIDResponse *pid = fEvent->GetPIDResponse();
200             fComputedValue = pid->NumberOfSigmasTOF(track, AliPID::kPion);
201             return kTRUE;
202          } else {
203             AliWarning("Cannot get TOF signal for non-track object");
204             fComputedValue = 0.0;
205             return kFALSE;
206          }
207       case kTOFnsigmaK:
208          if (track) {
209             AliPIDResponse *pid = fEvent->GetPIDResponse();
210             fComputedValue = pid->NumberOfSigmasTOF(track, AliPID::kKaon);
211             return kTRUE;
212          } else {
213             AliWarning("Cannot get TOF signal for non-track object");
214             fComputedValue = 0.0;
215             return kFALSE;
216          }
217       case kTOFnsigmaP:
218          if (track) {
219             AliPIDResponse *pid = fEvent->GetPIDResponse();
220             fComputedValue = pid->NumberOfSigmasTOF(track, AliPID::kProton);
221             return kTRUE;
222          } else {
223             AliWarning("Cannot get TOF signal for non-track object");
224             fComputedValue = 0.0;
225             return kFALSE;
226          }
227       default:
228          AliError(Form("[%s] Invalid value type for this computation", GetName()));
229          return kFALSE;
230    }
231 }