]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnMiniMonitor.cxx
fixing psi in MC from the header
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnMiniMonitor.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 "TH1.h"
43 #include "TH2.h"
44 #include "TH3.h"
45
46 #include "AliLog.h"
47 #include "AliPID.h"
48 #include "AliPIDResponse.h"
49
50 #include "AliRsnDaughter.h"
51 #include "AliRsnEvent.h"
52
53 #include "AliRsnMiniMonitor.h"
54
55 ClassImp(AliRsnMiniMonitor)
56
57 //__________________________________________________________________________________________________
58 AliRsnMiniMonitor::AliRsnMiniMonitor() :
59    TNamed(),
60    fType(kTypes),
61    fCutID(-1),
62    fCharge(0),
63    fListID(-1),
64    fList(0x0)
65 {
66 //
67 // Dummy constructor
68 //
69 }
70
71 //__________________________________________________________________________________________________
72 AliRsnMiniMonitor::AliRsnMiniMonitor(const char *name, EType type, Int_t cutID) :
73    TNamed(name, ""),
74    fType(type),
75    fCutID(cutID),
76    fCharge(0),
77    fListID(-1),
78    fList(0x0)
79 {
80 //
81 // Default constructor
82 //
83 }
84
85 //__________________________________________________________________________________________________
86 AliRsnMiniMonitor::AliRsnMiniMonitor(const AliRsnMiniMonitor &copy) :
87    TNamed(copy),
88    fType(copy.fType),
89    fCutID(copy.fCutID),
90    fCharge(copy.fCharge),
91    fListID(copy.fListID),
92    fList(copy.fList)
93 {
94 //
95 // Copy constructor
96 //
97 }
98
99 //__________________________________________________________________________________________________
100 AliRsnMiniMonitor &AliRsnMiniMonitor::operator=(const AliRsnMiniMonitor &copy)
101 {
102 //
103 // Assignment operator
104 //
105
106    TNamed::operator=(copy);
107    if (this == &copy)
108       return *this;
109    fType = copy.fType;
110    fCutID = copy.fCutID;
111    fCharge = copy.fCharge;
112    fListID = copy.fListID;
113    fList = copy.fList;
114
115    return (*this);
116 }
117
118 //__________________________________________________________________________________________________
119 Bool_t AliRsnMiniMonitor::Init(const char *name, TList *list)
120 {
121 //
122 // Initialize this output histogram and put into the passed list
123 //
124
125    // name
126    TString sname(name);
127    sname += '_';
128    sname += GetName();
129
130    // check list
131    fList = list;
132    if (!list) {
133       AliError("No list!");
134       return kFALSE;
135    }
136
137    // reset histogram
138    TH1 *histogram = 0x0;
139
140    switch (fType) {
141       case kTrackPt:
142          sname += "_TrackPt_cut";
143          sname += fCutID;
144          histogram = new TH1F(sname.Data(), "", 100, 0.0, 10.0);
145          break;
146       case kdEdxTPCvsP:
147          sname += "_TPCsignal";
148          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, 0.0, 1000.0);
149          break;
150       case ktimeTOFvsPPion:
151          sname += "_TOFsignalPi";
152          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, 0.0, 200000.0);
153       case ktimeTOFvsPKaon:
154          sname += "_TOFsignalK";
155          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, 0.0, 200000.0);
156          break;
157       case ktimeTOFvsPProton:
158          sname += "_TOFsignalP";
159          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, 0.0, 200000.0);
160          break;
161       default:
162          AliError("Wrong enum type");
163          return kFALSE;
164    }
165
166    // add to list
167    if (histogram && fList) {
168       histogram->Sumw2();
169       fList->Add(histogram);
170       fListID = fList->IndexOf(histogram);
171       AliInfo(Form("Histogram '%s' added to list in slot #%d", histogram->GetName(), fListID));
172       return kTRUE;
173    }
174
175    return kFALSE;
176 }
177
178 //_____________________________________________________________________________
179 Bool_t AliRsnMiniMonitor::Fill(AliRsnDaughter *track, AliRsnEvent *event)
180 {
181 //
182 // Fill the histogram
183 //
184
185    // retrieve object from list
186    if (!fList) {
187       AliError("List pointer is NULL");
188       return kFALSE;
189    }
190    TObject *obj = fList->At(fListID);
191    if (!obj) {
192       AliError("List object is NULL");
193       return kFALSE;
194    }
195
196    Double_t valueX, valueY;
197    AliVTrack *vtrack = track->Ref2Vtrack();
198
199    AliPIDResponse *pid = event->GetPIDResponse();
200
201    switch (fType) {
202       case kTrackPt:
203          if (!vtrack) {
204             AliWarning("Required vtrack for this value");
205             return kFALSE;
206          }
207          if (fCharge == '+' && vtrack->Charge() <= 0) return kFALSE;
208          if (fCharge == '-' && vtrack->Charge() >= 0) return kFALSE;
209          if (fCharge == '0' && vtrack->Charge() != 0) return kFALSE;
210          valueX = vtrack->Pt();
211          ((TH1F *)obj)->Fill(valueX);
212          return kTRUE;
213       case kdEdxTPCvsP:
214          if (!vtrack) {
215             AliWarning("Required vtrack for this value");
216             return kFALSE;
217          }
218          valueX = vtrack->GetTPCmomentum();
219          valueY = vtrack->GetTPCsignal();
220          ((TH2F *)obj)->Fill(valueX, valueY);
221          return kTRUE;
222       case ktimeTOFvsPPion:
223          if (!vtrack) {
224             AliWarning("Required vtrack for this value");
225             return kFALSE;
226          }
227          valueX = vtrack->P();
228          //valueY = vtrack->GetTOFsignal();
229          valueY = 1E20;
230          if (pid) valueY = pid->NumberOfSigmasTOF(vtrack, AliPID::kPion);
231          ((TH2F *)obj)->Fill(valueX, valueY);
232          return kTRUE;
233       case ktimeTOFvsPKaon:
234          if (!vtrack) {
235             AliWarning("Required vtrack for this value");
236             return kFALSE;
237          }
238          valueX = vtrack->P();
239          //valueY = vtrack->GetTOFsignal();
240          valueY = 1E20;
241          if (pid) valueY = pid->NumberOfSigmasTOF(vtrack, AliPID::kKaon);
242          ((TH2F *)obj)->Fill(valueX, valueY);
243          return kTRUE;
244       case ktimeTOFvsPProton:
245          if (!vtrack) {
246             AliWarning("Required vtrack for this value");
247             return kFALSE;
248          }
249          valueX = vtrack->P();
250          //valueY = vtrack->GetTOFsignal();
251          valueY = 1E20;
252          if (pid) valueY = pid->NumberOfSigmasTOF(vtrack, AliPID::kProton);
253          ((TH2F *)obj)->Fill(valueX, valueY);
254          return kTRUE;
255       default:
256          AliError("Invalid value type");
257          return kFALSE;
258    }
259 }