]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDEnergyFitter.h
Improvements
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDEnergyFitter.h
1 //
2 // Class to fit the energy distribution.  
3 //
4 #ifndef ALIFMDENERGYFITTER_H
5 #define ALIFMDENERGYFITTER_H
6 /**
7  * @file   AliFMDEnergyFitter.h
8  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
9  * @date   Wed Mar 23 14:02:23 2011
10  * 
11  * @brief  
12  * 
13  * 
14  * @ingroup pwglf_forward_eloss
15  */
16 #include <TNamed.h>
17 // #include <TH1D.h>
18 #include <TAxis.h>
19 #include <TList.h>
20 #include <TObjArray.h>
21 #include <TClonesArray.h>
22 #include "AliFMDCorrELossFit.h"
23 #include "AliForwardUtil.h"
24 class TH1;
25 class TH2;
26 class AliESDFMD;
27 class TFitResult;
28 class TF1;
29 class TArrayD;
30
31 /**
32  * Class to fit the energy distribution.  
33  *
34  * @par Input: 
35  *    - AliESDFMD object  - from reconstruction
36  *
37  * @par Output: 
38  *    - Lists of histogram - one per ring.  Each list has a number of 
39  *      histograms corresponding to the number of eta bins defined.  
40  *
41  * @par Corrections used: 
42  *    - None
43  *
44  * @image html alice-int-2012-040-eloss_fits.png "Summary of fits"
45  * 
46  * @ingroup pwglf_forward_algo
47  * @ingroup pwglf_forward_eloss
48  */
49 class AliFMDEnergyFitter : public TNamed
50 {
51 public: 
52   /** 
53    * Enumeration of parameters 
54    */
55   enum { 
56     /** Index of pre-constant @f$ C@f$ */
57     kC          = AliForwardUtil::ELossFitter::kC,
58     /** Index of most probable value @f$ \Delta_p@f$ */
59     kDelta      = AliForwardUtil::ELossFitter::kDelta, 
60     /** Index of Landau width @f$ \xi@f$ */
61     kXi         = AliForwardUtil::ELossFitter::kXi, 
62     /** Index of Gaussian width @f$ \sigma@f$ */
63     kSigma      = AliForwardUtil::ELossFitter::kSigma, 
64     /** Index of Gaussian additional width @f$ \sigma_n@f$ */
65     kSigmaN     = AliForwardUtil::ELossFitter::kSigmaN,
66     /** Index of Number of particles @f$ N@f$ */
67     kN          = AliForwardUtil::ELossFitter::kN, 
68     /** Base index of particle strengths @f$ a_i@f$ for 
69         @f$i=2,\ldots,N@f$ */
70     kA          = AliForwardUtil::ELossFitter::kA
71   };
72   /** 
73    * Enumeration of residual methods 
74    */
75   enum EResidualMethod {
76     /** Do not calculate residuals */
77     kNoResiduals = 0, 
78     /** The residuals stored are the difference, and the errors are
79         stored in the error bars of the histogram. */
80     kResidualDifference, 
81     /** The residuals stored are the differences scaled to the error
82         on the data */ 
83     kResidualScaledDifference, 
84     /** The residuals stored are the square difference scale to the
85         square error on the data. */
86     kResidualSquareDifference
87   };
88
89   /**
90    * FMD ring bits for skipping 
91    */
92    enum FMDRingBits { 
93      /** FMD1i */
94      kFMD1I=0x01,
95      /** All of FMD1 */
96      kFMD1 =kFMD1I,
97      /** FMD2i */
98      kFMD2I=0x02,
99      /** FMD2o */
100      kFMD2O=0x04,
101      /** All of FMD2 */
102      kFMD2 =kFMD2I|kFMD2O,
103      /** FMD3i */
104      kFMD3I=0x08,
105      /** FMD3o */
106      kFMD3O=0x10,
107      /** All of FMD3 */
108      kFMD3 =kFMD3I|kFMD3O
109   };
110   /** 
111    * Destructor
112    */
113   virtual ~AliFMDEnergyFitter();
114   /** 
115    * Default Constructor - do not use 
116    */
117   AliFMDEnergyFitter();
118   /** 
119    * Constructor 
120    * 
121    * @param title Title of object  - not significant 
122    */
123   AliFMDEnergyFitter(const char* title);
124
125   // -----------------------------------------------------------------
126   /** 
127    * @{ 
128    * @name Setters of options and parameters 
129    */
130   /** 
131    * Set the eta axis to use.  This will force the code to use this
132    * eta axis definition - irrespective of whatever axis is passed to
133    * the Init member function.  Therefore, this member function can be
134    * used to force another eta axis than one found in the correction
135    * objects. 
136    * 
137    * @param nBins  Number of bins 
138    * @param etaMin Minimum of the eta axis 
139    * @param etaMax Maximum of the eta axis 
140    */
141   void SetEtaAxis(Int_t nBins, Double_t etaMin, Double_t etaMax);
142   /** 
143    * Set the eta axis to use.  This will force the code to use this
144    * eta axis definition - irrespective of whatever axis is passed to
145    * the Init member function.  Therefore, this member function can be
146    * used to force another eta axis than one found in the correction
147    * objects. 
148    * 
149    * @param etaAxis Eta axis to use 
150    */
151   void SetEtaAxis(const TAxis& etaAxis);
152   /** 
153    * Set the centrality bins.  E.g., 
154    * @code 
155    * UShort_t n = 12;
156    * Double_t bins[] = {  0.,  5., 10., 15., 20., 30., 
157    *                     40., 50., 60., 70., 80., 100. };
158    * task->GetFitter().SetCentralityBins(n, bins);
159    * @endcode
160    * 
161    * @param nBins Size of @a bins
162    * @param bins  Bin limits. 
163    */
164   void SetCentralityAxis(UShort_t nBins, Double_t* bins);
165   /** 
166    * Set the low cut used for energy 
167    * 
168    * @param lowCut Low cut
169    */
170   void SetLowCut(Double_t lowCut=0.3) { fLowCut = lowCut; }
171   /** 
172    * Set the number of bins to subtract 
173    * 
174    * @param n 
175    */
176   void SetFitRangeBinWidth(UShort_t n=4) { fFitRangeBinWidth = n; }
177   /** 
178    * Whether or not to enable fitting of the final merged result.  
179    * Note, fitting takes quite a while and one should be careful not to do 
180    * this needlessly 
181    * 
182    * @param doFit Whether to do the fits or not 
183    */
184   void SetDoFits(Bool_t doFit=kTRUE) { fDoFits = doFit; }
185   /** 
186    * Set whether to make the corrections object on the output.  Note,
187    * fits should be enable for this to have any effect.
188    * 
189    * @param doMake If true (false is default), do make the corrections object. 
190    */
191   void SetDoMakeObject(Bool_t doMake=kTRUE) { fDoMakeObject = doMake; }
192   /** 
193    * Set how many particles we will try to fit at most to the data
194    * 
195    * @param n Max number of particle to try to fit 
196    */
197   void SetNParticles(UShort_t n) { fNParticles = (n<1 ? 1 : (n>5 ? 5 : n)); }
198   /** 
199    * Set the minimum number of entries each histogram must have 
200    * before we try to fit our response function to it
201    * 
202    * @param n Minimum number of entries
203    */
204   void SetMinEntries(UShort_t n) { fMinEntries = (n < 1 ? 1 : n); }
205   /**
206    * Set maximum energy loss to consider 
207    *
208    * @param x Maximum energy loss to consider 
209    */
210   void SetMaxE(Double_t x) { fMaxE = x; }
211   /**
212    * Set number of energy loss bins 
213    *
214    * @param x Number of energy loss bins 
215    */
216   void SetNEbins(Int_t x) { fNEbins = x; }
217   /** 
218    * Set the maximum relative error 
219    * 
220    * @param e Maximum relative error 
221    */
222   void SetMaxRelativeParameterError(Double_t e=0.2) { fMaxRelParError = e; }
223   /** 
224    * Set the maximum @f$ \chi^2/\nu@f$ 
225    * 
226    * @param c Maximum @f$ \chi^2/\nu@f$ 
227    */
228   void SetMaxChi2PerNDF(Double_t c=10) { fMaxChi2PerNDF = c; }
229   /** 
230    * Set the least weight
231    * 
232    * @param c Least weight
233    */
234   void SetMinWeight(Double_t c=1e-7) { fMinWeight = c; }
235   /**
236    * Set wheter to use increasing bin sizes 
237    *
238    * @param x Wheter to use increasing bin sizes 
239    */
240   void SetUseIncreasingBins(Bool_t x) { fUseIncreasingBins = x; }
241   /** 
242    * Set whether to make residuals, and in that case how. 
243    *
244    * - Square difference: @f$chi_i^2=(h_i - f(x_i))^2/\delta_i^2@f$ 
245    * - Scaled difference: @f$(h_i - f(x_i))/\delta_i@f$ 
246    * - Difference: @f$(h_i - f(x_i)) \pm\delta_i@f$ 
247    *
248    * where @f$h_i, x_i, \delta_i@f$ is the bin content, bin center,
249    * and bin error for bin @f$i@f$ respectively, and @f$ f@f$ is the
250    * fitted function.
251    * 
252    * @param x Residual method 
253    */
254   void SetStoreResiduals(EResidualMethod x=kResidualDifference) 
255   { 
256     fResidualMethod = x; 
257   }
258   /** 
259    * Set the regularization cut @f$c_{R}@f$.  If a @f$\Delta@f$
260    * distribution has more entries @f$ N_{dist}@f$ than @f$c_{R}@f$,
261    * then we modify the errors of the the distribution with the factor
262    * 
263    * @f[
264    * \sqrt{N_{dist}/c_{R}}
265    * @f]
266    *
267    * to keep the @f$\chi^2/\nu@f$ within resonable limits. 
268    *
269    * The large residuals @f$chi_i^2=(h_i - f(x_i))^2/\delta_i^2@f$
270    * (see also SetStoreResiduals) comes about on the boundary between
271    * the @f$N@f$ and @f$N+1@f$ particle contributions, and seems to
272    * fall off for larger @f$N@f$. This may indicate that there's a
273    * component in the distributions that the function
274    *
275    * @f[
276    *   f(\Delta;\Delta_p,\xi,\sigma,\mathbf{a}) = \sum_i=1^{n} a_i\int
277    *   d\Delta' L(\Delta;\Delta',\xi) G(\Delta';\Delta_p,\sigma)
278    * @f]
279    * 
280    * does not capture.   
281    *
282    * @param cut
283    */
284   void SetRegularizationCut(Double_t cut=3e6) 
285   {
286     fRegularizationCut = cut;
287   }
288   void SetSkips(UShort_t skip) { fSkips = skip; }
289   /** 
290    * Set the debug level.  The higher the value the more output 
291    * 
292    * @param dbg Debug level 
293    */
294   void SetDebug(Int_t dbg=1);
295   /* @} */
296   // -----------------------------------------------------------------
297   /** 
298    * @{ 
299    * @name Processing 
300    */
301   void Init();
302   /** 
303    * Define the output histograms.  These are put in a sub list of the
304    * passed list.   The histograms are merged before the parent task calls 
305    * AliAnalysisTaskSE::Terminate 
306    * 
307    * @param dir Directory to add to 
308    */
309   virtual void CreateOutputObjects(TList* dir);
310   /** 
311    * Initialise the task
312    * 
313    * @param etaAxis The eta axis to use.  Note, that if the eta axis
314    * has already been set (using SetEtaAxis), then this parameter will be 
315    * ignored
316    */
317  virtual void SetupForData(const TAxis& etaAxis);
318   /** 
319    * Fitter the input AliESDFMD object
320    * 
321    * @param input     Input 
322    * @param cent      Event centrality (or < 0 if not valid)
323    * @param empty     Whether the event is 'empty'
324    * 
325    * @return True on success, false otherwise 
326    */
327   virtual Bool_t Accumulate(const AliESDFMD& input, 
328                             Double_t         cent,
329                             Bool_t           empty);
330   /** 
331    * Scale the histograms to the total number of events 
332    * 
333    * @param dir Where the histograms are  
334    */
335   virtual void Fit(const TList* dir);
336   /** 
337    * Generate the corrections object 
338    * 
339    * @param dir List to analyse 
340    */
341   void MakeCorrectionsObject(TList* dir);
342   /** @} */
343   /** 
344    * Print information
345    * 
346    * @param option Not used 
347    */
348   void Print(Option_t* option="") const;
349   /** 
350    * Read the parameters from a list - used when re-running the code 
351    * 
352    * @param list Input list 
353    * 
354    * @return true if the parameter where read 
355    */
356   Bool_t ReadParameters(const TCollection* list);
357 protected:
358   /** 
359    * Copy constructor 
360    * 
361    * @param o Object to copy from 
362    */
363   AliFMDEnergyFitter(const AliFMDEnergyFitter& o);
364   /** 
365    * Assignment operator 
366    * 
367    * @param o Object to assign from 
368    * 
369    * @return Reference to this 
370    */
371   AliFMDEnergyFitter& operator=(const AliFMDEnergyFitter& o);
372
373   /** 
374    * Internal data structure to keep track of the histograms
375    */
376   struct RingHistos : public AliForwardUtil::RingHistos
377   { 
378     typedef AliFMDCorrELossFit::ELossFit ELossFit_t;
379     /** 
380      * Default CTOR
381      */
382     RingHistos();
383     /** 
384      * Constructor
385      * 
386      * @param d detector
387      * @param r ring 
388      */
389     RingHistos(UShort_t d, Char_t r);
390     /** 
391      * Copy constructor - not defined
392      * 
393      * @param o Object to copy from 
394      */
395     RingHistos(const RingHistos& o);
396     /** 
397      * Assignment operator  - not defined
398      * 
399      * @param o Object to assign from 
400      * 
401      * @return Reference to this 
402      */
403     RingHistos& operator=(const RingHistos& o);
404     /** 
405      * Destructor 
406      */
407     ~RingHistos();
408     /** 
409      * Make an axis with increasing bins 
410      * 
411      * @param n    Number of bins 
412      * @param min  Minimum 
413      * @param max  Maximum
414      * 
415      * @return An axis with quadratically increasing bin size 
416      */
417     virtual TArrayD MakeIncreasingAxis(Int_t    n, 
418                                        Double_t min, 
419                                        Double_t max) const;
420     /** 
421      * Make E/E_mip histogram 
422      * 
423      * @param name    Name of histogram
424      * @param title   Title of histogram
425      * @param eAxis   @f$\eta@f$ axis
426      * @param deMax   Maximum energy loss 
427      * @param nDeBins Number energy loss bins 
428      * @param incr    Whether to make bins of increasing size
429      */
430     TH2* Make(const char*  name, 
431               const char*  title, 
432               const TAxis& eAxis, 
433               Double_t     deMax=12, 
434               Int_t        nDeBins=300, 
435               Bool_t       incr=true);
436     /** 
437      * Define outputs
438      * 
439      * @param dir 
440      */
441     virtual void CreateOutputObjects(TList* dir);
442     /** 
443      * Initialise object 
444      * 
445      * @param eAxis      Eta axis
446      * @param cAxis      Centrality axis 
447      * @param maxDE      Max energy loss to consider 
448      * @param nDEbins    Number of bins 
449      * @param useIncrBin Whether to use an increasing bin size 
450      */
451     virtual void SetupForData(const TAxis& eAxis, 
452                               const TAxis& cAxis,
453                               Double_t     maxDE=10, 
454                               Int_t        nDEbins=300, 
455                               Bool_t       useIncrBin=true);
456     /** 
457      * Fill histogram 
458      * 
459      * @param empty  True if event is empty
460      * @param eta    @f$ Eta@f$
461      * @param icent  Centrality bin (1 based)
462      * @param mult   Signal 
463      */
464     virtual void Fill(Bool_t empty, Double_t eta, Int_t icent, Double_t mult);
465     /** 
466      * Get the the 2D histogram eloss name from our sub-list of @a dir
467      * and call the Fit function described below (with &fBest) as last
468      * argument.
469      * 
470      * @param dir         Output list 
471      * @param lowCut      Lower cut 
472      * @param nParticles  Max number of convolved landaus to fit
473      * @param minEntries  Minimum number of entries 
474      * @param minusBins   Number of bins from peak to subtract to 
475      *                    get the fit range 
476      * @param relErrorCut Cut applied to relative error of parameter. 
477      *                    Note, for multi-particle weights, the cut 
478      *                    is loosend by a factor of 2 
479      * @param chi2nuCut   Cut on @f$ \chi^2/\nu@f$ - 
480      *                    the reduced @f$\chi^2@f$ 
481      * @param minWeight   Least weight ot consider
482      * @param regCut      Regularization cut-off
483      * @param residuals   Mode for residual plots
484      *
485      * @return List of fit parameters 
486      */
487     virtual TObjArray* Fit(TList*          dir, 
488                            Double_t        lowCut, 
489                            UShort_t        nParticles,
490                            UShort_t        minEntries,
491                            UShort_t        minusBins,
492                            Double_t        relErrorCut, 
493                            Double_t        chi2nuCut,
494                            Double_t        minWeight,
495                            Double_t        regCut,
496                            EResidualMethod residuals) const;
497     /** 
498      * Get the the 2D histogram @a name from our sub-list of @a
499      * dir. Then for each eta slice, try to fit the energu loss
500      * distribution up to @a nParticles particle responses.
501      *
502      * The fitted distributions (along with the functions fitted) are
503      * stored in a newly created sublist (<i>name</i>Dists).
504      *
505      * The fit parameters are also recorded in the newly created sub-list 
506      * <i>name</i>Results.  
507      *
508      * If @a residuals is not equal to kNoResiduals, then the
509      * residuals of the fits will be stored in the newly created
510      * sub-list <i>name</i>Residuals.
511      *
512      * A histogram named <i>name</i>Status is also generated and
513      * stored in the output list.
514      * 
515      * @param dir         Output list 
516      * @param name        Name of 2D base histogram in list
517      * @param lowCut      Lower cut 
518      * @param nParticles  Max number of convolved landaus to fit
519      * @param minEntries  Minimum number of entries 
520      * @param minusBins   Number of bins from peak to subtract to 
521      *                    get the fit range 
522      * @param relErrorCut Cut applied to relative error of parameter. 
523      *                    Note, for multi-particle weights, the cut 
524      *                    is loosend by a factor of 2 
525      * @param chi2nuCut   Cut on @f$ \chi^2/\nu@f$ - 
526      *                    the reduced @f$\chi^2@f$ 
527      * @param minWeight   Least weight ot consider
528      * @param regCut      Regularization cut-off
529      * @param residuals   Mode for residual plots
530      * @param scaleToPeak If true, scale distribution to peak value
531      * @param best        Optional array to store fits in
532      *
533      * @return List of fit parameters 
534      */
535     virtual TObjArray* FitSlices(TList*          dir, 
536                                  const char*     name,
537                                  Double_t        lowCut, 
538                                  UShort_t        nParticles,
539                                  UShort_t        minEntries,
540                                  UShort_t        minusBins,
541                                  Double_t        relErrorCut, 
542                                  Double_t        chi2nuCut,
543                                  Double_t        minWeight,
544                                  Double_t        regCut,
545                                  EResidualMethod residuals,
546                                  Bool_t          scaleToPeak=true,
547                                  TObjArray*      best=0) const;
548      
549     /** 
550      * Fit a signal histogram.  First, the bin @f$ b_{min}@f$ with
551      * maximum bin content in the range @f$ [E_{min},\infty]@f$ is
552      * found.  Then the fit range is set to the bin range 
553      * @f$ [b_{min}-\Delta b,b_{min}+2\Delta b]@f$, and a 1 
554      * particle signal is fitted to that.  The parameters of that fit 
555      * is then used as seeds for a fit of the @f$ N@f$ particle response 
556      * to the data in the range 
557      * @f$ [b_{min}-\Delta b,N(\Delta_1+\xi_1\log(N))+2N\xi@f$
558      * 
559      * @param dist        Histogram to fit 
560      * @param lowCut      Lower cut @f$ E_{min}@f$ on signal 
561      * @param minEntries  Least number of entries required
562      * @param nParticles  Max number @f$ N@f$ of convolved landaus to fit
563      * @param minusBins   Number of bins @f$ \Delta b@f$ from peak to 
564      *                    subtract to get the fit range 
565      * @param relErrorCut Cut applied to relative error of parameter. 
566      *                    Note, for multi-particle weights, the cut 
567      *                    is loosend by a factor of 2 
568      * @param chi2nuCut   Cut on @f$ \chi^2/\nu@f$ - 
569      *                    the reduced @f$\chi^2@f$ 
570      * @param minWeight   Least weight ot consider
571      * @param regCut      Regularization cut-off
572      * @param scaleToPeak If true, scale distribution to peak value
573      * @param status      On return, contain the status code (0: OK, 1:
574      *                    empty, 2: low statistics, 3: fit failed)
575      * 
576      * @return The best fit function 
577      */
578     virtual ELossFit_t* FitHist(TH1*      dist,
579                                 Double_t  lowCut, 
580                                 UShort_t  nParticles,
581                                 UShort_t  minEntries,
582                                 UShort_t  minusBins,
583                                 Double_t  relErrorCut, 
584                                 Double_t  chi2nuCut,
585                                 Double_t  minWeight,
586                                 Double_t  regCut,
587                                 Bool_t    scaleToPeak,
588                                 UShort_t& status) const;
589     /** 
590      * Find the best fit 
591      * 
592      * @param dist           Histogram 
593      * @param relErrorCut    Cut applied to relative error of parameter. 
594      *                       Note, for multi-particle weights, the cut 
595      *                       is loosend by a factor of 2 
596      * @param chi2nuCut      Cut on @f$ \chi^2/\nu@f$ - 
597      *                       the reduced @f$\chi^2@f$ 
598      * @param minWeightCut   Least valid @f$ a_i@f$ 
599      * 
600      * @return Best fit 
601      */
602     virtual ELossFit_t* FindBestFit(const TH1* dist,
603                                     Double_t   relErrorCut, 
604                                     Double_t   chi2nuCut,
605                                     Double_t   minWeightCut) const;
606     /** 
607      * Calculate residuals of the fit 
608      * 
609      * @param mode   How to calculate 
610      * @param lowCut Lower cut 
611      * @param dist   Distribution 
612      * @param fit    Function fitted to distribution
613      * @param out    Output list to store residual histogram in
614      */
615     virtual void CalculateResiduals(EResidualMethod  mode,
616                                     Double_t         lowCut,
617                                     TH1*             dist, 
618                                     ELossFit_t*      fit, 
619                                     TCollection*     out) const;
620     /** 
621      * Find the best fits.  This assumes that the array fBest has been
622      * filled with the best possible fits for each eta bin, and that
623      * the fits are placed according to the bin number of the eta bin.
624      *
625      * This is called by the parent class when generating the corretion 
626      * object. 
627      * 
628      * @param d    Parent list
629      * @param obj  Object to add fits to
630      * @param eta  Eta axis 
631      */
632     virtual void FindBestFits(const TList*        d, 
633                               AliFMDCorrELossFit& obj,
634                               const TAxis&        eta);
635     /** 
636      * Make a parameter histogram
637      * 
638      * @param name   Name of histogram.
639      * @param title  Title of histogram. 
640      * @param eta    Eta axis 
641      * 
642      * @return 
643      */
644     TH1* MakePar(const char* name, const char* title, const TAxis& eta) const;
645     /** 
646      * Make a histogram that contains the results of the fit over the
647      * full ring
648      * 
649      * @param name  Name 
650      * @param title Title
651      * @param eta   Eta axis 
652      * @param low   Least bin
653      * @param high  Largest bin
654      * @param val   Value of parameter 
655      * @param err   Error on parameter 
656      * 
657      * @return The newly allocated histogram 
658      */
659     TH1* MakeTotal(const char*  name, 
660                     const char*  title, 
661                     const TAxis& eta, 
662                     Int_t        low, 
663                     Int_t        high, 
664                     Double_t     val, 
665                     Double_t     err) const;
666     TH1*                 fEDist;     // Ring energy distribution 
667     TH1*                 fEmpty;     // Ring energy dist for empty events
668     TH2*                 fHist;      // Two dimension Delta distribution
669     // TList*               fEtaEDists; // Energy distributions per eta bin. 
670     TList*               fList;
671     mutable TObjArray    fBest;
672     mutable TClonesArray fFits;
673     Int_t                fDebug;
674     ClassDef(RingHistos,4);
675   };
676   virtual RingHistos* CreateRingHistos(UShort_t d, Char_t r) const;
677   /** 
678    * Get the ring histogram container 
679    * 
680    * @param d Detector
681    * @param r Ring 
682    * 
683    * @return Ring histogram container 
684    */
685   RingHistos* GetRingHistos(UShort_t d, Char_t r) const;
686   /** 
687    * Check if the detector @a d, ring @a r is listed <i>in</i> the @a
688    * skips bit mask.  If the detector/ring is in the mask, return true.
689    * 
690    * That is, use case is 
691    * @code 
692    *  for (UShort_t d=1. d<=3, d++) {
693    *    UShort_t nr = (d == 1 ? 1 : 2);
694    *    for (UShort_t q = 0; q < nr; q++) { 
695    *      Char_t r = (q == 0 ? 'I' : 'O');
696    *      if (CheckSkips(d, r, skips)) continue; 
697    *      // Process detector/ring 
698    *    }
699    *  }
700    * @endcode
701    *
702    * @param d      Detector
703    * @param r      Ring 
704    * @param skips  Mask of detector/rings to skip
705    * 
706    * @return True if detector @a d, ring @a r is in the mask @a skips 
707    */
708   static Bool_t CheckSkip(UShort_t d, Char_t r, UShort_t skips);
709
710   TList           fRingHistos;        // List of histogram containers
711   Double_t        fLowCut;            // Low cut on energy
712   UShort_t        fNParticles;        // Number of landaus to try to fit 
713   UShort_t        fMinEntries;        // Minimum number of entries
714   UShort_t        fFitRangeBinWidth;  // N-bins to subtract from found max
715   Bool_t          fDoFits;            // Whether to actually do the fits 
716   Bool_t          fDoMakeObject;      // Whether to make corrections object
717   TAxis           fEtaAxis;           // Eta axis 
718   TAxis           fCentralityAxis;    // Centrality axis 
719   Double_t        fMaxE;              // Maximum energy loss to consider 
720   Int_t           fNEbins;            // Number of energy loss bins 
721   Bool_t          fUseIncreasingBins; // Wheter to use increasing bin sizes 
722   Double_t        fMaxRelParError;    // Relative error cut
723   Double_t        fMaxChi2PerNDF;     // chi^2/nu cit
724   Double_t        fMinWeight;         // Minimum weight value 
725   Int_t           fDebug;             // Debug level 
726   EResidualMethod fResidualMethod;    // Whether to store residuals (debugging)
727   UShort_t        fSkips;             // Rings to skip when fitting 
728   Double_t        fRegularizationCut; // When to regularize the chi^2
729
730   ClassDef(AliFMDEnergyFitter,6); //
731 };
732
733 #endif
734 // Local Variables:
735 //  mode: C++ 
736 // End: