]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnMiniMonitor.cxx
two track study as function of pt
[u/mrichter/AliRoot.git] / PWG2 / 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
48 #include "AliRsnDaughter.h"
49 #include "AliRsnEvent.h"
50
51 #include "AliRsnMiniMonitor.h"
52
53 ClassImp(AliRsnMiniMonitor)
54
55 //__________________________________________________________________________________________________
56 AliRsnMiniMonitor::AliRsnMiniMonitor() :
57    TNamed(), 
58    fType(kTypes), 
59    fCutID(-1), 
60    fListID(-1), 
61    fList(0x0)
62 {
63 //
64 // Dummy constructor
65 //
66 }
67
68 //__________________________________________________________________________________________________
69 AliRsnMiniMonitor::AliRsnMiniMonitor(const char *name, EType type, Int_t cutID) :
70    TNamed(name, ""), 
71    fType(type), 
72    fCutID(cutID), 
73    fListID(-1), 
74    fList(0x0)
75 {
76 //
77 // Default constructor
78 //
79 }
80    
81 //__________________________________________________________________________________________________
82 AliRsnMiniMonitor::AliRsnMiniMonitor(const AliRsnMiniMonitor& copy) :
83    TNamed(copy), 
84    fType(copy.fType), 
85    fCutID(copy.fCutID), 
86    fListID(copy.fListID), 
87    fList(copy.fList)
88 {
89 //
90 // Copy constructor
91 //
92 }
93    
94 //__________________________________________________________________________________________________
95 AliRsnMiniMonitor& AliRsnMiniMonitor::operator=(const AliRsnMiniMonitor& copy) 
96 {
97 //
98 // Assignment operator
99 //
100    
101    TNamed::operator=(copy); 
102    fType = copy.fType; 
103    fCutID = copy.fCutID; 
104    fListID = copy.fListID; 
105    fList = copy.fList; 
106    
107    return (*this); 
108 }
109
110 //__________________________________________________________________________________________________
111 Bool_t AliRsnMiniMonitor::Init(const char *name, TList *list)
112 {
113 //
114 // Initialize this output histogram and put into the passed list
115 //
116
117    // name
118    TString sname(name);
119    sname += '_';
120    sname += GetName();
121    
122    // check list
123    fList = list;
124    if (!list) {
125       AliError("No list!");
126       return kFALSE;
127    }
128    
129    // reset histogram
130    TH1 *histogram = 0x0;
131
132    switch (fType) {
133       case kdEdxTPCvsP: 
134          sname += "_TPCsignal";
135          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, 0.0, 1000.0);
136          break;
137       case ktimeTOFvsP:
138          sname += "_TOFsignal";
139          histogram = new TH2F(sname.Data(), "", 500, 0.0, 5.0, 1000, -5.0, 5.0);
140          break;
141       default:
142          AliError("Wrong enum type");
143          return kFALSE;
144    }
145    
146    // add to list
147    if (histogram && fList) {
148       histogram->Sumw2();
149       fList->Add(histogram);
150       fListID = fList->IndexOf(histogram);
151       AliInfo(Form("Histogram '%s' added to list in slot #%d", histogram->GetName(), fListID));
152       return kTRUE;
153    }
154    
155    return kFALSE;
156 }
157
158 //_____________________________________________________________________________
159 Bool_t AliRsnMiniMonitor::Fill(AliRsnDaughter *track, AliRsnEvent *)
160 {
161 //
162 // Fill the histogram
163 //
164
165    // retrieve object from list
166    if (!fList) {
167       AliError("List pointer is NULL");
168       return kFALSE;
169    }
170    TObject *obj = fList->At(fListID);
171    if (!obj) {
172       AliError("List object is NULL");
173       return kFALSE;
174    }
175
176    Double_t valueX, valueY;
177    AliVTrack *vtrack = track->Ref2Vtrack();
178
179    switch (fType) {
180       case kdEdxTPCvsP: 
181          if (!vtrack) {
182             AliWarning("Required vtrack for this value");
183             return kFALSE;
184          }
185          valueX = vtrack->GetTPCmomentum();
186          valueY = vtrack->GetTPCsignal();
187          ((TH2F*)obj)->Fill(valueX, valueY);
188          return kTRUE;
189       case ktimeTOFvsP:
190          if (!vtrack) {
191             AliWarning("Required vtrack for this value");
192             return kFALSE;
193          }
194          valueX = vtrack->P();
195          valueY = vtrack->GetTOFsignal();
196          ((TH2F*)obj)->Fill(valueX, valueY);
197          return kTRUE;
198       default:
199          AliError("Invalid value type");
200          return kFALSE;
201    }
202 }