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.
25 * TFile* file = TFile::Open("AliAODs.root","READ");
26 * TTree* tree = static_cast<TTree*>(file->Get("aodTree"));
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
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"));
50 * // Other trigger/event requirements could be defined
51 * if (!mult->IsTriggerBits(mask)) continue;
54 * // Check if we have vertex
55 * if (!mult->HasIpZ()) continue;
58 * // Select vertex range (in centimeters)
59 * if (!mult->InRange(vzLow, vzHigh) continue;
62 * // Add contribution from this event
63 * sum->Add(&(mult->GetHistogram()));
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
73 * dndeta->Scale(Double_t(nWithVertex)/nTriggered, "width");
74 * // And draw the result
79 * The above code will draw the final @f$ dN_{ch}/d\eta@f$ for the
80 * selected event class and vertex range
82 * The histogram can be used as input for other kinds of analysis too,
83 * like flow, event-plane, centrality, and so on.
85 * @ingroup pwg2_forward_analysis
87 class AliAODForwardMult : public TObject
91 * Bits of the trigger pattern
94 /** In-elastic collision */
96 /** In-elastic collision with at least one SPD tracklet */
98 /** Non-single diffractive collision */
100 /** Empty bunch crossing */
102 /** A-side trigger */
104 /** B(arrel) trigger */
106 /** C-side trigger */
112 * Default constructor
114 * Used by ROOT I/O sub-system - do not use
121 AliAODForwardMult(Bool_t);
125 ~AliAODForwardMult() {}
129 * @param etaAxis Pseudo-rapidity axis
131 void Init(const TAxis& etaAxis);
133 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
135 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
137 const TH2D& GetHistogram() const { return fHist; }
139 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
141 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
143 TH2D& GetHistogram() { return fHist; }
145 * Get the trigger mask
147 * @return Trigger mask
149 UInt_t GetTriggerMask() const { return fTriggers; }
151 * Set the trigger mask
153 * @param trg Trigger mask
155 void SetTriggerMask(UInt_t trg) { fTriggers = trg; }
157 * Set bit(s) in trigger mask
159 * @param bits bit(s) to set
161 void SetTriggerBits(UInt_t bits) { fTriggers |= bits; }
163 * Check if bit(s) are set in the trigger mask
165 * @param bits Bits to test for
169 Bool_t IsTriggerBits(UInt_t bits) const;
171 * Whether we have any trigger bits
173 Bool_t HasTrigger() const { return fTriggers != 0; }
177 * @param option Passed on to TH2::Reset verbatim
179 void Clear(Option_t* option="");
185 void Browse(TBrowser* b);
189 * @return Always true
191 Bool_t IsFolder() const { return kTRUE; }
195 * @param option Passed verbatim to TH2::Print
197 void Print(Option_t* option="") const;
199 * Set the z coordinate of the interaction point
201 * @param ipZ Interaction point z coordinate
203 void SetIpZ(Float_t ipZ) { fIpZ = ipZ; }
205 * Set the z coordinate of the interaction point
207 * @return Interaction point z coordinate
209 Float_t GetIpZ() const { return fIpZ; }
211 * Check if we have a valid z coordinate of the interaction point
213 * @return True if we have a valid interaction point z coordinate
215 Bool_t HasIpZ() const;
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.
221 * @param low Lower cut (inclusive)
222 * @param high Upper cut (exclusive)
224 * @return true if @f$ low \ge ipz < high@f$
226 Bool_t InRange(Float_t low, Float_t high) const;
228 * Get the name of the object
230 * @return Name of object
232 const Char_t* GetName() const { return "Forward"; }
234 * Get a string correspondig to the trigger mask
236 * @param mask Trigger mask
238 * @return Static string (copy before use)
240 static const Char_t* GetTriggerString(UInt_t mask);
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
246 static const Float_t fgkInvalidIpZ; // Invalid IpZ value
247 ClassDef(AliAODForwardMult,1); // AOD forward multiplicity
250 //____________________________________________________________________
252 AliAODForwardMult::InRange(Float_t low, Float_t high) const
254 return HasIpZ() && fIpZ >= low && fIpZ < high;
257 //____________________________________________________________________
259 AliAODForwardMult::IsTriggerBits(UInt_t bits) const
261 return HasTrigger() && ((fTriggers & bits) != 0);