]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/AliAODForwardMult.h
Small updates
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliAODForwardMult.h
CommitLineData
7e4038b5 1#ifndef ALIROOT_PWG2_FORWARD_ANALYSIS_ALIAODFORWARDMULT_H
2#define ALIROOT_PWG2_FORWARD_ANALYSIS_ALIAODFORWARDMULT_H
3#include <TObject.h>
4#include <TH2D.h>
5class TBrowser;
6/**
7 * Class that contains the forward multiplicity data per event
8 *
9 * This class contains a histogram of
10 * @f[
11 * \frac{d^2N_{ch}}{d\eta d\phi}\quad,
12 * @f]
13 * as well as a trigger mask for each analysed event.
14 *
15 * The eta acceptance of the event is stored in the underflow bins of
16 * the histogram. So to build the final histogram, one needs to
17 * correct for this acceptance (properly weighted by the events), and
18 * the vertex efficiency. This simply boils down to defining a 2D
19 * histogram and summing the event histograms in that histogram. One
20 * should of course also do proper book-keeping of the accepted event.
21 *
22 * @code
fa4236ed 23 * TTree* GetAODTree()
24 * {
25 * TFile* file = TFile::Open("AliAODs.root","READ");
26 * TTree* tree = static_cast<TTree*>(file->Get("aodTree"));
27 * return tree;
28 * }
29 *
30 * void Analyse()
31 * {
32 * TH2D* sum = 0; // Summed hist
33 * TTree* tree = GetAODTree(); // AOD tree
34 * AliAODForwardMult* mult = 0; // AOD object
35 * Int_t nTriggered = 0; // # of triggered ev.
36 * Int_t nWithVertex= 0; // # of ev. w/vertex
37 * Int_t nAccepted = 0; // # of ev. used
38 * Int_t nAvailable = tree->GetEntries(); // How many entries
39 * Float_t vzLow = -10; // Lower ip cut
40 * Float_t vzHigh = 10; // Upper ip cut
41 * Int_t mask = AliAODForwardMult::kInel;// Trigger mask
42 * tree->SetBranchAddress("forward", &forward); // Set the address
7e4038b5 43 *
44 * for (int i = 0; i < nAvailable; i++) {
45 * // Create sum histogram on first event - to match binning to input
46 * if (!sum) sum = static_cast<TH2D*>(mult->Clone("d2ndetadphi"));
47 *
48 * tree->GetEntry(i);
49 *
50 * // Other trigger/event requirements could be defined
51 * if (!mult->IsTriggerBits(mask)) continue;
52 * nTriggered++;
fa4236ed 53 *
54 * // Check if we have vertex
55 * if (!mult->HasIpZ()) continue;
56 * nWithVertex++;
7e4038b5 57 *
58 * // Select vertex range (in centimeters)
59 * if (!mult->InRange(vzLow, vzHigh) continue;
60 * nAccepted++;
61 *
62 * // Add contribution from this event
63 * sum->Add(&(mult->GetHistogram()));
64 * }
65 *
66 * // Get acceptance normalisation from underflow bins
67 * TH1D* norm = sum->Projection("norm", 0, 1, "");
68 * // Project onto eta axis - _ignoring_underflow_bins_!
69 * TH1D* dndeta = sum->Projection("dndeta", 1, -1, "e");
70 * // Normalize to the acceptance
71 * dndeta->Divide(norm);
72 * // Scale by the vertex efficiency
fa4236ed 73 * dndeta->Scale(Double_t(nWithVertex)/nTriggered, "width");
7e4038b5 74 * // And draw the result
75 * dndeta->Draw();
fa4236ed 76 * }
7e4038b5 77 * @endcode
78 *
79 * The above code will draw the final @f$ dN_{ch}/d\eta@f$ for the
80 * selected event class and vertex range
81 *
82 * The histogram can be used as input for other kinds of analysis too,
83 * like flow, event-plane, centrality, and so on.
84 *
7c1a1f1d 85 * @ingroup pwg2_forward
7e4038b5 86 */
87class AliAODForwardMult : public TObject
88{
89public:
90 /**
91 * Bits of the trigger pattern
92 */
93 enum {
94 /** In-elastic collision */
95 kInel = 0x001,
96 /** In-elastic collision with at least one SPD tracklet */
97 kInelGt0 = 0x002,
98 /** Non-single diffractive collision */
99 kNSD = 0x004,
100 /** Empty bunch crossing */
101 kEmpty = 0x008,
102 /** A-side trigger */
103 kA = 0x010,
104 /** B(arrel) trigger */
105 kB = 0x020,
106 /** C-side trigger */
107 kC = 0x080,
108 /** Empty trigger */
109 kE = 0x100
110 };
111 /**
112 * Default constructor
113 *
114 * Used by ROOT I/O sub-system - do not use
115 */
116 AliAODForwardMult();
117 /**
118 * Constructor
119 *
120 */
121 AliAODForwardMult(Bool_t);
122 /**
123 * Destructor
124 */
125 ~AliAODForwardMult() {}
126 /**
127 * Initialize
128 *
129 * @param etaAxis Pseudo-rapidity axis
130 */
131 void Init(const TAxis& etaAxis);
132 /**
133 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
134 *
135 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
136 */
137 const TH2D& GetHistogram() const { return fHist; }
138 /**
139 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
140 *
141 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
142 */
143 TH2D& GetHistogram() { return fHist; }
144 /**
145 * Get the trigger mask
146 *
147 * @return Trigger mask
148 */
149 UInt_t GetTriggerMask() const { return fTriggers; }
150 /**
151 * Set the trigger mask
152 *
153 * @param trg Trigger mask
154 */
155 void SetTriggerMask(UInt_t trg) { fTriggers = trg; }
156 /**
157 * Set bit(s) in trigger mask
158 *
159 * @param bits bit(s) to set
160 */
161 void SetTriggerBits(UInt_t bits) { fTriggers |= bits; }
162 /**
163 * Check if bit(s) are set in the trigger mask
164 *
165 * @param bits Bits to test for
166 *
167 * @return
168 */
169 Bool_t IsTriggerBits(UInt_t bits) const;
170 /**
171 * Whether we have any trigger bits
172 */
173 Bool_t HasTrigger() const { return fTriggers != 0; }
174 /**
175 * Clear all data
176 *
177 * @param option Passed on to TH2::Reset verbatim
178 */
179 void Clear(Option_t* option="");
180 /**
181 * browse this object
182 *
183 * @param b Browser
184 */
185 void Browse(TBrowser* b);
186 /**
187 * This is a folder
188 *
189 * @return Always true
190 */
191 Bool_t IsFolder() const { return kTRUE; }
192 /**
193 * Print content
194 *
195 * @param option Passed verbatim to TH2::Print
196 */
197 void Print(Option_t* option="") const;
198 /**
199 * Set the z coordinate of the interaction point
200 *
201 * @param ipZ Interaction point z coordinate
202 */
203 void SetIpZ(Float_t ipZ) { fIpZ = ipZ; }
204 /**
205 * Set the z coordinate of the interaction point
206 *
207 * @return Interaction point z coordinate
208 */
209 Float_t GetIpZ() const { return fIpZ; }
210 /**
211 * Check if we have a valid z coordinate of the interaction point
212 *
213 * @return True if we have a valid interaction point z coordinate
214 */
215 Bool_t HasIpZ() const;
216 /**
217 * Check if the z coordinate of the interaction point is within the
218 * given limits. Note that the convention used corresponds to the
219 * convention used in ROOTs TAxis.
220 *
221 * @param low Lower cut (inclusive)
222 * @param high Upper cut (exclusive)
223 *
224 * @return true if @f$ low \ge ipz < high@f$
225 */
226 Bool_t InRange(Float_t low, Float_t high) const;
227 /**
228 * Get the name of the object
229 *
230 * @return Name of object
231 */
232 const Char_t* GetName() const { return "Forward"; }
233 /**
234 * Get a string correspondig to the trigger mask
235 *
236 * @param mask Trigger mask
237 *
238 * @return Static string (copy before use)
239 */
240 static const Char_t* GetTriggerString(UInt_t mask);
241protected:
242 TH2D fHist; // Histogram of d^2N_{ch}/(deta dphi) for this event
243 UInt_t fTriggers; // Trigger bit mask
244 Float_t fIpZ; // Z coordinate of the interaction point
245
246 static const Float_t fgkInvalidIpZ; // Invalid IpZ value
247 ClassDef(AliAODForwardMult,1); // AOD forward multiplicity
248};
249
250//____________________________________________________________________
251inline Bool_t
252AliAODForwardMult::InRange(Float_t low, Float_t high) const
253{
254 return HasIpZ() && fIpZ >= low && fIpZ < high;
255}
256
257//____________________________________________________________________
258inline Bool_t
259AliAODForwardMult::IsTriggerBits(UInt_t bits) const
260{
261 return HasTrigger() && ((fTriggers & bits) != 0);
262}
263
264
265#endif
266// Local Variables:
267// mode: C++
268// End:
269