]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliForwardMultiplicityBase.h
Fixed references from PWG2 -> PWGLF - very efficiently done using ETags.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardMultiplicityBase.h
1 // 
2 // Base class for classes that calculate the multiplicity in the
3 // forward regions event-by-event
4 // 
5 #ifndef ALIFORWARDMULTIPLICITYBASE_H
6 #define ALIFORWARDMULTIPLICITYBASE_H
7 /**
8  * @file   AliForwardMultiplicityBase.h
9  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
10  * @date   Wed Mar 23 14:06:29 2011
11  * 
12  * @brief  
13  * 
14  * 
15  * @ingroup pwglf_forward_aod
16  */
17 #include <AliAnalysisTaskSE.h>
18 class AliFMDEventInspector;
19 class AliFMDEnergyFitter;
20 class AliFMDSharingFilter;
21 class AliFMDDensityCalculator;
22 class AliFMDCorrector;
23 class AliFMDHistCollector;
24 class AliForwardCorrectionManager;
25 class AliESDEvent;
26 class TH2D;
27 class TList;
28 class TTree;
29 class TAxis;
30
31 /** 
32  * @mainpage ALICE PWGLF Forward Multiplcity Analysis 
33  * 
34  * This is the analysis code for analysis of the Forward data. 
35  * 
36  * @par Code overview 
37  * 
38  * See the <a href="modules.html">Modules</a> page 
39  * 
40  * @par Run.sh script 
41  * 
42  * @verbatim 
43  * Usage: Run.sh [OPTIONS]
44  * 
45  * Do Pass1 and Pass2 on ESD files in current directory.  
46  * 
47  * Options:
48  *      -h,--help               This help                  
49  *      -n,--events N           Number of events            (-1)
50  *      -1,--pass1              Run pass 1, only AOD        (0)
51  *      -2,--pass2              Run pass 2, only Hists      (0)
52  *      -3,--pass3              Draw results                (0)
53  *      -v,--vz-min CM          Minimum value of vz         (-10)
54  *      -V,--vz-max CM          Maximum value of vz         (10)
55  *      -t,--trigger TYPE       Select trigger TYPE         (INEL)
56  *      -b,--batch              Do batch processing         (0)
57  *      -P,--proof NWORKERS     Run in PROOF(Lite) mode     (0)
58  *      -M,--mc                 Run over MC data            (0)
59  *      -g,--gdb                Run in GDB mode             (0)
60  *      -E,--eloss              Run energy loss script      
61  *      -r,--rebin              Rebin factor                (1)
62  *      -C,--use-centrality     Run centrality task         (0)
63  *      -O,--show-older         Show older data             (0)
64  *      -J,--show-published     Show ALICE published data   (1)
65  *      -R,--show-ratios        Show ratios to other data   (1)
66  *      -Z,--show-asymmetry     Show asymmetry              (1)
67  *      -S,--scheme SCHEME      Normalisation scheme        (full)
68  *      -T,--title STRING       Title on plots              ()
69  * 
70  * TYPE is a comma or space separated list of 
71  *  
72  *   INEL             Inelastic triggers (V0A|V0C|SPD)
73  *   INEL>0      As above + N_ch > 0 in -0.5<eta<+0.5
74  *   NSD         Non-single diffractive ((VOA&VOC)|N_ch > 5 -1.9<eta<+1.9)
75  * 
76  * SCHEME is a comma or space separated list of 
77  * 
78  *   NONE          No event-level normalization except trivial one 
79  *   EVENTLEVEL    Event-level normalization 
80  *   ALTEVENTLEVEL Event-level normalization (alternative version)
81  *   BACKGROUND    Not implemented yet 
82  *   SHAPE         Shape correction 
83  *   FULL          Same as EVENTLEVEL,BACKGROUND,SHAPE
84  *   ALTFULL       Same as ALTEVENTLEVEL,BACKGROUND,SHAPE
85  * 
86  * If NWORKERS is 0, then the analysis will be run in local mode.
87  * @endverbatim
88  */
89 /** 
90  * @defgroup pwglf_forward PWGLF Forward analysis
91  *
92  * Code to do the multiplicity analysis in the forward psuedo-rapidity
93  * regions
94  *
95  */
96 /** 
97  * @defgroup pwglf_forward_tasks Tasks
98  *
99  * Code to do the multiplicity analysis in the forward psuedo-rapidity
100  * regions
101  *
102  * @ingroup pwglf_forward 
103  */
104 /** 
105  * @defgroup pwglf_forward_topical Topical
106  */
107 /** 
108  * @defgroup pwglf_forward_aod AOD
109  * @ingroup pwglf_forward_topical
110  */
111 /** 
112  * Base class for classes that calculate the multiplicity in the
113  * forward regions event-by-event
114  * 
115  * @par Inputs: 
116  *   - AliESDEvent 
117  *
118  * @par Outputs: 
119  *   - AliAODForwardMult 
120  * 
121  * @par Histograms 
122  *   
123  * @par Corrections used 
124  * 
125  * @ingroup pwglf_forward_tasks
126  * @ingroup pwglf_forward_aod
127  * 
128  */
129 class AliForwardMultiplicityBase : public AliAnalysisTaskSE
130 {
131 public:
132   /** 
133    * @{ 
134    * @name Interface methods 
135    */
136   /** 
137    * Initialize the task 
138    * 
139    */
140   virtual void Init() { fFirstEvent = true; }
141   /** 
142    * Create output objects 
143    * 
144    */
145   virtual void UserCreateOutputObjects() = 0;
146   /** 
147    * Process each event 
148    *
149    * @param option Not used
150    */  
151   virtual void UserExec(Option_t* option) = 0;
152   /** 
153    * End of job
154    * 
155    * @param option Not used 
156    */
157   virtual void Terminate(Option_t* option) = 0;
158   /** 
159    * @} 
160    */
161   /** 
162    * Print information 
163    * 
164    * @param option Not used
165    */
166   virtual void Print(Option_t* option="") const;
167   /** 
168    * Whether to enable low-flux code 
169    * 
170    * @param use IF true, enable low-flux code 
171    */
172   virtual void SetEnableLowFlux(Bool_t use=true) { fEnableLowFlux = use; }
173   /** 
174    * @{ 
175    * @name Access to sub-algorithms 
176    */
177   /**
178    * Get reference to the EventInspector algorithm 
179    * 
180    * @return Reference to AliFMDEventInspector object 
181    */
182   virtual AliFMDEventInspector& GetEventInspector() = 0;
183   /**
184    * Get reference to the SharingFilter algorithm 
185    * 
186    * @return Reference to AliFMDSharingFilter object 
187    */
188   virtual AliFMDSharingFilter& GetSharingFilter() = 0;
189   /**
190    * Get reference to the DensityCalculator algorithm 
191    * 
192    * @return Reference to AliFMDDensityCalculator object 
193    */
194   virtual AliFMDDensityCalculator& GetDensityCalculator() = 0;
195   /**
196    * Get reference to the Corrections algorithm 
197    * 
198    * @return Reference to AliFMDCorrector object 
199    */
200   virtual AliFMDCorrector& GetCorrections() = 0;
201   /**
202    * Get reference to the HistCollector algorithm 
203    * 
204    * @return Reference to AliFMDHistCollector object 
205    */
206   virtual AliFMDHistCollector& GetHistCollector() = 0;
207   /**
208    * Get reference to the EventInspector algorithm 
209    * 
210    * @return Reference to AliFMDEventInspector object 
211    */
212   virtual const AliFMDEventInspector& GetEventInspector() const = 0;
213   /**
214    * Get reference to the SharingFilter algorithm 
215    * 
216    * @return Reference to AliFMDSharingFilter object 
217    */
218   virtual const AliFMDSharingFilter& GetSharingFilter() const = 0;
219   /**
220    * Get reference to the DensityCalculator algorithm 
221    * 
222    * @return Reference to AliFMDDensityCalculator object 
223    */
224   virtual const AliFMDDensityCalculator& GetDensityCalculator() const = 0;
225   /**
226    * Get reference to the Corrections algorithm 
227    * 
228    * @return Reference to AliFMDCorrector object 
229    */
230   virtual const AliFMDCorrector& GetCorrections() const = 0;
231   /**
232    * Get reference to the HistCollector algorithm 
233    * 
234    * @return Reference to AliFMDHistCollector object 
235    */
236   virtual const AliFMDHistCollector& GetHistCollector() const = 0;
237   /** 
238    * @} 
239    */
240   virtual void SetDebug(Int_t dbg) = 0;
241 protected: 
242   /** 
243    * Constructor 
244    * 
245    * @param name Name of task 
246    */
247   AliForwardMultiplicityBase(const char* name); 
248   /** 
249    * Constructor
250    */
251   AliForwardMultiplicityBase() 
252   : AliAnalysisTaskSE(), 
253     fEnableLowFlux(true), 
254     fFirstEvent(true),
255     fCorrManager(0)
256   {}
257   /** 
258    * Copy constructor 
259    * 
260    * @param o Object to copy from 
261    */
262   AliForwardMultiplicityBase(const AliForwardMultiplicityBase& o)
263     : AliAnalysisTaskSE(o),
264       fEnableLowFlux(o.fEnableLowFlux), 
265       fFirstEvent(o.fFirstEvent),
266       fCorrManager(o.fCorrManager)
267   {}
268   /** 
269    * Assignment operator 
270    * 
271    * @param o Object to assign from 
272    * 
273    * @return Reference to this object 
274    */
275   AliForwardMultiplicityBase& operator=(const AliForwardMultiplicityBase& o);
276   /** 
277    * Check if all needed corrections are there and accounted for.  If not,
278    * do a Fatal exit 
279    * 
280    * @param what Which corrections is needed
281    * 
282    * @return true if all present, false otherwise
283    */  
284   Bool_t CheckCorrections(UInt_t what) const;
285   /**
286    * Read corrections
287    *
288    */
289   virtual Bool_t ReadCorrections(const TAxis*& pe, 
290                                  const TAxis*& pv,
291                                  Bool_t mc=false);
292   /**
293    * Get the ESD event. IF this is the first event, initialise
294    */
295   virtual AliESDEvent* GetESDEvent();
296   /** 
297    * Initialise the sub objects and stuff.  Called on first event
298    *
299    */
300   virtual void InitializeSubs() = 0;
301   /**
302    * Mark this event as one to store in the AOD 
303    * 
304    */
305   virtual void MarkEventForStore() const;
306   /** 
307    * Make Ring @f$ dN/d\eta @f$ histogram and a stack 
308    * 
309    * @param input      List with summed signals 
310    * @param output     Output list 
311    * @param inName     Input name 
312    * @param outName    Output name
313    * @param style      Style 
314    */
315   virtual void MakeRingdNdeta(const TList* input, 
316                               const char*  inName,
317                               TList*       output,
318                               const char*  outName,
319                               Int_t        style=20) const;
320   Bool_t fEnableLowFlux;// Whether to use low-flux specific code
321   Bool_t fFirstEvent;   // Whether the event is the first seen 
322 private:
323   /**
324    * A pointer to the corrections manager.  This is here to make the
325    * corrections manager persistent - that is, when we write the
326    * analysis train to a file (as done in PROOF) we should also write
327    * down the corrections mananger.   This pointer ensures that. 
328    * 
329    */
330   AliForwardCorrectionManager* fCorrManager; // Pointer to corrections manager
331
332   ClassDef(AliForwardMultiplicityBase,2) // Forward multiplicity class
333 };
334
335 #endif
336 // Local Variables:
337 //  mode: C++
338 // End:
339