2 // See implementation or Doxygen comments for more information
4 #ifndef ALIAODFORWARDMULT_H
5 #define ALIAODFORWARDMULT_H
10 * Class that contains the forward multiplicity data per event
12 * This class contains a histogram of
14 * \frac{d^2N_{ch}}{d\eta d\phi}\quad,
16 * as well as a trigger mask for each analysed event.
18 * The eta acceptance of the event is stored in the underflow bins of
19 * the histogram. So to build the final histogram, one needs to
20 * correct for this acceptance (properly weighted by the events), and
21 * the vertex efficiency. This simply boils down to defining a 2D
22 * histogram and summing the event histograms in that histogram. One
23 * should of course also do proper book-keeping of the accepted event.
28 * TFile* file = TFile::Open("AliAODs.root","READ");
29 * TTree* tree = static_cast<TTree*>(file->Get("aodTree"));
35 * TH2D* sum = 0; // Summed hist
36 * TTree* tree = GetAODTree(); // AOD tree
37 * AliAODForwardMult* mult = 0; // AOD object
38 * Int_t nTriggered = 0; // # of triggered ev.
39 * Int_t nWithVertex= 0; // # of ev. w/vertex
40 * Int_t nAccepted = 0; // # of ev. used
41 * Int_t nAvailable = tree->GetEntries(); // How many entries
42 * Float_t vzLow = -10; // Lower ip cut
43 * Float_t vzHigh = 10; // Upper ip cut
44 * Int_t mask = AliAODForwardMult::kInel;// Trigger mask
45 * tree->SetBranchAddress("forward", &forward); // Set the address
47 * for (int i = 0; i < nAvailable; i++) {
48 * // Create sum histogram on first event - to match binning to input
49 * if (!sum) sum = static_cast<TH2D*>(mult->Clone("d2ndetadphi"));
53 * // Other trigger/event requirements could be defined
54 * if (!mult->IsTriggerBits(mask)) continue;
57 * // Check if we have vertex
58 * if (!mult->HasIpZ()) continue;
61 * // Select vertex range (in centimeters)
62 * if (!mult->InRange(vzLow, vzHigh) continue;
65 * // Add contribution from this event
66 * sum->Add(&(mult->GetHistogram()));
69 * // Get acceptance normalisation from underflow bins
70 * TH1D* norm = sum->Projection("norm", 0, 1, "");
71 * // Project onto eta axis - _ignoring_underflow_bins_!
72 * TH1D* dndeta = sum->Projection("dndeta", 1, -1, "e");
73 * // Normalize to the acceptance
74 * dndeta->Divide(norm);
75 * // Scale by the vertex efficiency
76 * dndeta->Scale(Double_t(nWithVertex)/nTriggered, "width");
77 * // And draw the result
82 * The above code will draw the final @f$ dN_{ch}/d\eta@f$ for the
83 * selected event class and vertex range
85 * The histogram can be used as input for other kinds of analysis too,
86 * like flow, event-plane, centrality, and so on.
88 * @ingroup pwg2_forward
90 class AliAODForwardMult : public TObject
94 * Bits of the trigger pattern
97 /** In-elastic collision */
99 /** In-elastic collision with at least one SPD tracklet */
101 /** Non-single diffractive collision */
103 /** Empty bunch crossing */
105 /** A-side trigger */
107 /** B(arrel) trigger */
109 /** C-side trigger */
115 * Default constructor
117 * Used by ROOT I/O sub-system - do not use
123 * @param isMC Whether this was from MC or not
125 AliAODForwardMult(Bool_t isMC);
129 virtual ~AliAODForwardMult() {} // Destructor
133 * @param etaAxis Pseudo-rapidity axis
135 void Init(const TAxis& etaAxis);
137 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
139 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
141 const TH2D& GetHistogram() const { return fHist; } // Get histogram
143 * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
145 * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram,
147 TH2D& GetHistogram() { return fHist; } // Get histogram
149 * Get the trigger mask
151 * @return Trigger mask
153 UInt_t GetTriggerMask() const { return fTriggers; } // Get triggers
155 * Set the trigger mask
157 * @param trg Trigger mask
159 void SetTriggerMask(UInt_t trg) { fTriggers = trg; } // Set triggers
161 * Set bit(s) in trigger mask
163 * @param bits bit(s) to set
165 void SetTriggerBits(UInt_t bits) { fTriggers |= bits; } // Set trigger bits
167 * Check if bit(s) are set in the trigger mask
169 * @param bits Bits to test for
173 Bool_t IsTriggerBits(UInt_t bits) const;
175 * Whether we have any trigger bits
177 Bool_t HasTrigger() const { return fTriggers != 0; } // Check for triggers
181 * @param option Passed on to TH2::Reset verbatim
183 void Clear(Option_t* option="");
189 void Browse(TBrowser* b);
193 * @return Always true
195 Bool_t IsFolder() const { return kTRUE; } // Always true
199 * @param option Passed verbatim to TH2::Print
201 void Print(Option_t* option="") const;
203 * Set the z coordinate of the interaction point
205 * @param ipZ Interaction point z coordinate
207 void SetIpZ(Float_t ipZ) { fIpZ = ipZ; } // Set Ip's Z coordinate
209 * Set the z coordinate of the interaction point
211 * @return Interaction point z coordinate
213 Float_t GetIpZ() const { return fIpZ; } // Get Ip's Z coordinate
215 * Check if we have a valid z coordinate of the interaction point
217 * @return True if we have a valid interaction point z coordinate
219 Bool_t HasIpZ() const;
221 * Check if the z coordinate of the interaction point is within the
222 * given limits. Note that the convention used corresponds to the
223 * convention used in ROOTs TAxis.
225 * @param low Lower cut (inclusive)
226 * @param high Upper cut (exclusive)
228 * @return true if @f$ low \ge ipz < high@f$
230 Bool_t InRange(Float_t low, Float_t high) const;
232 * Get the name of the object
234 * @return Name of object
236 const Char_t* GetName() const { return (fIsMC ? "ForwardMC" : "Forward"); }
238 * Get a string correspondig to the trigger mask
240 * @param mask Trigger mask
242 * @return Static string (copy before use)
244 static const Char_t* GetTriggerString(UInt_t mask);
246 Bool_t fIsMC; // Whether this is from MC
247 TH2D fHist; // Histogram of d^2N_{ch}/(deta dphi) for this event
248 UInt_t fTriggers; // Trigger bit mask
249 Float_t fIpZ; // Z coordinate of the interaction point
251 static const Float_t fgkInvalidIpZ; // Invalid IpZ value
252 ClassDef(AliAODForwardMult,1); // AOD forward multiplicity
255 //____________________________________________________________________
257 AliAODForwardMult::InRange(Float_t low, Float_t high) const
259 return HasIpZ() && fIpZ >= low && fIpZ < high;
262 //____________________________________________________________________
264 AliAODForwardMult::IsTriggerBits(UInt_t bits) const
266 return HasTrigger() && ((fTriggers & bits) != 0);