]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliAODForwardMult.h
Fix compilation warning.
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliAODForwardMult.h
1 //
2 // See implementation or Doxygen comments for more information
3 //
4 #ifndef ALIAODFORWARDMULT_H
5 #define ALIAODFORWARDMULT_H
6 #include <TObject.h>
7 #include <TH2D.h>
8 class TBrowser;
9 /**
10  * Class that contains the forward multiplicity data per event 
11  *
12  * This class contains a histogram of 
13  * @f[
14  *   \frac{d^2N_{ch}}{d\eta d\phi}\quad,
15  * @f]
16  * as well as a trigger mask for each analysed event.  
17  * 
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.
24  *
25  * @code 
26  * TTree* GetAODTree()
27  * { 
28  *    TFile* file = TFile::Open("AliAODs.root","READ");
29  *    TTree* tree = static_cast<TTree*>(file->Get("aodTree"));
30  *    return tree;
31  * }
32  * 
33  * void Analyse()
34  * { 
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
46  * 
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"));
50  * 
51  *     tree->GetEntry(i);
52  * 
53  *     // Other trigger/event requirements could be defined 
54  *     if (!mult->IsTriggerBits(mask)) continue; 
55  *     nTriggered++;
56  *
57  *     // Check if we have vertex 
58  *     if (!mult->HasIpZ()) continue;
59  *     nWithVertex++;
60  * 
61  *     // Select vertex range (in centimeters) 
62  *     if (!mult->InRange(vzLow, vzHigh) continue; 
63  *     nAccepted++;
64  * 
65  *     // Add contribution from this event
66  *     sum->Add(&(mult->GetHistogram()));
67  *   }
68  * 
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
78  *   dndeta->Draw();
79  * }
80  * @endcode   
81  *     
82  * The above code will draw the final @f$ dN_{ch}/d\eta@f$ for the
83  * selected event class and vertex range
84  *
85  * The histogram can be used as input for other kinds of analysis too, 
86  * like flow, event-plane, centrality, and so on. 
87  *
88  * @ingroup pwg2_forward 
89  */
90 class AliAODForwardMult : public TObject
91 {
92 public:
93   /** 
94    * Bits of the trigger pattern
95    */
96   enum { 
97     /** In-elastic collision */
98     kInel     = 0x001, 
99     /** In-elastic collision with at least one SPD tracklet */
100     kInelGt0  = 0x002, 
101     /** Non-single diffractive collision */
102     kNSD      = 0x004, 
103     /** Empty bunch crossing */
104     kEmpty    = 0x008, 
105     /** A-side trigger */
106     kA        = 0x010, 
107     /** B(arrel) trigger */
108     kB        = 0x020, 
109     /** C-side trigger */
110     kC        = 0x080,  
111     /** Empty trigger */
112     kE        = 0x100
113   };
114   /** 
115    * Default constructor 
116    * 
117    * Used by ROOT I/O sub-system - do not use
118    */
119   AliAODForwardMult();
120   /** 
121    * Constructor 
122    * 
123    * @param isMC Whether this was from MC or not 
124    */
125   AliAODForwardMult(Bool_t isMC);
126   /** 
127    * Destructor 
128    */
129   virtual ~AliAODForwardMult() {} // Destructor 
130   /** 
131    * Initialize 
132    * 
133    * @param etaAxis  Pseudo-rapidity axis
134    */
135   void Init(const TAxis& etaAxis);
136   /** 
137    * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram, 
138    *
139    * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram, 
140    */  
141   const TH2D& GetHistogram() const { return fHist; } // Get histogram 
142   /** 
143    * Get the @f$ d^2N_{ch}/d\eta d\phi@f$ histogram, 
144    *
145    * @return @f$ d^2N_{ch}/d\eta d\phi@f$ histogram, 
146    */  
147   TH2D& GetHistogram() { return fHist; } // Get histogram 
148   /** 
149    * Get the trigger mask 
150    * 
151    * @return Trigger mask 
152    */
153   UInt_t GetTriggerMask() const { return fTriggers; } // Get triggers
154   /** 
155    * Set the trigger mask 
156    * 
157    * @param trg Trigger mask
158    */
159   void SetTriggerMask(UInt_t trg) { fTriggers = trg; } // Set triggers 
160   /** 
161    * Set bit(s) in trigger mask 
162    * 
163    * @param bits bit(s) to set 
164    */
165   void SetTriggerBits(UInt_t bits) { fTriggers |= bits; } // Set trigger bits
166   /** 
167    * Check if bit(s) are set in the trigger mask 
168    * 
169    * @param bits Bits to test for 
170    * 
171    * @return 
172    */
173   Bool_t IsTriggerBits(UInt_t bits) const;
174   /** 
175    * Whether we have any trigger bits 
176    */
177   Bool_t HasTrigger() const { return fTriggers != 0; } // Check for triggers
178   /** 
179    * Clear all data 
180    * 
181    * @param option  Passed on to TH2::Reset verbatim
182    */
183   void Clear(Option_t* option="");
184   /** 
185    * browse this object 
186    * 
187    * @param b Browser 
188    */
189   void Browse(TBrowser* b);
190   /** 
191    * This is a folder 
192    * 
193    * @return Always true
194    */
195   Bool_t IsFolder() const { return kTRUE; } // Always true 
196   /** 
197    * Print content 
198    * 
199    * @param option Passed verbatim to TH2::Print 
200    */
201   void Print(Option_t* option="") const;
202   /** 
203    * Set the z coordinate of the interaction point
204    * 
205    * @param ipZ Interaction point z coordinate
206    */
207   void SetIpZ(Float_t ipZ) { fIpZ = ipZ; } // Set Ip's Z coordinate
208   /** 
209    * Set the z coordinate of the interaction point
210    * 
211    * @return Interaction point z coordinate
212    */
213   Float_t GetIpZ() const { return fIpZ; } // Get Ip's Z coordinate 
214   /** 
215    * Check if we have a valid z coordinate of the interaction point
216    *
217    * @return True if we have a valid interaction point z coordinate
218    */
219   Bool_t HasIpZ() const;
220   /** 
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.
224    * 
225    * @param low  Lower cut (inclusive)
226    * @param high Upper cut (exclusive)
227    * 
228    * @return true if @f$ low \ge ipz < high@f$ 
229    */
230   Bool_t InRange(Float_t low, Float_t high) const;
231   /** 
232    * Get the name of the object 
233    * 
234    * @return Name of object 
235    */
236   const Char_t* GetName() const { return (fIsMC ? "ForwardMC" : "Forward"); }
237   /** 
238    * Get a string correspondig to the trigger mask
239    * 
240    * @param mask Trigger mask 
241    * 
242    * @return Static string (copy before use)
243    */
244   static const Char_t* GetTriggerString(UInt_t mask);
245 protected: 
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
250
251   static const Float_t fgkInvalidIpZ; // Invalid IpZ value 
252   ClassDef(AliAODForwardMult,1); // AOD forward multiplicity 
253 };
254
255 //____________________________________________________________________
256 inline Bool_t
257 AliAODForwardMult::InRange(Float_t low, Float_t high) const 
258 {
259   return HasIpZ() && fIpZ >= low && fIpZ < high;
260 }
261
262 //____________________________________________________________________
263 inline Bool_t 
264 AliAODForwardMult::IsTriggerBits(UInt_t bits) const 
265
266   return HasTrigger() && ((fTriggers & bits) != 0); 
267 }
268
269 #endif
270 // Local Variables:
271 //  mode: C++
272 // End:
273