1 /**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //____________________________________________________________________
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
26 #include "AliESDFMD.h" // ALIFMDESD_H
27 #include "AliLog.h" // ALILOG_H
28 #include "Riostream.h" // ROOT_Riostream
31 //____________________________________________________________________
34 ; // This is here to keep Emacs for indenting the next line
38 //____________________________________________________________________
40 // Private implementation of a AliFMDMap::ForOne to use in
41 // forwarding to AliESDFMD::ForOne
42 class ForMultiplicity : public AliFMDMap::ForOne
45 ForMultiplicity(const AliESDFMD& o, AliESDFMD::ForOne& a)
46 : fObject(o), fAlgo(a)
48 Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
51 Float_t e = fObject.Eta(d, r, 0, t);
52 return fAlgo.operator()(d, r, s, t, m, e);
54 Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Int_t)
58 Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, UShort_t)
62 Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Bool_t)
67 const AliESDFMD& fObject;
68 AliESDFMD::ForOne& fAlgo;
71 // Private implementation of AliESDFMD::ForOne to print an
73 class Printer : public AliESDFMD::ForOne
76 Printer() : fOldD(0), fOldR('-'), fOldS(1024) {}
77 Bool_t operator()(UShort_t d, Char_t r, UShort_t s, UShort_t t,
81 if (fOldD != 0) printf("\n");
84 printf("FMD%d", fOldD);
89 printf("\n %s ring", (r == 'I' ? "Inner" : "Outer"));
93 printf("\n Sector %d", fOldS);
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);
110 //____________________________________________________________________
111 AliESDFMD::AliESDFMD()
112 : fMultiplicity(0, 0, 0, 0),
113 fEta(AliFMDFloatMap::kMaxDetectors,
114 AliFMDFloatMap::kMaxRings,
116 AliFMDFloatMap::kMaxStrips),
118 fAngleCorrected(kFALSE)
123 //____________________________________________________________________
124 AliESDFMD::AliESDFMD(const AliESDFMD& other)
126 fMultiplicity(other.fMultiplicity),
128 fNoiseFactor(other.fNoiseFactor),
129 fAngleCorrected(other.fAngleCorrected)
134 //____________________________________________________________________
136 AliESDFMD::operator=(const AliESDFMD& other)
140 TObject::operator=(other);
141 fMultiplicity = other.fMultiplicity;
147 //____________________________________________________________________
149 AliESDFMD::Copy(TObject &obj) const
151 // this overwrites the virtual TOBject::Copy()
152 // to allow run time copying without casting
155 if(this==&obj)return;
156 AliESDFMD *robj = dynamic_cast<AliESDFMD*>(&obj);
157 if(!robj)return; // not an AliESDFMD
161 //____________________________________________________________________
163 AliESDFMD::CheckNeedUShort(TFile* file)
165 fMultiplicity.CheckNeedUShort(file);
166 fEta.CheckNeedUShort(file);
169 //____________________________________________________________________
171 AliESDFMD::Clear(Option_t* )
173 fMultiplicity.Reset(kInvalidMult);
174 fEta.Reset(kInvalidEta);
178 //____________________________________________________________________
180 AliESDFMD::Multiplicity(UShort_t detector, Char_t ring, UShort_t sector,
181 UShort_t strip) const
183 // Return rough estimate of charged particle multiplicity in the
184 // strip FMD<detector><ring>[<sector>,<strip>].
186 // Note, that this should at most be interpreted as the sum
187 // multiplicity of secondaries and primaries.
188 return fMultiplicity(detector, ring, sector, strip);
191 //____________________________________________________________________
193 AliESDFMD::Eta(UShort_t detector, Char_t ring, UShort_t /* sector */,
194 UShort_t strip) const
196 // Return pseudo-rapidity of the strip
197 // FMD<detector><ring>[<sector>,<strip>]. (actually, the sector
198 // argument is ignored, as it is assumed that the primary vertex is
199 // a (x,y) = (0,0), and that the modules are aligned with a
200 // precision better than 2 degrees in the azimuthal angle).
202 return fEta(detector, ring, 0, strip);
205 //____________________________________________________________________
207 AliESDFMD::Phi(UShort_t detector, Char_t ring, UShort_t sector, UShort_t) const
209 // Return azimuthal angle (in degrees) of the strip
210 // FMD<detector><ring>[<sector>,<strip>].
212 Float_t baseAng = (detector == 1 ? 90 :
213 detector == 2 ? 0 : 180);
214 Float_t dAng = ((detector == 3 ? -1 : 1) * 360 /
215 (ring == 'I' || ring == 'i' ?
216 AliFMDMap::kNSectorInner :
217 AliFMDMap::kNSectorOuter));
218 Float_t ret = baseAng + dAng * (sector + .5);
219 if (ret > 360) ret -= 360;
220 if (ret < 0) ret += 360;
225 //____________________________________________________________________
227 AliESDFMD::Theta(UShort_t detector, Char_t ring, UShort_t, UShort_t strip) const
229 // Return polar angle from beam line (in degrees) of the strip
230 // FMD<detector><ring>[<sector>,<strip>].
232 // This value is calculated from eta and therefor takes into account
233 // the Z position of the interaction point.
234 Float_t eta = Eta(detector, ring, 0, strip);
235 Float_t theta = TMath::ATan(2 * TMath::Exp(-eta));
236 if (theta < 0) theta += TMath::Pi();
237 theta *= 180. / TMath::Pi();
241 //____________________________________________________________________
243 AliESDFMD::R(UShort_t, Char_t ring, UShort_t, UShort_t strip) const
245 // Return radial distance from beam line (in cm) of the strip
246 // FMD<detector><ring>[<sector>,<strip>].
249 // Numbers are from AliFMDRing
250 Float_t lR = (ring == 'I' || ring == 'i' ? 4.522 : 15.4);
251 Float_t hR = (ring == 'I' || ring == 'i' ? 17.2 : 28.0);
252 UShort_t nS = (ring == 'I' || ring == 'i' ?
253 AliFMDMap::kNStripInner :
254 AliFMDMap::kNStripOuter);
255 Float_t dR = (hR - lR) / nS;
256 Float_t ret = lR + dR * (strip + .5);
261 //____________________________________________________________________
263 AliESDFMD::SetMultiplicity(UShort_t detector, Char_t ring, UShort_t sector,
264 UShort_t strip, Float_t mult)
266 // Return rough estimate of charged particle multiplicity in the
267 // strip FMD<detector><ring>[<sector>,<strip>].
269 // Note, that this should at most be interpreted as the sum
270 // multiplicity of secondaries and primaries.
271 fMultiplicity(detector, ring, sector, strip) = mult;
274 //____________________________________________________________________
276 AliESDFMD::SetEta(UShort_t detector, Char_t ring, UShort_t /* sector */,
277 UShort_t strip, Float_t eta)
279 // Set pseudo-rapidity of the strip
280 // FMD<detector><ring>[<sector>,<strip>]. (actually, the sector
281 // argument is ignored, as it is assumed that the primary vertex is
282 // a (x,y) = (0,0), and that the modules are aligned with a
283 // precision better than 2 degrees in the azimuthal angle).
285 fEta(detector, ring, 0, strip) = eta;
288 //____________________________________________________________________
290 AliESDFMD::ForEach(AliESDFMD::ForOne& a) const
292 ForMultiplicity i(*this, a);
293 return fMultiplicity.ForEach(i);
296 //____________________________________________________________________
298 AliESDFMD::Print(Option_t* /* option*/) const
300 // Print all information to standard output.
301 std::cout << "AliESDFMD:" << std::endl;
306 for (UShort_t det = 1; det <= fMultiplicity.MaxDetectors(); det++) {
307 for (UShort_t ir = 0; ir < fMultiplicity.MaxRings(); ir++) {
308 Char_t ring = (ir == 0 ? 'I' : 'O');
309 std::cout << "FMD" << det << ring << ":" << std::endl;
310 for (UShort_t sec = 0; sec < fMultiplicity.MaxSectors(); sec++) {
311 std::cout << " Sector # " << sec << ":" << std::flush;
312 for (UShort_t str = 0; str < fMultiplicity.MaxStrips(); str++) {
313 if (str % 6 == 0) std::cout << "\n " << std::flush;
314 Float_t m = fMultiplicity(det, ring, sec, str);
315 Float_t e = fEta(det, ring, 0, str);
316 if (m == kInvalidMult && e == kInvalidEta) break;
317 if (m == kInvalidMult) std::cout << " ---- ";
318 else std::cout << Form("%6.3f", m);
320 if (e == kInvalidEta) std::cout << " ---- ";
321 else std::cout << Form("%-6.3f", e);
322 std::cout << " " << std::flush;
324 std::cout << std::endl;
332 //____________________________________________________________________