AliRoot Core  3dc7879 (3dc7879)
AliESDFMD.h
Go to the documentation of this file.
1 #ifndef ALIESDFMD_H
2 #define ALIESDFMD_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved.
5  *
6  * See cxx source for full Copyright notice
7  */
8 //___________________________________________________________________
9 //
10 // AliESDFMD is the Event Summary Data entry for the FMD. It contains
11 // a rough estimate of the charged particle multiplicity in each strip
12 // of the FMD. It also contains the psuedo-rapidity of each strip.
13 // This is important, as it varies from event to event, due to a
14 // finite interaction point probability distribution.
15 //
16 #ifndef ROOT_TObject
17 # include <TObject.h>
18 #endif
19 #ifndef ALIFMDFLOATMAP_H
20 # include "AliFMDFloatMap.h"
21 #endif
22 
23 //___________________________________________________________________
30 class AliESDFMD : public TObject
31 {
32 public:
50  class ForOne
51  {
52  public:
56  virtual ~ForOne() {}
70  virtual bool operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
71  Float_t m, Float_t e) = 0;
72  };
73  enum {
74  kNeedNoiseFix = (1 << 14)
75  };
79  AliESDFMD();
85  AliESDFMD(const AliESDFMD& other);
93  AliESDFMD& operator=(const AliESDFMD& other);
97  virtual ~AliESDFMD() {}
104  virtual void Copy(TObject &obj) const;
105 
109  void Clear(Option_t *option="");
122  Float_t Multiplicity(UShort_t detector, Char_t ring,
123  UShort_t sector, UShort_t strip) const;
136  Float_t Eta(UShort_t detector, Char_t ring,
137  UShort_t sector, UShort_t strip) const;
150  Float_t Phi(UShort_t detector, Char_t ring,
151  UShort_t sector, UShort_t strip) const;
164  Float_t Theta(UShort_t detector, Char_t ring,
165  UShort_t sector, UShort_t strip) const;
178  Float_t R(UShort_t detector, Char_t ring,
179  UShort_t sector, UShort_t strip) const;
191  void SetMultiplicity(UShort_t detector, Char_t ring,
192  UShort_t sector, UShort_t strip,
193  Float_t mult);
205  void SetEta(UShort_t detector, Char_t ring,
206  UShort_t sector, UShort_t strip,
207  Float_t eta);
211  void SetNoiseFactor(Float_t f) { fNoiseFactor = f; }
215  void SetAngleCorrected(Bool_t done) { fAngleCorrected = done; }
219  Bool_t IsAngleCorrected() const { return fAngleCorrected; }
223  Float_t GetNoiseFactor() const { return fNoiseFactor; }
227  UShort_t MaxDetectors() const { return fMultiplicity.MaxDetectors(); }
231  UShort_t MaxRings() const { return fMultiplicity.MaxRings(); }
235  UShort_t MaxSectors() const { return fMultiplicity.MaxSectors(); }
239  UShort_t MaxStrips() const { return fMultiplicity.MaxStrips(); }
245  void Print(Option_t* option="") const;
251  void CheckNeedUShort(TFile* file);
257  Bool_t NeedNoiseFix() const { return TestBit(kNeedNoiseFix); }
267  Bool_t ForEach(ForOne& algo) const;
268  enum {
271  };
272  enum {
274  kInvalidEta = 1024
275  };
279  const AliFMDFloatMap& MultiplicityMap() const { return fMultiplicity; }
283  const AliFMDFloatMap& EtaMap() const { return fEta; }
284 protected:
285  AliFMDFloatMap fMultiplicity; // Psuedo multplicity per strip
286  AliFMDFloatMap fEta; // Psuedo-rapidity per strip
287  Float_t fNoiseFactor; // Factor used for noise suppresion
288  Bool_t fAngleCorrected; // Whether we've done angle correction
289  ClassDef(AliESDFMD,4) // ESD info from FMD
290 };
291 #endif
292 //____________________________________________________________________
293 //
294 // Local Variables:
295 // mode: C++
296 // End:
297 //
298 // EOF
299 //
UShort_t MaxStrips() const
Definition: AliFMDMap.h:226
Bool_t NeedNoiseFix() const
Definition: AliESDFMD.h:257
void SetNoiseFactor(Float_t f)
Definition: AliESDFMD.h:211
virtual void Copy(TObject &obj) const
Definition: AliESDFMD.cxx:154
virtual ~AliESDFMD()
Definition: AliESDFMD.h:97
Float_t Phi(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:212
UShort_t MaxRings() const
Definition: AliESDFMD.h:231
Float_t GetNoiseFactor() const
Definition: AliESDFMD.h:223
Bool_t IsAngleCorrected() const
Definition: AliESDFMD.h:219
Float_t R(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:248
Bool_t ForEach(ForOne &algo) const
Definition: AliESDFMD.cxx:295
UShort_t MaxStrips() const
Definition: AliESDFMD.h:239
Float_t fNoiseFactor
Definition: AliESDFMD.h:287
AliESDFMD & operator=(const AliESDFMD &other)
Definition: AliESDFMD.cxx:136
void SetMultiplicity(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip, Float_t mult)
Definition: AliESDFMD.cxx:268
const AliFMDFloatMap & MultiplicityMap() const
Definition: AliESDFMD.h:279
void SetEta(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip, Float_t eta)
Definition: AliESDFMD.cxx:281
Float_t Multiplicity(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:185
Float_t Eta(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:198
UShort_t MaxRings() const
Definition: AliFMDMap.h:222
UShort_t MaxSectors() const
Definition: AliESDFMD.h:235
UShort_t MaxDetectors() const
Definition: AliFMDMap.h:220
virtual ~ForOne()
Definition: AliESDFMD.h:56
AliFMDFloatMap fEta
Definition: AliESDFMD.h:286
Event Summary Data for the Forward Multiplicity Detector.This stores the psuedo-multiplicity and -rap...
Definition: AliESDFMD.h:30
TF1 * f
Definition: interpolTest.C:21
void Clear(Option_t *option="")
Definition: AliESDFMD.cxx:176
const AliFMDFloatMap & EtaMap() const
Definition: AliESDFMD.h:283
Float_t Theta(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:232
void CheckNeedUShort(TFile *file)
Definition: AliESDFMD.cxx:168
Bool_t fAngleCorrected
Definition: AliESDFMD.h:288
UShort_t MaxSectors() const
Definition: AliFMDMap.h:224
AliFMDFloatMap fMultiplicity
Definition: AliESDFMD.h:285
void Print(Option_t *option="") const
Definition: AliESDFMD.cxx:302
virtual bool operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t, Float_t m, Float_t e)=0
void SetAngleCorrected(Bool_t done)
Definition: AliESDFMD.h:215
UShort_t MaxDetectors() const
Definition: AliESDFMD.h:227