]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/FORWARD/analysis2/AliFMDMCTrackELoss.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDMCTrackELoss.h
CommitLineData
7a3e34f4 1#ifndef ALIFMDMCTRACKELOSS_MC
2#define ALIFMDMCTRACKELOSS_MC
3#include "AliForwardUtil.h"
4#include "AliBaseMCTrackDensity.h"
5#include "AliFMDFloatMap.h"
6class TH1D;
7class TH2;
8class TH2D;
9class AliESDFMD;
10class TClonesArray;
11class TTree;
12
13/**
14 * A class to calculate the particle eloss from track references.
15 * This code is used both in AliForwardMCCorrectionsTask and
16 * AliFMDMCEloss calculator.
17 *
18 * @par Input:
19 * - AliESDFMD object - from reconstruction
20 * - Kinematics
21 * - Track-References
22 *
23 * @par Output:
24 * - AliESDFMD object - content is # of track references/strip
25 *
26 * @par Corrections used:
27 * - None
28 *
29 * @par Histograms:
30 * - Incident angle vs number of track references
31 * - Incident angle vs number of strips/cluster
32 *
33 * @ingroup pwglf_forward_algo
34 * @ingroup pwglf_forward_mc
35 * @ingroup pwglf_forward_aod
36 */
37class AliFMDMCTrackELoss : public AliBaseMCTrackDensity
38{
39public:
40 /**
41 * Default constructor. Do not use - for ROOT I/O system use only
42 */
43 AliFMDMCTrackELoss();
44 /**
45 * Normal constructor
46 *
47 * @param name Not used
48 */
49 AliFMDMCTrackELoss(const char* name);
50 /**
51 * Copy constructor
52 *
53 * @param o Object to copy from
54 */
ee83b849 55 AliFMDMCTrackELoss(const AliFMDMCTrackELoss& o){;}
7a3e34f4 56 /**
57 * Assignment operator
58 *
59 * @param o Object to assign from
60 *
61 * @return Reference to this
62 */
ee83b849 63 AliFMDMCTrackELoss& operator=(const AliFMDMCTrackELoss& o){return *this;}
7a3e34f4 64 /**
65 * Destructor.
66 */
67 virtual ~AliFMDMCTrackELoss() {}
68
69 /**
70 * @{
71 * @name Options
72 */
73 /**
74 * Whether to make an output nTree
75 *
76 * @param use If true, make an nTree of hits
77 */
78 void SetUseTree(Bool_t use=true) { fUseTree = use; }
79 /**
80 * Set maximum number of strips per 'cluster'
81 *
82 * @param n Maximum number of strips per 'cluster'
83 */
84 void SetMaxConsequtiveStrips(UShort_t n) { fMaxConsequtiveStrips = n; }
85 /* @} */
86
87 /**
88 * Loops over all the particles in the passed event. If @a primary
89 * is not null, then that histogram is filled with the primary
90 * particle information - irrespective of whether the particle
91 * actually hits the FMD or not. For each track (primary or
92 * secondary, unless only primary information is requested - see
93 * SetUseOnlyPrimary) loop over all track references to that
94 * particle and check if they come from the FMD. In that case,
95 * figure out which strip(s) to assign the track to, and fill the @a
96 * hits structure.
97 *
98 * @param esd FMD ESD structure
99 * @param event MC event
100 * @param vz IP z-coordinate
101 * @param cent Event centrality
102 *
103 * @return true
104 */
105 Bool_t Calculate(const AliESDFMD& esd,
106 const AliMCEvent& event,
107 Double_t vz,
108 Double_t cent);
109 /**
110 * Define ouputs
111 *
112 * @param list List to add outputs to
113 */
114 void CreateOutputObjects(TList* list);
115
116 /**
117 * @{
118 * @name Access to cache
119 */
120 /**
121 * Clear caches
122 */
123 void Clear(Option_t* opt="");
124 /**
125 * Get reference to cache of eloss per primary
126 *
127 * @return Reference to cache of eloss per primary
128 */
129 AliFMDFloatMap& GetPrimaries() { return fPrimaries; }
130 /**
131 * Get constant reference to cache of eloss per primary
132 *
133 * @return Constant reference to cache of eloss per primary
134 */
135 const AliFMDFloatMap& GetPrimaries() const { return fPrimaries; }
136 /**
137 * Get reference to cache of eloss per secondary
138 *
139 * @return Reference to cache of eloss per secondary
140 */
141 AliFMDFloatMap& GetSecondaries() { return fSecondaries; }
142 /**
143 * Get constant reference to cache of eloss per secondary
144 *
145 * @return Constant reference to cache of eloss per secondary
146 */
147 const AliFMDFloatMap& GetSecondaries() const { return fSecondaries; }
148 /**
149 * Get reference to cache of eloss for all
150 *
151 * @return Reference to cache of eloss for all
152 */
153 AliFMDFloatMap& GetAll() { return fAll; }
154 /**
155 * Get constant reference to cache of eloss for all
156 *
157 * @return Constant reference to cache of eloss for all
158 */
159 const AliFMDFloatMap& GetAll() const { return fAll; }
160 /**
161 * Get reference to cache of psuedo-rapidity (@f$ \eta@f$)
162 *
163 * @return Reference to cache of psuedo-rapidity (@f$ \eta@f$)
164 */
165 AliFMDFloatMap& GetEta() { return fEta; }
166 /**
167 * Get constant reference to cache of psuedo-rapidity (@f$ \eta@f$)
168 *
169 * @return Constant reference to cache of psuedo-rapidity (@f$ \eta@f$)
170 */
171 const AliFMDFloatMap& GetEta() const { return fEta; }
172
173 /**
174 * Get pointer to NTutple if defined
175 *
176 * @return Pointer to nTree or null if not enabled
177 */
178 TTree* GetTree() const { return fTree; }
179 TClonesArray* GetHits() const { return fHits; }
180
181 TH2* GetBetaGammadEdx() const { return fBetaGammadEdx; }
182 TH2* GetBetaGammaEta() const { return fBetaGammaEta; }
183 TH2* GetDEdxEta() const { return fDEdxEta; }
184 /* @} */
185
186 /**
187 * Print this task
188 *
189 * @param option Not used
190 */
191 void Print(Option_t* option="") const;
192
193
194 /**
195 * Structure to hold hit information
196 */
197 struct Hit : public TObject {
198 Double_t fGamma;
199 Double_t fBeta;
200 Double_t fEta;
201 Double_t fDe;
202 Double_t fDx;
203 UInt_t fDetId;
204 Int_t fPdg;
205 Bool_t fPrimary;
206 mutable UShort_t fDetector; //!
207 mutable Char_t fRing; //!
208 mutable UShort_t fSector; //!
209 mutable UShort_t fStrip; //!
210
211 Hit();
212 void Decode() const;
213 UShort_t D() const { Decode(); return fDetector; }
214 Char_t R() const { Decode(); return fRing; }
215 UShort_t S() const { Decode(); return fSector; }
216 UShort_t T() const { Decode(); return fStrip; }
217 Double_t DeDx() const { return (fDx > 0 ? fDe/fDx : 0); }
218 Double_t Eta() const { return fEta; }
219 Double_t BetaGamma() const { return fBeta*fGamma; }
220 UInt_t AbsPdg() const;
221 Bool_t IsPrimary() const { return fPrimary; }
222 Bool_t IsPion() const { return AbsPdg() == 211; }
223 Bool_t IsKaon() const { return AbsPdg() == 321; }
224 Bool_t IsProton() const { return AbsPdg() == 2212; }
225 Bool_t IsElectron() const { return AbsPdg() == 11 || AbsPdg() == 13; }
226
227 ClassDef(Hit,1);
228 };
229 /**
230 * Structure to hold event information
231 */
232 struct Event {
233 Double_t fIpZ;
234 Double_t fCent;
235 };
236
237protected:
238 /**
239 * Must be defined to return the track-reference ID for this detector
240 *
241 * @return Detector id set on track references
242 */
243 Int_t GetDetectorId() const;
244 /**
245 * Process a track reference
246 *
247 * @param particle Particle
248 * @param mother Ultimate mother (if not primary)
249 * @param ref Reference
250 *
251 * @return 0 if no output should be generated for this reference, or
252 * pointer to track-reference to produce output for.
253 */
254 AliTrackReference* ProcessRef(AliMCParticle* particle,
255 const AliMCParticle* mother,
256 AliTrackReference* ref);
257 /**
258 * Called at before loop over track references
259 *
260 */
261 void BeginTrackRefs();
262 /**
263 * Called at before loop over track references
264 *
265 * @param nRefs Number of references
266 */
267 void EndTrackRefs(Int_t nRefs);
268 /**
269 * Store a particle hit in Base<i>dr</i>[<i>s,t</i>] in @a output
270 *
271 *
272 * @param particle Particle to store
273 * @param mother Ultimate mother of particle
274 * @param ref Longest track reference
275 *
276 * @return weight
277 */
278 Double_t StoreParticle(AliMCParticle* particle,
279 const AliMCParticle* mother,
280 AliTrackReference* ref) const;
281 //__________________________________________________________________
282 /**
283 * Structure holding the state of the `tracker'
284 *
285 */
286 mutable struct State
287 {
288 Double_t angle; // Angle
289 UShort_t oldDetector; // Last detector
290 Char_t oldRing; // Last ring
291 UShort_t oldSector; // Last sector
292 UShort_t oldStrip; // Last strip
293 UShort_t startStrip; // First strip
294 UShort_t nRefs; // Number of references
295 UShort_t nStrips; // Number of strips
296 UShort_t count; // Count of hit strips
297 Double_t de; // Summed energy loss
298 Double_t dx; // Summed dx
299 AliTrackReference* longest; //! Longest track through
300 /**
301 * Clear this state
302 *
303 * @param alsoCount If true, also clear count
304 */
305 void Clear(Bool_t alsoCount=false);
306 /**
307 * Assignment operator
308 *
309 * @param o Object to assign from
310 *
311 * @return Reference to this object
312 */
313 State& operator=(const State& o);
314 } fState; //! State
315 Event fEvent;
316
317
318 UShort_t fMaxConsequtiveStrips; // Max 'cluster' size
319 Bool_t fUseTree; // Whether to make an nTree
320 TClonesArray* fHits;
321 TTree* fTree;
322 TH1D* fNr; // Number of track-refs per cluster
323 TH1D* fNt; // Size of cluster in strips
324 TH1D* fNc; // Number of clusters per track
325 TH2D* fNcr; // Number of clusters per track
326 TH2* fBetaGammadEdx;
327 TH2* fBetaGammaEta;
328 TH2* fDEdxEta;
329 mutable AliFMDFloatMap fPrimaries; // Cache of eloss per primary
330 mutable AliFMDFloatMap fSecondaries; // Cache of eloss per secondary
331 mutable AliFMDFloatMap fAll; // Cache of eloss for all
332 AliFMDFloatMap fEta; // Cache of pseudo-rapidity
333
334
335 ClassDef(AliFMDMCTrackELoss,1); // Calculate track-ref density
336};
337
338#endif
339// Local Variables:
340// mode: C++
341// End: