1 #ifndef ALIROOT_PWG2_FORWARD_ANALYSIS_ALIAODFORWARDMULT_H
2 #define ALIROOT_PWG2_FORWARD_ANALYSIS_ALIAODFORWARDMULT_H
7 * Class that contains the forward multiplicity data per event
9 * This class contains a histogram of
11 * \frac{d^2N_{ch}}{d\eta d\phi}\quad,
13 * as well as a trigger mask for each analysed event.
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.
23 * TH2D* sum = 0; // Summed hist
24 * TTree* tree = GetAODTree(); // AOD tree
25 * AliAODForwardMult* mult = 0; // AOD object
26 * Int_t nTriggered = 0; // # of triggered ev.
27 * Int_t nAccepted = 0; // # of ev. w/vertex
28 * Int_t nAvailable = tree->GetEntries(); // How many entries
29 * Float_t vzLow = -10; // Lower ip cut
30 * Float_t vzHigh = 10; // Upper ip cut
31 * Int_t mask = AliAODForward::kINEL;// Trigger mask
32 * tree->SetBranchAddress("forward", &forward); // Set the address
34 * for (int i = 0; i < nAvailable; i++) {
35 * // Create sum histogram on first event - to match binning to input
36 * if (!sum) sum = static_cast<TH2D*>(mult->Clone("d2ndetadphi"));
40 * // Other trigger/event requirements could be defined
41 * if (!mult->IsTriggerBits(mask)) continue;
44 * // Select vertex range (in centimeters)
45 * if (!mult->InRange(vzLow, vzHigh) continue;
48 * // Add contribution from this event
49 * sum->Add(&(mult->GetHistogram()));
52 * // Get acceptance normalisation from underflow bins
53 * TH1D* norm = sum->Projection("norm", 0, 1, "");
54 * // Project onto eta axis - _ignoring_underflow_bins_!
55 * TH1D* dndeta = sum->Projection("dndeta", 1, -1, "e");
56 * // Normalize to the acceptance
57 * dndeta->Divide(norm);
58 * // Scale by the vertex efficiency
59 * dndeta->Scale(Double_t(nAccepted)/nTriggered, "width");
60 * // And draw the result
64 * The above code will draw the final @f$ dN_{ch}/d\eta@f$ for the
65 * selected event class and vertex range
67 * The histogram can be used as input for other kinds of analysis too,
68 * like flow, event-plane, centrality, and so on.
70 * @ingroup pwg2_forward_analysis
72 class AliAODForwardMult : public TObject
76 * Bits of the trigger pattern
79 /** In-elastic collision */
81 /** In-elastic collision with at least one SPD tracklet */
83 /** Non-single diffractive collision */
85 /** Empty bunch crossing */
89 /** B(arrel) trigger */
99 * Used by ROOT I/O sub-system - do not use
106 AliAODForwardMult(Bool_t);
110 ~AliAODForwardMult() {}
114 * @param etaAxis Pseudo-rapidity axis
116 void Init(const TAxis& etaAxis);
118 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
120 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
122 const TH2D& GetHistogram() const { return fHist; }
124 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
126 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
128 TH2D& GetHistogram() { return fHist; }
130 * Get the trigger mask
132 * @return Trigger mask
134 UInt_t GetTriggerMask() const { return fTriggers; }
136 * Set the trigger mask
138 * @param trg Trigger mask
140 void SetTriggerMask(UInt_t trg) { fTriggers = trg; }
142 * Set bit(s) in trigger mask
144 * @param bits bit(s) to set
146 void SetTriggerBits(UInt_t bits) { fTriggers |= bits; }
148 * Check if bit(s) are set in the trigger mask
150 * @param bits Bits to test for
154 Bool_t IsTriggerBits(UInt_t bits) const;
156 * Whether we have any trigger bits
158 Bool_t HasTrigger() const { return fTriggers != 0; }
162 * @param option Passed on to TH2::Reset verbatim
164 void Clear(Option_t* option="");
170 void Browse(TBrowser* b);
174 * @return Always true
176 Bool_t IsFolder() const { return kTRUE; }
180 * @param option Passed verbatim to TH2::Print
182 void Print(Option_t* option="") const;
184 * Set the z coordinate of the interaction point
186 * @param ipZ Interaction point z coordinate
188 void SetIpZ(Float_t ipZ) { fIpZ = ipZ; }
190 * Set the z coordinate of the interaction point
192 * @return Interaction point z coordinate
194 Float_t GetIpZ() const { return fIpZ; }
196 * Check if we have a valid z coordinate of the interaction point
198 * @return True if we have a valid interaction point z coordinate
200 Bool_t HasIpZ() const;
202 * Check if the z coordinate of the interaction point is within the
203 * given limits. Note that the convention used corresponds to the
204 * convention used in ROOTs TAxis.
206 * @param low Lower cut (inclusive)
207 * @param high Upper cut (exclusive)
209 * @return true if @f$ low \ge ipz < high@f$
211 Bool_t InRange(Float_t low, Float_t high) const;
213 * Get the name of the object
215 * @return Name of object
217 const Char_t* GetName() const { return "Forward"; }
219 * Get a string correspondig to the trigger mask
221 * @param mask Trigger mask
223 * @return Static string (copy before use)
225 static const Char_t* GetTriggerString(UInt_t mask);
227 TH2D fHist; // Histogram of d^2N_{ch}/(deta dphi) for this event
228 UInt_t fTriggers; // Trigger bit mask
229 Float_t fIpZ; // Z coordinate of the interaction point
231 static const Float_t fgkInvalidIpZ; // Invalid IpZ value
232 ClassDef(AliAODForwardMult,1); // AOD forward multiplicity
235 //____________________________________________________________________
237 AliAODForwardMult::InRange(Float_t low, Float_t high) const
239 return HasIpZ() && fIpZ >= low && fIpZ < high;
242 //____________________________________________________________________
244 AliAODForwardMult::IsTriggerBits(UInt_t bits) const
246 return HasTrigger() && ((fTriggers & bits) != 0);