]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/JET/AliHLTJETReader.h
Removed depricted histograms ... cleanup
[u/mrichter/AliRoot.git] / HLT / JET / AliHLTJETReader.h
CommitLineData
0734d112 1//-*- Mode: C++ -*-
2
3// $Id: $
4
5#ifndef ALIHLTJETREADER_H
6#define ALIHLTJETREADER_H
7
8/* This file is property of and copyright by the ALICE HLT Project *
9 * ALICE Experiment at CERN, All rights reserved. *
10 * See cxx source for full Copyright notice */
11
12/** @file AliHLTJETReader.h
13 @author Jochen Thaeder
14 @date
15 @brief Reader for jet finder
16*/
17
6ce099ba 18#ifdef HAVE_FASTJET
19#include "fastjet/PseudoJet.hh"
20#endif
0734d112 21
22#include "AliJetReader.h"
23#include "AliJetReaderHeader.h"
24
6ce099ba 25#include "AliESDEvent.h"
26#include "AliAODEvent.h"
27
6ce099ba 28#include "AliHLTLogging.h"
0734d112 29#include "AliHLTMCEvent.h"
30
6ce099ba 31#include "AliHLTJETBase.h"
32#include "AliHLTJETTrackCuts.h"
33#include "AliHLTJETReaderHeader.h"
0734d112 34
6ce099ba 35#include "AliHLTJETConeSeedCuts.h"
36#include "AliHLTJETConeGrid.h"
0734d112 37
38/**
39 * @class AliHLTJETReader
1f9fec4a 40 * This class is the reader class for the JetFinder in the HLT
41 * It implements the reading of ESDs and MCs. AOD reading is
42 * not yet implemented
43 * <br>
44 * Usage :<br>
45 * - Initilize() // Initializes the reader dependent of the algorithm
0734d112 46 *
47 * @ingroup alihlt_jet
48 */
49
50class AliHLTJETReader : public AliJetReader, public AliHLTLogging {
51
52public:
53
54 /*
55 * ---------------------------------------------------------------------------------
56 * Constructor / Destructor
57 * ---------------------------------------------------------------------------------
58 */
59
60 /** standard constructor */
61 AliHLTJETReader();
62
63 /** destructor */
64 virtual ~AliHLTJETReader();
65
66 /*
67 * ---------------------------------------------------------------------------------
6ce099ba 68 * Initialize / Reset
69 * ---------------------------------------------------------------------------------
70 */
71
1f9fec4a 72 /** Initialize reader
6ce099ba 73 * Calls AliHLTJETReaderHeader::Initialize
1f9fec4a 74 * and the private Initialize methods
6ce099ba 75 * @return 0 on success, otherwise <0
76 */
77 Int_t Initialize();
78
79 /** Reset the event */
80 void ResetEvent();
81
1f9fec4a 82 /*
83 * ---------------------------------------------------------------------------------
84 * Setter
85 * ---------------------------------------------------------------------------------
86 */
87
88 /** Set pointer to input event
89 * Needs "useMC" flag for running in analysis task only
90 * @param esd an AliESDEvent
91 * @param aod an AliAODEvent
92 * @param mc an AliHLTMCEvent
93 */
94 void SetInputEvent(const TObject* esd, const TObject* aod, const TObject* mc);
95
96 /** Set number of jet candates = seeds */
97 void SetNJetCandidates( Int_t i ) { fNJetCandidates = i; }
98
6ce099ba 99 /*
100 * ---------------------------------------------------------------------------------
101 * Fastjet Reader functionality
0734d112 102 * ---------------------------------------------------------------------------------
103 */
104
6ce099ba 105#ifdef HAVE_FASTJET
106 /** Fill tracks in fastjet momemtum vector
107 * @return kTRUE on success, otherwise kFALSE
108 */
1f9fec4a 109 Bool_t FillVector();
6ce099ba 110
111 /** Fill MC tracks in fastjet momemtum vector
112 * @return kTRUE on success, otherwise kFALSE
113 */
1f9fec4a 114 Bool_t FillVectorMC();
115
116 /** Fill HLT MC tracks in fastjet momemtum vector
117 * @return kTRUE on success, otherwise kFALSE
118 */
119 Bool_t FillVectorHLTMC();
6ce099ba 120
121 /** Fill ESD tracks in fastjet momemtum vector
122 * @return kTRUE on success, otherwise kFALSE
123 */
1f9fec4a 124 Bool_t FillVectorESD();
6ce099ba 125
126 /** Fill AOD tracks in fastjet momemtum vector
127 * @return kTRUE on success, otherwise kFALSE
128 */
1f9fec4a 129 Bool_t FillVectorAOD();
6ce099ba 130#endif
131
132 /*
133 * ---------------------------------------------------------------------------------
134 * Grid functionality
135 * ---------------------------------------------------------------------------------
136 */
137
0734d112 138 /** Fill tracks in momentum array
139 * @return kTRUE on success, otherwise kFALSE
140 */
6ce099ba 141 Bool_t FillGrid();
0734d112 142
143 /** Fill MC tracks in momentum array
144 * @return kTRUE on success, otherwise kFALSE
145 */
6ce099ba 146 Bool_t FillGridMC();
0734d112 147
1f9fec4a 148 /** Fill HLT MC tracks in momentum array
149 * @return kTRUE on success, otherwise kFALSE
150 */
151 Bool_t FillGridHLTMC();
152
0734d112 153 /** Fill ESD tracks in momentum array
154 * @return kTRUE on success, otherwise kFALSE
155 */
6ce099ba 156 Bool_t FillGridESD();
0734d112 157
158 /** Fill AOD tracks in momentum array
159 * @return kTRUE on success, otherwise kFALSE
160 */
6ce099ba 161 Bool_t FillGridAOD();
0734d112 162
0734d112 163 /*
164 * ---------------------------------------------------------------------------------
165 * Getter
166 * ---------------------------------------------------------------------------------
167 */
168
169 /** Get Ptr to AliHLTJETReaderHeader
170 * @return ptr to AliHLTJETReaderHeader
171 */
77c2a39b 172 AliHLTJETReaderHeader* GetReaderHeader() const { return dynamic_cast<AliHLTJETReaderHeader*>(fReaderHeader);}
6ce099ba 173
174#ifdef HAVE_FASTJET
175 /** Get Ptr to input vector of Fastjet
176 * @return ptr to input vector of Fastjet
177 */
1f9fec4a 178 vector<fastjet::PseudoJet>* GetVector() { return fMomentumVector; }
6ce099ba 179#endif
180
181 /** Get Ptr to grid of cone finder
182 * @return ptr to grid of cone finder
183 */
1f9fec4a 184 AliHLTJETConeGrid* GetGrid() { return fGrid; }
185
6ce099ba 186 /** Get number of jet candates = seeds */
1f9fec4a 187 Int_t GetNJetCandidates() { return fNJetCandidates; }
6ce099ba 188
189 /** Get ptr to jet candiates = seeds for cone finder */
1f9fec4a 190 TClonesArray* GetJetCandidates() { return fJetCandidates; }
6ce099ba 191
192 /*
193 * ---------------------------------------------------------------------------------
194 * Seeds
195 * ---------------------------------------------------------------------------------
196 */
197
198 /** Add new seed
199 * @param aEtaPhi eta and phi of the seed
200 * @param aGridIdx indeces in the grid
201 * @param coneRadius coneRadius
202 */
203 void AddSeed( const Float_t* aEtaPhi, const Int_t* aGridIdx,
204 const Float_t coneRadius );
0734d112 205
206 ///////////////////////////////////////////////////////////////////////////////////
207
208private:
209
210 /** copy constructor prohibited */
211 AliHLTJETReader (const AliHLTJETReader&);
212
213 /** assignment operator prohibited */
214 AliHLTJETReader& operator= (const AliHLTJETReader&);
215
1f9fec4a 216 /*
217 * ---------------------------------------------------------------------------------
218 * Initialize - private
219 * ---------------------------------------------------------------------------------
220 */
221
222 /** Initialize reader for the FFSC cone jet finder
223 * @return 0 on success, otherwise <0
224 */
225 Int_t InitializeFFSC();
226
227#ifdef HAVE_FASTJET
228 /** Initialize reader for the fastjet jet finders
229 * @return 0 on success, otherwise <0
230 */
231 Int_t InitializeFastjet();
232#endif
233
0734d112 234 /*
235 * ---------------------------------------------------------------------------------
236 * Members - private
237 * ---------------------------------------------------------------------------------
238 */
239
1f9fec4a 240 // -- Input
241 // ----------
242
0734d112 243 /** ESD event */
6ce099ba 244 AliESDEvent *fESD; //! transient
0734d112 245
1f9fec4a 246 /** off-line MC event */
247 AliMCEvent *fMC; //! transient
248
249 /** on-line MC event */
250 AliHLTMCEvent *fHLTMC; //! transient
0734d112 251
252 /** AOD event */
6ce099ba 253 AliAODEvent *fAOD; //! transient
254
1f9fec4a 255 // -- Particle structures
256 // ------------------------
257
6ce099ba 258#ifdef HAVE_FASTJET
259 /** Vector of fastjet momemtum entries */
260 vector<fastjet::PseudoJet> *fMomentumVector; //! transient
261#endif
262
263 /** Grid for cone finder */
264 AliHLTJETConeGrid *fGrid; //! transient
265
1f9fec4a 266 // -- Output
267 // -----------
268
6ce099ba 269 /** Number of jet candates = seeds */
270 Int_t fNJetCandidates; // see above
271
272 /** Jet candiates = seeds for cone finder */
273 TClonesArray *fJetCandidates; //! transient
274
1f9fec4a 275 // -- Cuts
276 // ---------
277
6ce099ba 278 /** Ptr to seed cuts */
279 AliHLTJETConeSeedCuts *fSeedCuts; //! transient
280
281 /** Ptr to track cuts */
282 AliHLTJETTrackCuts *fTrackCuts; //! transient
0734d112 283
284 ClassDef(AliHLTJETReader, 1)
285
286};
287#endif
288