]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliBaseMCTrackDensity.h
Updates to DebugGuard
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliBaseMCTrackDensity.h
1 #ifndef ALIBaseMCTRACKDENSITY_MC
2 #define ALIBaseMCTRACKDENSITY_MC
3 #include <TNamed.h>
4 #include "AliForwardFlowWeights.h"
5 class TList;
6 class TH1D;
7 class TH2D;
8 class TF1;
9 class TGraph;
10 class AliMCEvent;
11 class AliMCParticle;
12 class AliTrackReference;
13 class 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  */
35 class AliBaseMCTrackDensity : public TNamed
36 {
37 public:
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;
98 protected:
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    * 
110    * @param esd      Base ESD structure 
111    * @param event    MC event 
112    * @param vz       IP z-coordinate
113    * @param output   Output of Base hits
114    * @param primary  Primary information, if available. 
115    * 
116    * @return true 
117    */
118   Bool_t ProcessTracks(const AliMCEvent&   event, 
119                        Double_t            vz,
120                        TH2D*               primary);
121   /** 
122    * Process a single track 
123    * 
124    * @param particle   Particle 
125    * @param mother     Ultimate mother
126    * 
127    * @return true on succsess 
128    */
129   Bool_t ProcessTrack(AliMCParticle* particle, 
130                       const AliMCParticle* mother);
131   /** 
132    * Must be defined to return the track-reference ID for this detector
133    * 
134    * @return Detector id set on track references
135    */
136   virtual Int_t GetDetectorId() const = 0;
137   /** 
138    * Process a track reference 
139    * 
140    * @param particle Particle 
141    * @param mother   Ultimate mother (if not primary)
142    * @param ref      Reference 
143    * 
144    * @return 0 if no output should be generated for this reference, or
145    * pointer to track-reference to produce output for.
146    */
147   virtual AliTrackReference* ProcessRef(AliMCParticle* particle, 
148                                         const AliMCParticle* mother, 
149                                         AliTrackReference* ref) = 0;
150   /** 
151    * Called at before loop over track references
152    * 
153    */
154   virtual void BeginTrackRefs() {}
155   virtual Bool_t CheckTrackRef(AliTrackReference* /*ref*/) const { return true; }
156   /** 
157    * Called at before loop over track references
158    * 
159    */
160   virtual void EndTrackRefs(Int_t /*nRefs*/) {}
161   /** 
162    * Store a particle hit in Base<i>dr</i>[<i>s,t</i>] in @a output
163    * 
164    * 
165    * @param particle  Particle to store
166    * @param mother    Ultimate mother of particle 
167    * @param longest   Longest track reference
168    * @param vz        Z coordinate of IP
169    * @param nC        Total number of track-references in this sector  
170    * @param nT        Number of distint strips hit in this sector
171    * @param output    Output structure 
172    */  
173   virtual Double_t StoreParticle(AliMCParticle*       particle, 
174                                  const AliMCParticle* mother,
175                                  AliTrackReference*   ref) const;
176   /** 
177    * Get the parameters of this event 
178    * 
179    * @param event Event
180    * 
181    * @return true if found, false otherwise 
182    */
183   Bool_t GetCollisionParameters(const AliMCEvent& event);
184   /** 
185    * Get incident angle of this track reference
186    * 
187    * @param ref Track reference
188    * @param vz  Z coordinate of the IP
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
221   AliForwardFlowWeights fWeights; // Flow weights
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
227   ClassDef(AliBaseMCTrackDensity,2); // Calculate track-ref density
228 };
229
230 #endif
231 // Local Variables:
232 //  mode: C++ 
233 // End: