]>
Commit | Line | Data |
---|---|---|
1c5acb87 | 1 | #ifndef ALICALOTRACKREADER_H |
2 | #define ALICALOTRACKREADER_H | |
3 | /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
4 | * See cxx source for full Copyright notice */ | |
5 | /* $Id: $ */ | |
6 | ||
7 | //_________________________________________________________________________ | |
8 | // Base class for reading data: MonteCarlo, ESD or AOD, of PHOS EMCAL and | |
9 | // Central Barrel Tracking detectors. | |
10 | // Not all MC particles/tracks/clusters are kept, some kinematical restrictions are done. | |
591cc579 | 11 | // Mother class of : AliCaloTrackESDReader: Fills ESD data in 3 TObjArrays (PHOS, EMCAL, CTS) |
12 | // : AliCaloTrackMCReader: Fills Kinematics data in 3 TObjArrays (PHOS, EMCAL, CTS) | |
13 | // : AliCaloTrackAODReader: Fills AOD data in 3 TObjArrays (PHOS, EMCAL, CTS) | |
1c5acb87 | 14 | // |
15 | // -- Author: Gustavo Conesa (INFN-LNF) | |
16 | ||
17 | // --- ROOT system --- | |
18 | #include "TObject.h" | |
591cc579 | 19 | class TObjArray ; |
1c5acb87 | 20 | class TLorentzVector ; |
a97cd470 | 21 | #include "TString.h" |
8dacfd76 | 22 | #include "TObjArray.h" |
477d6cee | 23 | class TArrayF; |
591cc579 | 24 | class TTree ; |
1c5acb87 | 25 | |
26 | //--- ANALYSIS system --- | |
27 | class AliStack ; | |
28 | class AliHeader ; | |
29 | class AliGenEventHeader ; | |
477d6cee | 30 | class AliVEvent; |
8dacfd76 | 31 | class AliAODEvent; |
477d6cee | 32 | class AliMCEvent; |
ff45398a | 33 | class AliFiducialCut; |
591cc579 | 34 | class AliAODMCHeader; |
afabc52f | 35 | #include "AliPHOSGeoUtils.h" |
36 | #include "AliEMCALGeoUtils.h" | |
1c5acb87 | 37 | |
38 | class AliCaloTrackReader : public TObject { | |
39 | ||
477d6cee | 40 | public: |
41 | ||
1c5acb87 | 42 | AliCaloTrackReader() ; // ctor |
43 | AliCaloTrackReader(const AliCaloTrackReader & g) ; // cpy ctor | |
44 | AliCaloTrackReader & operator = (const AliCaloTrackReader & g) ;//cpy assignment | |
45 | virtual ~AliCaloTrackReader() ;//virtual dtor | |
46 | ||
47 | enum inputDataType {kESD, kAOD, kMC}; | |
48 | ||
591cc579 | 49 | //Select generated events, depending on comparison of pT hard and jets. |
50 | virtual Bool_t ComparePtHardAndJetPt() ; | |
51 | virtual Bool_t IsPtHardAndJetPtComparisonSet() const {return fComparePtHardAndJetPt ;} | |
52 | virtual void SetPtHardAndJetPtComparison(Bool_t compare) { fComparePtHardAndJetPt = compare ;} | |
53 | virtual Float_t GetPtHardAndJetFactor() const {return fPtHardAndJetPtFactor ;} | |
54 | virtual void SetPtHardAndJetPtFactor(Float_t factor) { fPtHardAndJetPtFactor = factor ;} | |
55 | ||
1c5acb87 | 56 | virtual void InitParameters(); |
57 | virtual void Print(const Option_t * opt) const; | |
58 | ||
59 | virtual Int_t GetDebug() const { return fDebug ; } | |
60 | virtual void SetDebug(Int_t d) { fDebug = d ; } | |
61 | virtual Int_t GetDataType() const { return fDataType ; } | |
62 | virtual void SetDataType(Int_t data ) { fDataType = data ; } | |
63 | ||
a79a2424 | 64 | virtual Int_t GetEventNumber() const {return fEventNumber ; } |
65 | virtual TString GetCurrentFileName() const {return fCurrentFileName ; } | |
66 | ||
1c5acb87 | 67 | //Minimum pt setters and getters |
68 | virtual Float_t GetEMCALPtMin() const { return fEMCALPtMin ; } | |
69 | virtual Float_t GetPHOSPtMin() const { return fPHOSPtMin ; } | |
70 | virtual Float_t GetCTSPtMin() const { return fCTSPtMin ; } | |
71 | ||
72 | virtual void SetEMCALPtMin(Float_t pt) { fEMCALPtMin = pt ; } | |
73 | virtual void SetPHOSPtMin(Float_t pt) { fPHOSPtMin = pt ; } | |
74 | virtual void SetCTSPtMin(Float_t pt) { fCTSPtMin = pt ; } | |
75 | ||
76 | //Input setters and getters | |
77 | ||
78 | Bool_t IsCTSSwitchedOn() const { return fFillCTS ; } | |
79 | void SwitchOnCTS() {fFillCTS = kTRUE ; } | |
80 | void SwitchOffCTS() {fFillCTS = kFALSE ; } | |
81 | ||
82 | Bool_t IsEMCALSwitchedOn() const { return fFillEMCAL ; } | |
83 | void SwitchOnEMCAL() {fFillEMCAL = kTRUE ; } | |
84 | void SwitchOffEMCAL() {fFillEMCAL = kFALSE ; } | |
85 | ||
86 | Bool_t IsPHOSSwitchedOn() const { return fFillPHOS ; } | |
87 | void SwitchOnPHOS() {fFillPHOS = kTRUE ; } | |
88 | void SwitchOffPHOS() {fFillPHOS = kFALSE ; } | |
89 | ||
90 | Bool_t IsEMCALCellsSwitchedOn() const { return fFillEMCALCells ; } | |
91 | void SwitchOnEMCALCells() {fFillEMCALCells = kTRUE ; } | |
92 | void SwitchOffEMCALCells() {fFillEMCALCells = kFALSE ; } | |
93 | ||
94 | Bool_t IsPHOSCellsSwitchedOn() const { return fFillPHOSCells ; } | |
95 | void SwitchOnPHOSCells() {fFillPHOSCells = kTRUE ; } | |
96 | void SwitchOffPHOSCells() {fFillPHOSCells = kFALSE ; } | |
97 | ||
29b2ceec | 98 | virtual Bool_t FillInputEvent(const Int_t iEntry, const char *currentFileName) ; |
1c5acb87 | 99 | virtual void FillInputCTS() {;} |
100 | virtual void FillInputEMCAL() {;} | |
101 | virtual void FillInputPHOS() {;} | |
102 | virtual void FillInputEMCALCells() {;} | |
103 | virtual void FillInputPHOSCells() {;} | |
104 | ||
591cc579 | 105 | virtual TObjArray* GetAODCTS() const {return fAODCTS ;} |
106 | virtual TObjArray* GetAODEMCAL() const {return fAODEMCAL ;} | |
107 | virtual TObjArray* GetAODPHOS() const {return fAODPHOS ;} | |
477d6cee | 108 | virtual TNamed* GetEMCALCells() const {return fEMCALCells ;} |
109 | virtual TNamed* GetPHOSCells() const {return fPHOSCells ;} | |
1c5acb87 | 110 | |
591cc579 | 111 | //Get MC informatio |
112 | //Kinematics and galice.root available | |
477d6cee | 113 | virtual AliStack* GetStack() const ; |
114 | virtual AliHeader* GetHeader() const ; | |
1c5acb87 | 115 | virtual AliGenEventHeader* GetGenEventHeader() const ; |
591cc579 | 116 | //Filtered kinematics in AOD |
117 | virtual TClonesArray* GetAODMCParticles(Int_t input = 0) const ; | |
118 | virtual AliAODMCHeader* GetAODMCHeader(Int_t input = 0) const ; | |
119 | ||
477d6cee | 120 | virtual AliVEvent* GetInputEvent() const {return fInputEvent;} |
121 | virtual AliAODEvent* GetOutputEvent() const {return fOutputEvent;} | |
122 | virtual AliMCEvent* GetMC() const {return fMC;} | |
233e0df8 | 123 | virtual void GetVertex(Double_t *) const {;} |
124 | virtual void GetSecondInputAODVertex(Double_t *) const {;} | |
8dacfd76 | 125 | virtual Double_t GetBField() const { return 0.;} |
591cc579 | 126 | |
127 | virtual void Init(); | |
128 | ||
8dacfd76 | 129 | virtual void SetInputEvent(AliVEvent* const input) {fInputEvent = input;} |
130 | virtual void SetOutputEvent(AliAODEvent* const aod) {fOutputEvent = aod;} | |
131 | virtual void SetMC(AliMCEvent* const mc) {fMC = mc;} | |
1c5acb87 | 132 | |
133 | virtual void ResetLists(); | |
134 | ||
ff45398a | 135 | virtual AliFiducialCut * GetFiducialCut() const {return fFiducialCut ;} |
136 | virtual void SetFiducialCut(AliFiducialCut * const fc) { fFiducialCut = fc ;} | |
29b2ceec | 137 | |
477d6cee | 138 | virtual void SetInputOutputMCEvent(AliVEvent* /*esd*/, AliAODEvent* /*aod*/, AliMCEvent* /*mc*/) {;} |
591cc579 | 139 | |
140 | //Methods for mixing with external input file (AOD) | |
141 | virtual TTree* GetSecondInputAODTree() const {return fSecondInputAODTree ; } | |
142 | //virtual void SetSecondInputAODTree(TTree * tree) {fSecondInputAODTree = tree ; | |
143 | // fSecondInputAODEvent->ReadFromTree(tree);}//Connect tree and AOD event. | |
144 | ||
145 | virtual AliAODEvent* GetSecondInputAODEvent() const { return fSecondInputAODEvent ; } | |
146 | ||
147 | TString GetSecondInputFileName() const {return fSecondInputFileName ; } | |
148 | void SetSecondInputFileName(TString name) { fSecondInputFileName = name ; } | |
1c5acb87 | 149 | |
591cc579 | 150 | Int_t GetSecondInputFirstEvent() const {return fSecondInputFirstEvent ; } |
151 | void SetSecondInputFirstEvent(Int_t iEvent0) { fSecondInputFirstEvent = iEvent0 ; } | |
152 | ||
153 | Int_t GetAODCTSNormalInputEntries() {if(!fSecondInputAODTree) { fAODCTSNormalInputEntries = fAODCTS->GetEntriesFast() ;} | |
154 | return fAODCTSNormalInputEntries ; } | |
155 | Int_t GetAODEMCALNormalInputEntries() {if(!fSecondInputAODTree) { fAODEMCALNormalInputEntries = fAODEMCAL->GetEntriesFast();} | |
156 | return fAODEMCALNormalInputEntries ; } | |
157 | Int_t GetAODPHOSNormalInputEntries() {if(!fSecondInputAODTree) { fAODPHOSNormalInputEntries = fAODPHOS->GetEntriesFast() ;} | |
158 | return fAODPHOSNormalInputEntries ; } | |
159 | ||
160 | ULong_t GetTrackStatus() const {return fTrackStatus ; } | |
161 | void SetTrackStatus(ULong_t bit) { fTrackStatus = bit ; } | |
162 | ||
163 | void SwitchOnStack() { fReadStack = kTRUE ; } | |
164 | void SwitchOffStack() { fReadStack = kFALSE ; } | |
165 | void SwitchOnAODMCParticles() { fReadAODMCParticles = kTRUE ; } | |
166 | void SwitchOffAODMCParticles() { fReadAODMCParticles = kFALSE ; } | |
167 | Bool_t ReadStack() const { return fReadStack ; } | |
168 | Bool_t ReadAODMCParticles() const { return fReadAODMCParticles ; } | |
169 | ||
08a064bc | 170 | void SwitchOnCleanStdAOD() {fCleanOutputStdAOD = kTRUE;} |
171 | void SwitchOffCleanStdAOD() {fCleanOutputStdAOD = kFALSE;} | |
172 | ||
42dc8e7d | 173 | void SetDeltaAODFileName(TString name ) {fDeltaAODFileName = name ; } |
174 | TString GetDeltaAODFileName() const {return fDeltaAODFileName ; } | |
72d2488e | 175 | |
176 | void SetFiredTriggerClassName(TString name ) {fFiredTriggerClassName = name ; } | |
177 | TString GetFiredTriggerClassName() const {return fFiredTriggerClassName ; } | |
178 | virtual TString GetFiredTriggerClasses() {return "";} | |
42dc8e7d | 179 | |
0866d83a | 180 | //Calorimeters Geometry Methods |
afabc52f | 181 | void SetEMCALGeometryName(TString name) { fEMCALGeoName = name ; } |
182 | TString EMCALGeometryName() const { return fEMCALGeoName ; } | |
4e2b43d8 | 183 | void InitEMCALGeometry() ; |
afabc52f | 184 | AliEMCALGeoUtils * GetEMCALGeometry() const {return fEMCALGeo;} |
0866d83a | 185 | Bool_t IsEMCALGeoMatrixSet() {return fEMCALGeoMatrixSet; } |
186 | ||
afabc52f | 187 | void SetPHOSGeometryName(TString name) { fPHOSGeoName = name ; } |
188 | TString PHOSGeometryName() const { return fPHOSGeoName ; } | |
4e2b43d8 | 189 | void InitPHOSGeometry() ; |
0866d83a | 190 | AliPHOSGeoUtils * GetPHOSGeometry() const {return fPHOSGeo;} |
191 | Bool_t IsPHOSGeoMatrixSet() {return fPHOSGeoMatrixSet ; } | |
192 | ||
afabc52f | 193 | |
1c5acb87 | 194 | protected: |
477d6cee | 195 | Int_t fEventNumber; // Event number |
a79a2424 | 196 | TString fCurrentFileName; // Current file name under analysis |
1c5acb87 | 197 | Int_t fDataType ; // Select MC:Kinematics, Data:ESD/AOD, MCData:Both |
198 | Int_t fDebug; // Debugging level | |
614701c6 | 199 | AliFiducialCut * fFiducialCut; //! Acceptance cuts |
29b2ceec | 200 | |
591cc579 | 201 | Bool_t fComparePtHardAndJetPt; // In MonteCarlo, jet events, reject fake events with wrong jet energy. |
202 | Float_t fPtHardAndJetPtFactor; // Factor between ptHard and jet pT to reject/accept event. | |
203 | ||
477d6cee | 204 | Float_t fCTSPtMin; // pT Threshold on charged particles |
1c5acb87 | 205 | Float_t fEMCALPtMin; // pT Threshold on emcal clusters |
477d6cee | 206 | Float_t fPHOSPtMin; // pT Threshold on phos clusters |
1c5acb87 | 207 | |
591cc579 | 208 | TObjArray * fAODCTS ; //! temporal referenced array with tracks |
209 | TObjArray * fAODEMCAL ; //! temporal referenced array with EMCAL CaloClusters | |
210 | TObjArray * fAODPHOS ; //! temporal referenced array with PHOS CaloClusters | |
1c5acb87 | 211 | TNamed * fEMCALCells ; //! temporal array with EMCAL CaloCells, ESD or AOD |
212 | TNamed * fPHOSCells ; //! temporal array with PHOS CaloCells, ESD or AOD | |
213 | ||
477d6cee | 214 | AliVEvent * fInputEvent; //! pointer to esd or aod input |
215 | AliAODEvent * fOutputEvent; //! pointer to aod output | |
1c5acb87 | 216 | AliMCEvent * fMC; //! Monte Carlo Event Handler |
217 | ||
218 | Bool_t fFillCTS; // use data from CTS | |
219 | Bool_t fFillEMCAL; // use data from EMCAL | |
220 | Bool_t fFillPHOS; // use data from PHOS | |
221 | Bool_t fFillEMCALCells; // use data from EMCAL | |
222 | Bool_t fFillPHOSCells; // use data from PHOS | |
223 | ||
591cc579 | 224 | TTree * fSecondInputAODTree; // Tree with second input AOD, for mixing analysis. |
225 | AliAODEvent* fSecondInputAODEvent; //! pointer to second input AOD event. | |
226 | TString fSecondInputFileName; // File with AOD data to mix with normal stream of data. | |
227 | Int_t fSecondInputFirstEvent; // First event to be considered in the mixing. | |
228 | ||
229 | Int_t fAODCTSNormalInputEntries; // Number of entries in CTS in case of standard input, larger with mixing. | |
230 | Int_t fAODEMCALNormalInputEntries; // Number of entries in EMCAL in case of standard input, larger with mixing. | |
231 | Int_t fAODPHOSNormalInputEntries; // Number of entries in PHOS in case of standard input, larger with mixing. | |
232 | ||
233 | ULong_t fTrackStatus ; // Track selection bit, select tracks refitted in TPC, ITS ... | |
234 | Bool_t fReadStack ; // Access kine information from stack | |
235 | Bool_t fReadAODMCParticles ; // Access kine information from filtered AOD MC particles | |
236 | ||
08a064bc | 237 | Bool_t fCleanOutputStdAOD; // clean the written standard tracks and caloclusters in output AOD |
42dc8e7d | 238 | TString fDeltaAODFileName ; // Delta AOD file name |
72d2488e | 239 | TString fFiredTriggerClassName ; // Name of trigger event type used to do the analysis |
afabc52f | 240 | |
241 | TString fEMCALGeoName; // Name of geometry to use for EMCAL. | |
242 | TString fPHOSGeoName; // Name of geometry to use for PHOS. | |
243 | AliEMCALGeoUtils * fEMCALGeo ; //! EMCAL geometry pointer | |
0866d83a | 244 | AliPHOSGeoUtils * fPHOSGeo ; //! PHOS geometry pointer |
245 | Bool_t fEMCALGeoMatrixSet; // Check if the transformation matrix is set for EMCAL | |
246 | Bool_t fPHOSGeoMatrixSet ; // Check if the transformation matrix is set for PHOS | |
afabc52f | 247 | |
0866d83a | 248 | ClassDef(AliCaloTrackReader,11) |
1c5acb87 | 249 | } ; |
250 | ||
251 | ||
252 | #endif //ALICALOTRACKREADER_H | |
253 | ||
254 | ||
255 |