]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/FORWARD/analysis2/AliBaseMCTrackDensity.h
Maker scripts now use TrainSetup exclusively.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliBaseMCTrackDensity.h
CommitLineData
4bcdcbc1 1#ifndef ALIBaseMCTRACKDENSITY_MC
2#define ALIBaseMCTRACKDENSITY_MC
3#include <TNamed.h>
936b0a6c 4#include "AliForwardFlowWeights.h"
4bcdcbc1 5class TList;
6class TH1D;
7class TH2D;
8class TF1;
9class TGraph;
10class AliMCEvent;
11class AliMCParticle;
12class AliTrackReference;
13class AliStack;
14
15/**
16 * A class to calculate the particle density from track references.
17 * This code is used both in AliForwardMCCorrectionsTask and
18 * AliBaseMCDensity calculator.
19 *
20 * @par Input:
21 * - Kinematics
22 * - Track-References
23 *
24 * @par Corrections used:
25 * - None
26 *
27 * @par Histograms:
28 * - Incident angle vs number of track references
29 * - Incident angle vs number of strips/cluster
30 *
31 * @ingroup pwglf_forward_algo
32 * @ingroup pwglf_forward_mc
33 * @ingroup pwglf_forward_aod
34 */
35class AliBaseMCTrackDensity : public TNamed
36{
37public:
38 /**
39 * Default constructor. Do not use - for ROOT I/O system use only
40 */
41 AliBaseMCTrackDensity();
42 /**
43 * Normal constructor
44 *
45 * @param name Not used
46 */
47 AliBaseMCTrackDensity(const char* name);
48 /**
49 * Copy constructor
50 *
51 * @param o Object to copy from
52 */
53 AliBaseMCTrackDensity(const AliBaseMCTrackDensity& o);
54 /**
55 * Assignment operator
56 *
57 * @param o Object to assign from
58 *
59 * @return Reference to this
60 */
61 AliBaseMCTrackDensity& operator=(const AliBaseMCTrackDensity& o);
62 /**
63 * Destructor.
64 */
65 virtual ~AliBaseMCTrackDensity() {}
66
67 /**
68 * Set whether to only consider primaries
69 *
70 * @param use If true, consider only primaries
71 */
72 void SetUseOnlyPrimary(Bool_t use) { fUseOnlyPrimary = use; }
73 /**
74 * Whether to `after-burn' the particles with flow
75 *
76 * @param use
77 */
78 void SetUseFlowWeights(Bool_t use) { fUseFlowWeights = use; }
79 /**
80 * Set whether to print debug messages. Please note this will
81 * produce a lot of output.
82 *
83 * @param debug Whether to enable debug messages or not
84 */
85 void SetDebug(Bool_t debug=true) { fDebug = debug; }
86 /**
87 * Define ouputs
88 *
89 * @param list List to add outputs to
90 */
91 virtual void DefineOutput(TList* list);
92 /**
93 * Print information to standard out
94 *
95 * @param option Not used
96 */
97 virtual void Print(Option_t* option="") const;
98protected:
4bcdcbc1 99 /**
100 * Loops over all the particles in the passed event. If @a primary
101 * is not null, then that histogram is filled with the primary
102 * particle information - irrespective of whether the particle
103 * actually hits the Base or not. For each track (primary or
104 * secondary, unless only primary information is requested - see
105 * SetUseOnlyPrimary) loop over all track references to that
106 * particle and check if they come from the Base. In that case,
107 * figure out which strip(s) to assign the track to, and fill the @a
108 * hits structure.
109 *
4bcdcbc1 110 * @param event MC event
111 * @param vz IP z-coordinate
4bcdcbc1 112 * @param primary Primary information, if available.
113 *
114 * @return true
115 */
116 Bool_t ProcessTracks(const AliMCEvent& event,
117 Double_t vz,
118 TH2D* primary);
119 /**
120 * Process a single track
121 *
122 * @param particle Particle
123 * @param mother Ultimate mother
124 *
125 * @return true on succsess
126 */
127 Bool_t ProcessTrack(AliMCParticle* particle,
128 const AliMCParticle* mother);
129 /**
130 * Must be defined to return the track-reference ID for this detector
131 *
132 * @return Detector id set on track references
133 */
134 virtual Int_t GetDetectorId() const = 0;
135 /**
136 * Process a track reference
137 *
138 * @param particle Particle
139 * @param mother Ultimate mother (if not primary)
140 * @param ref Reference
141 *
142 * @return 0 if no output should be generated for this reference, or
143 * pointer to track-reference to produce output for.
144 */
145 virtual AliTrackReference* ProcessRef(AliMCParticle* particle,
146 const AliMCParticle* mother,
147 AliTrackReference* ref) = 0;
148 /**
149 * Called at before loop over track references
150 *
151 */
152 virtual void BeginTrackRefs() {}
290052e7 153 /**
154 * Check a track reference
155 *
156 * @return true if the track reference should be used
157 */
158 virtual Bool_t CheckTrackRef(AliTrackReference*) const { return true; }
4bcdcbc1 159 /**
160 * Called at before loop over track references
161 *
162 */
163 virtual void EndTrackRefs(Int_t /*nRefs*/) {}
164 /**
165 * Store a particle hit in Base<i>dr</i>[<i>s,t</i>] in @a output
166 *
167 *
168 * @param particle Particle to store
169 * @param mother Ultimate mother of particle
290052e7 170 * @param ref Longest track reference
171 *
172 * @return Weight factor
4bcdcbc1 173 */
174 virtual Double_t StoreParticle(AliMCParticle* particle,
175 const AliMCParticle* mother,
176 AliTrackReference* ref) const;
177 /**
178 * Get the parameters of this event
179 *
180 * @param event Event
181 *
182 * @return true if found, false otherwise
183 */
184 Bool_t GetCollisionParameters(const AliMCEvent& event);
185 /**
186 * Get incident angle of this track reference
187 *
188 * @param ref Track reference
4bcdcbc1 189 *
190 * @return incident angle (in radians)
191 */
192 Double_t GetTrackRefTheta(const AliTrackReference* ref) const;
193 /**
194 * Get ultimate mother of a track
195 *
196 * @param iTr Track number
197 * @param event Event
198 *
199 * @return Pointer to mother or null
200 */
201 const AliMCParticle* GetMother(Int_t iTr, const AliMCEvent& event) const;
202 /**
203 * Calculate flow weight
204 *
205 * @param eta Pseudo rapidity
206 * @param pt Transverse momemtum
207 * @param phi Azimuthal angle
208 * @param id Particle PDG code
209 *
210 * @return Flow weight for the particle
211 */
212 Double_t CalculateWeight(Double_t eta, Double_t pt,
213 Double_t phi, Int_t id) const;
214
215 Bool_t fUseOnlyPrimary; // Only use primaries
216 Bool_t fUseFlowWeights; // Wether to use flow weights
217 TH2D* fBinFlow; // eta,phi bin flow
218 TH2D* fEtaBinFlow; // dEta vs eta of strip
219 TH2D* fPhiBinFlow; // dPhi vs phi of strip
220 TH1D* fNRefs; // Number of track-references per track
936b0a6c 221 AliForwardFlowWeights fWeights; // Flow weights
4bcdcbc1 222 Double_t fVz; // IP z-coordinate of this event
223 Double_t fB; // Impact parameter of this event
224 Double_t fPhiR; // Reaction plane of this event
225 Bool_t fDebug; // Debug flag
226
936b0a6c 227 ClassDef(AliBaseMCTrackDensity,2); // Calculate track-ref density
4bcdcbc1 228};
229
230#endif
231// Local Variables:
232// mode: C++
233// End: