]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDFMD.cxx
Exchange information only through a data member
[u/mrichter/AliRoot.git] / STEER / AliESDFMD.cxx
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
30 //____________________________________________________________________
31 ClassImp(AliESDFMD)
32 #if 0
33   ; // This is here to keep Emacs for indenting the next line
34 #endif
35
36
37 //____________________________________________________________________
38 AliESDFMD::AliESDFMD()
39   : fMultiplicity(AliFMDFloatMap::kMaxDetectors, 
40                   AliFMDFloatMap::kMaxRings, 
41                   AliFMDFloatMap::kMaxSectors, 
42                   AliFMDFloatMap::kMaxStrips), 
43     fEta(AliFMDFloatMap::kMaxDetectors, 
44          AliFMDFloatMap::kMaxRings, 
45          1,
46          AliFMDFloatMap::kMaxStrips), 
47     fNoiseFactor(0),
48     fAngleCorrected(kFALSE)
49 {
50   // Default CTOR
51 }
52   
53 //____________________________________________________________________
54 AliESDFMD::AliESDFMD(const AliESDFMD& other)
55   : TObject(other), 
56     fMultiplicity(other.fMultiplicity),
57     fEta(other.fEta),
58     fNoiseFactor(other.fNoiseFactor),
59     fAngleCorrected(other.fAngleCorrected)
60 {
61   // Default CTOR
62 }
63
64 //____________________________________________________________________
65 AliESDFMD& 
66 AliESDFMD::operator=(const AliESDFMD& other)
67 {
68   // Default CTOR
69   if(this!=&other){
70     TObject::operator=(other);
71     fMultiplicity = other.fMultiplicity;
72     fEta          = other.fEta;
73   }
74   return *this;
75 }
76
77 void AliESDFMD::Copy(TObject &obj) const{
78   // this overwrites the virtual TOBject::Copy()
79   // to allow run time copying without casting
80   // in AliESDEvent
81
82   if(this==&obj)return;
83   AliESDFMD *robj = dynamic_cast<AliESDFMD*>(&obj);
84   if(!robj)return; // not an AliESDFMD
85   *robj = *this;
86 }
87
88 //____________________________________________________________________
89 void
90 AliESDFMD::CheckNeedUShort(TFile* file) 
91 {
92   fMultiplicity.CheckNeedUShort(file);
93   fEta.CheckNeedUShort(file);
94 }
95
96 //____________________________________________________________________
97 void
98 AliESDFMD::Clear(Option_t* )
99 {
100   fMultiplicity.Reset(kInvalidMult);
101   fEta.Reset(kInvalidEta);
102 }
103
104
105 //____________________________________________________________________
106 Float_t
107 AliESDFMD::Multiplicity(UShort_t detector, Char_t ring, UShort_t sector, 
108                         UShort_t strip) const
109 {
110   // Return rough estimate of charged particle multiplicity in the
111   // strip FMD<detector><ring>[<sector>,<strip>]. 
112   // 
113   // Note, that this should at most be interpreted as the sum
114   // multiplicity of secondaries and primaries. 
115   return fMultiplicity(detector, ring, sector, strip);
116 }
117
118 //____________________________________________________________________
119 Float_t
120 AliESDFMD::Eta(UShort_t detector, Char_t ring, UShort_t /* sector */, 
121                UShort_t strip) const
122 {
123   // Return pseudo-rapidity of the strip
124   // FMD<detector><ring>[<sector>,<strip>].  (actually, the sector
125   // argument is ignored, as it is assumed that the primary vertex is
126   // a (x,y) = (0,0), and that the modules are aligned with a
127   // precision better than 2 degrees in the azimuthal angle). 
128   // 
129   return fEta(detector, ring, 0, strip);
130 }
131
132 //____________________________________________________________________
133 void
134 AliESDFMD::SetMultiplicity(UShort_t detector, Char_t ring, UShort_t sector, 
135                            UShort_t strip, Float_t mult)
136 {
137   // Return rough estimate of charged particle multiplicity in the
138   // strip FMD<detector><ring>[<sector>,<strip>]. 
139   // 
140   // Note, that this should at most be interpreted as the sum
141   // multiplicity of secondaries and primaries. 
142   fMultiplicity(detector, ring, sector, strip) = mult;
143 }
144
145 //____________________________________________________________________
146 void
147 AliESDFMD::SetEta(UShort_t detector, Char_t ring, UShort_t /* sector */, 
148                   UShort_t strip, Float_t eta)
149 {
150   // Set pseudo-rapidity of the strip
151   // FMD<detector><ring>[<sector>,<strip>].  (actually, the sector
152   // argument is ignored, as it is assumed that the primary vertex is
153   // a (x,y) = (0,0), and that the modules are aligned with a
154   // precision better than 2 degrees in the azimuthal angle). 
155   // 
156   fEta(detector, ring, 0, strip) = eta;
157 }
158
159 //____________________________________________________________________
160 void
161 AliESDFMD::Print(Option_t* /* option*/) const
162 {
163   // Print all information to standard output. 
164   std::cout << "AliESDFMD:" << std::endl;
165   for (UShort_t det = 1; det <= fMultiplicity.MaxDetectors(); det++) {
166     for (UShort_t ir = 0; ir < fMultiplicity.MaxRings(); ir++) {
167       Char_t ring = (ir == 0 ? 'I' : 'O');
168       std::cout << "FMD" << det << ring << ":" << std::endl;
169       for  (UShort_t sec = 0; sec < fMultiplicity.MaxSectors(); sec++) {
170         std::cout << " Sector # " << sec << ":" << std::flush;
171         for (UShort_t str = 0; str < fMultiplicity.MaxStrips(); str++) {
172           if (str % 6 == 0) std::cout << "\n  " << std::flush;
173           Float_t m = fMultiplicity(det, ring, sec, str);
174           Float_t e = fEta(det, ring, 0, str);
175           if (m == kInvalidMult && e == kInvalidEta) break;
176           if (m == kInvalidMult) std::cout << " ---- ";
177           else                   std::cout << Form("%6.3f", m);
178           std::cout << "/";
179           if (e == kInvalidEta)  std::cout << " ---- ";
180           else                   std::cout << Form("%-6.3f", e);
181           std::cout << " " << std::flush;
182         }
183         std::cout << std::endl;
184       }
185     }
186   }
187 }
188
189
190 //____________________________________________________________________
191 //
192 // EOF
193 //