AliRoot Core  edcc906 (edcc906)
AliESDFMD.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 /* $Id$ */
17 
18 //____________________________________________________________________
19 //
20 // ESD information from the FMD
21 // Contains information on:
22 // Charged particle multiplicty per strip (rough estimate)
23 // Psuedo-rapdity per strip
24 // Latest changes by Christian Holm Christensen
25 //
26 #include "AliESDFMD.h" // ALIFMDESD_H
27 #include "AliLog.h" // ALILOG_H
28 #include "Riostream.h" // ROOT_Riostream
29 #include <TMath.h>
30 
31 //____________________________________________________________________
32 ClassImp(AliESDFMD)
33 #if 0
34  ; // This is here to keep Emacs for indenting the next line
35 #endif
36 
37 
38 //____________________________________________________________________
39 namespace {
40  // Private implementation of a AliFMDMap::ForOne to use in
41  // forwarding to AliESDFMD::ForOne
42  class ForMultiplicity : public AliFMDMap::ForOne
43  {
44  public:
45  ForMultiplicity(const AliESDFMD& o, AliESDFMD::ForOne& a)
46  : fObject(o), fAlgo(a)
47  {}
48  Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
49  Float_t m)
50  {
51  Float_t e = fObject.Eta(d, r, 0, t);
52  return fAlgo.operator()(d, r, s, t, m, e);
53  }
54  Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Int_t)
55  {
56  return kTRUE;
57  }
58  Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, UShort_t)
59  {
60  return kTRUE;
61  }
62  Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Bool_t)
63  {
64  return kTRUE;
65  }
66  protected:
67  const AliESDFMD& fObject;
68  AliESDFMD::ForOne& fAlgo;
69  };
70 
71  // Private implementation of AliESDFMD::ForOne to print an
72  // object
73  class Printer : public AliESDFMD::ForOne
74  {
75  public:
76  Printer() : fOldD(0), fOldR('-'), fOldS(1024) {}
77  Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
78  Float_t m, Float_t e)
79  {
80  if (d != fOldD) {
81  if (fOldD != 0) printf("\n");
82  fOldD = d;
83  fOldR = '-';
84  printf("FMD%d", fOldD);
85  }
86  if (r != fOldR) {
87  fOldR = r;
88  fOldS = 1024;
89  printf("\n %s ring", (r == 'I' ? "Inner" : "Outer"));
90  }
91  if (s != fOldS) {
92  fOldS = s;
93  printf("\n Sector %d", fOldS);
94  }
95  if (t % 4 == 0) printf("\n %3d-%3d ", t, t+3);
96  if (m == AliESDFMD::kInvalidMult) printf("------/");
97  else printf("%6.3f/", m);
98  if (e == AliESDFMD::kInvalidEta) printf("------ ");
99  else printf("%6.3f ", e);
100 
101  return kTRUE;
102  }
103  private:
104  UShort_t fOldD;
105  Char_t fOldR;
106  UShort_t fOldS;
107  };
108 }
109 
110 //____________________________________________________________________
112  : fMultiplicity(0, 0, 0, 0),
113  fEta(AliFMDFloatMap::kMaxDetectors,
114  AliFMDFloatMap::kMaxRings,
115  1,
116  AliFMDFloatMap::kMaxStrips),
117  fNoiseFactor(0),
118  fAngleCorrected(kFALSE)
119 {
120  // Default CTOR
121 }
122 
123 //____________________________________________________________________
125  : TObject(other),
127  fEta(other.fEta),
128  fNoiseFactor(other.fNoiseFactor),
130 {
131  // Default CTOR
132 }
133 
134 //____________________________________________________________________
135 AliESDFMD&
137 {
138  // Default CTOR
139  if(this == &other) return *this;
140 
141  TObject::operator=(other);
143  fEta = other.fEta;
144 
145  // These two lines were missing prior to version 4 of this class
146  fNoiseFactor = other.fNoiseFactor;
148 
149  return *this;
150 }
151 
152 //____________________________________________________________________
153 void
154 AliESDFMD::Copy(TObject &obj) const
155 {
156  // this overwrites the virtual TOBject::Copy()
157  // to allow run time copying without casting
158  // in AliESDEvent
159 
160  if(this==&obj)return;
161  AliESDFMD *robj = dynamic_cast<AliESDFMD*>(&obj);
162  if(!robj)return; // not an AliESDFMD
163  *robj = *this;
164 }
165 
166 //____________________________________________________________________
167 void
169 {
171  fEta.CheckNeedUShort(file);
172 }
173 
174 //____________________________________________________________________
175 void
176 AliESDFMD::Clear(Option_t* )
177 {
180 }
181 
182 
183 //____________________________________________________________________
184 Float_t
185 AliESDFMD::Multiplicity(UShort_t detector, Char_t ring, UShort_t sector,
186  UShort_t strip) const
187 {
188  // Return rough estimate of charged particle multiplicity in the
189  // strip FMD<detector><ring>[<sector>,<strip>].
190  //
191  // Note, that this should at most be interpreted as the sum
192  // multiplicity of secondaries and primaries.
193  return fMultiplicity(detector, ring, sector, strip);
194 }
195 
196 //____________________________________________________________________
197 Float_t
198 AliESDFMD::Eta(UShort_t detector, Char_t ring, UShort_t /* sector */,
199  UShort_t strip) const
200 {
201  // Return pseudo-rapidity of the strip
202  // FMD<detector><ring>[<sector>,<strip>]. (actually, the sector
203  // argument is ignored, as it is assumed that the primary vertex is
204  // a (x,y) = (0,0), and that the modules are aligned with a
205  // precision better than 2 degrees in the azimuthal angle).
206  //
207  return fEta(detector, ring, 0, strip);
208 }
209 
210 //____________________________________________________________________
211 Float_t
212 AliESDFMD::Phi(UShort_t detector, Char_t ring, UShort_t sector, UShort_t) const
213 {
214  // Return azimuthal angle (in degrees) of the strip
215  // FMD<detector><ring>[<sector>,<strip>].
216  //
217  Float_t baseAng = (detector == 1 ? 90 :
218  detector == 2 ? 0 : 180);
219  Float_t dAng = ((detector == 3 ? -1 : 1) * 360 /
220  (ring == 'I' || ring == 'i' ?
223  Float_t ret = baseAng + dAng * (sector + .5);
224  if (ret > 360) ret -= 360;
225  if (ret < 0) ret += 360;
226  return ret;
227 
228 }
229 
230 //____________________________________________________________________
231 Float_t
232 AliESDFMD::Theta(UShort_t detector, Char_t ring, UShort_t, UShort_t strip) const
233 {
234  // Return polar angle from beam line (in degrees) of the strip
235  // FMD<detector><ring>[<sector>,<strip>].
236  //
237  // This value is calculated from eta and therefor takes into account
238  // the Z position of the interaction point.
239  Float_t eta = Eta(detector, ring, 0, strip);
240  Float_t theta = TMath::ATan(2 * TMath::Exp(-eta));
241  if (theta < 0) theta += TMath::Pi();
242  theta *= 180. / TMath::Pi();
243  return theta;
244 }
245 
246 //____________________________________________________________________
247 Float_t
248 AliESDFMD::R(UShort_t, Char_t ring, UShort_t, UShort_t strip) const
249 {
250  // Return radial distance from beam line (in cm) of the strip
251  // FMD<detector><ring>[<sector>,<strip>].
252  //
253 
254  // Numbers are from AliFMDRing
255  Float_t lR = (ring == 'I' || ring == 'i' ? 4.522 : 15.4);
256  Float_t hR = (ring == 'I' || ring == 'i' ? 17.2 : 28.0);
257  UShort_t nS = (ring == 'I' || ring == 'i' ?
260  Float_t dR = (hR - lR) / nS;
261  Float_t ret = lR + dR * (strip + .5);
262  return ret;
263 
264 }
265 
266 //____________________________________________________________________
267 void
268 AliESDFMD::SetMultiplicity(UShort_t detector, Char_t ring, UShort_t sector,
269  UShort_t strip, Float_t mult)
270 {
271  // Return rough estimate of charged particle multiplicity in the
272  // strip FMD<detector><ring>[<sector>,<strip>].
273  //
274  // Note, that this should at most be interpreted as the sum
275  // multiplicity of secondaries and primaries.
276  fMultiplicity(detector, ring, sector, strip) = mult;
277 }
278 
279 //____________________________________________________________________
280 void
281 AliESDFMD::SetEta(UShort_t detector, Char_t ring, UShort_t /* sector */,
282  UShort_t strip, Float_t eta)
283 {
284  // Set pseudo-rapidity of the strip
285  // FMD<detector><ring>[<sector>,<strip>]. (actually, the sector
286  // argument is ignored, as it is assumed that the primary vertex is
287  // a (x,y) = (0,0), and that the modules are aligned with a
288  // precision better than 2 degrees in the azimuthal angle).
289  //
290  fEta(detector, ring, 0, strip) = eta;
291 }
292 
293 //____________________________________________________________________
294 Bool_t
296 {
297  ForMultiplicity i(*this, a);
298  return fMultiplicity.ForEach(i);
299 }
300 //____________________________________________________________________
301 void
302 AliESDFMD::Print(Option_t* /* option*/) const
303 {
304  // Print all information to standard output.
305  std::cout << "AliESDFMD:" << std::endl;
306  Printer p;
307  ForEach(p);
308  printf("\n");
309 #if 0
310  for (UShort_t det = 1; det <= fMultiplicity.MaxDetectors(); det++) {
311  for (UShort_t ir = 0; ir < fMultiplicity.MaxRings(); ir++) {
312  Char_t ring = (ir == 0 ? 'I' : 'O');
313  std::cout << "FMD" << det << ring << ":" << std::endl;
314  for (UShort_t sec = 0; sec < fMultiplicity.MaxSectors(); sec++) {
315  std::cout << " Sector # " << sec << ":" << std::flush;
316  for (UShort_t str = 0; str < fMultiplicity.MaxStrips(); str++) {
317  if (str % 6 == 0) std::cout << "\n " << std::flush;
318  Float_t m = fMultiplicity(det, ring, sec, str);
319  Float_t e = fEta(det, ring, 0, str);
320  if (m == kInvalidMult && e == kInvalidEta) break;
321  if (m == kInvalidMult) std::cout << " ---- ";
322  else std::cout << Form("%6.3f", m);
323  std::cout << "/";
324  if (e == kInvalidEta) std::cout << " ---- ";
325  else std::cout << Form("%-6.3f", e);
326  std::cout << " " << std::flush;
327  }
328  std::cout << std::endl;
329  }
330  }
331  }
332 #endif
333 }
334 
335 
336 //____________________________________________________________________
337 //
338 // EOF
339 //
UShort_t MaxStrips() const
Definition: AliFMDMap.h:226
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
virtual void Copy(TObject &obj) const
Definition: AliESDFMD.cxx:154
Float_t Phi(UShort_t detector, Char_t ring, UShort_t sector, UShort_t strip) const
Definition: AliESDFMD.cxx:212
Float_t p[]
Definition: kNNTest.C:133
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
void CheckNeedUShort(TFile *file)
Definition: AliFMDMap.cxx:73
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
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 MaxDetectors() const
Definition: AliFMDMap.h:220
virtual Bool_t ForEach(ForOne &algo) const
Definition: AliFMDMap.cxx:335
AliFMDFloatMap fEta
Definition: AliESDFMD.h:286
virtual void Reset(const Float_t &v=Float_t())
Event Summary Data for the Forward Multiplicity Detector.This stores the psuedo-multiplicity and -rap...
Definition: AliESDFMD.h:30
virtual Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t, Float_t v)
Definition: AliFMDMap.h:591
void Clear(Option_t *option="")
Definition: AliESDFMD.cxx:176
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