]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnValueEvent.cxx
introduce the option to analyze with PROOF, plus cosmetics
[u/mrichter/AliRoot.git] / PWG2 / 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    fType = copy.fType;
88
89    return (*this);
90 }
91
92 //_____________________________________________________________________________
93 const char* AliRsnValueEvent::GetTypeName() const
94 {
95 //
96 // This method returns a string to give a name to each possible
97 // computation value.
98 //
99
100    switch (fType) {
101       case kLeadingPt:       return "EventLeadingPt";
102       case kMult:            return "EventMult";
103       case kMultMC:          return "EventMultMC";
104       case kMultESDCuts:     return "EventMultESDCuts";
105       case kMultSPD:         return "EventMultSPD";
106       case kVz:              return "EventVz";
107       case kCentralityV0:    return "EventCentralityV0";
108       case kCentralityTrack: return "EventCentralityTrack";
109       case kCentralityCL1:   return "EventCentralityCL1";
110       default:               return "Undefined";
111    }
112 }
113
114 //_____________________________________________________________________________
115 Bool_t AliRsnValueEvent::Eval(TObject *object)
116 {
117 //
118 // Evaluation of the required value.
119 // In this implementation, fills the member 4-vectors with data
120 // coming from the object passed as argument, and then returns the value
121 //
122    
123    // coherence check, which also casts object 
124    // to AliRsnTarget data members and returns kFALSE
125    // in case the object is NULL
126    if (!TargetOK(object)) return kFALSE;
127    if (!fEvent->GetRef()) {
128       AliWarning("NULL ref");
129       return kFALSE;
130    }
131    
132    // declare support variables
133    AliCentrality *centrality = fEvent->GetRef()->GetCentrality();
134    
135    // compute value depending on types in the enumeration
136    // if the type does not match any available choice, or if
137    // the computation is not doable due to any problem
138    // (not initialized support object, wrong values, risk of floating point errors)
139    // the method returng kFALSE and sets the computed value to a meaningless number
140    switch (fType) {
141       case kMult:
142          fComputedValue = (Double_t)fEvent->GetRef()->GetNumberOfTracks();
143          return (fComputedValue >= 0);
144       case kMultMC:
145          fComputedValue = -999.0;
146          if (fEvent->GetRefMC()) {
147             if (fEvent->IsESD()) 
148                fComputedValue = (Double_t)fEvent->GetRefMC()->GetNumberOfTracks();
149             else {
150                AliAODEvent *aod = (AliAODEvent*)fEvent->GetRefMC();
151                TClonesArray *mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
152                if (mcArray) fComputedValue = (Double_t)mcArray->GetEntries();
153             }
154          }
155          return (fComputedValue >= 0);
156       case kMultESDCuts:
157          fComputedValue = -999.0;
158          if (fEvent->IsESD()) {
159             fComputedValue = AliESDtrackCuts::GetReferenceMultiplicity(fEvent->GetRefESD(), kTRUE);
160          } else {
161             AliWarning("Cannot compute ESD cuts multiplicity in AOD");
162             return kFALSE;
163          }
164          return (fComputedValue >= 0);
165       case kMultSPD:
166          fComputedValue = -999.0;
167          if (fEvent->IsESD()) {
168             const AliMultiplicity *mult = fEvent->GetRefESD()->GetMultiplicity();
169             Float_t nClusters[6] = {0.0,0.0,0.0,0.0,0.0,0.0};
170             for(Int_t ilay = 0; ilay < 6; ilay++) nClusters[ilay] = (Float_t)mult->GetNumberOfITSClusters(ilay);
171             fComputedValue = AliESDUtils::GetCorrSPD2(nClusters[1], fEvent->GetRef()->GetPrimaryVertex()->GetZ());
172          } else {
173             AliWarning("Cannot compute SPD multiplicity with AOD");
174             return kFALSE;
175          }
176          return (fComputedValue >= 0);
177       case kLeadingPt: 
178          if (fEvent->GetLeadingIndex() >= 0) {
179             AliRsnDaughter leadingPart;
180             fEvent->SetLeadingParticle(leadingPart);
181             fComputedValue = leadingPart.GetRef()->Pt();
182             return kTRUE;
183          } else {
184             AliError("Not found good leading particle");
185             return kFALSE;
186          }
187       case kVz:
188          fComputedValue = fEvent->GetRef()->GetPrimaryVertex()->GetZ();
189          return kTRUE;
190       case kCentralityV0:
191          if (centrality) {
192             fComputedValue = centrality->GetCentralityPercentile("V0M");
193             return kTRUE;
194          } else {
195             AliError("Centrality undefined");
196             return kFALSE;
197          }
198       case kCentralityTrack:
199          if (centrality) {
200             fComputedValue = centrality->GetCentralityPercentile("TRK");
201             return kTRUE;
202          } else {
203             AliError("Centrality undefined");
204             return kFALSE;
205          }
206       case kCentralityCL1:
207          if (centrality) {
208             fComputedValue = centrality->GetCentralityPercentile("CL1");
209             return kTRUE;
210          } else {
211             AliError("Centrality undefined");
212             return kFALSE;
213          }
214       default:
215          AliError(Form("[%s] Invalid value type for this computation", GetName()));
216          return kFALSE;
217    }
218 }