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 */
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.
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)
15 // -- Author: Gustavo Conesa (INFN-LNF)
17 // --- ROOT system ---
20 class TLorentzVector ;
22 #include "TObjArray.h"
26 //--- ANALYSIS system ---
29 class AliGenEventHeader ;
36 class AliCaloTrackReader : public TObject {
40 AliCaloTrackReader() ; // ctor
41 AliCaloTrackReader(const AliCaloTrackReader & g) ; // cpy ctor
42 AliCaloTrackReader & operator = (const AliCaloTrackReader & g) ;//cpy assignment
43 virtual ~AliCaloTrackReader() ;//virtual dtor
45 enum inputDataType {kESD, kAOD, kMC};
47 //Select generated events, depending on comparison of pT hard and jets.
48 virtual Bool_t ComparePtHardAndJetPt() ;
49 virtual Bool_t IsPtHardAndJetPtComparisonSet() const {return fComparePtHardAndJetPt ;}
50 virtual void SetPtHardAndJetPtComparison(Bool_t compare) { fComparePtHardAndJetPt = compare ;}
51 virtual Float_t GetPtHardAndJetFactor() const {return fPtHardAndJetPtFactor ;}
52 virtual void SetPtHardAndJetPtFactor(Float_t factor) { fPtHardAndJetPtFactor = factor ;}
54 virtual void InitParameters();
55 virtual void Print(const Option_t * opt) const;
57 virtual Int_t GetDebug() const { return fDebug ; }
58 virtual void SetDebug(Int_t d) { fDebug = d ; }
59 virtual Int_t GetDataType() const { return fDataType ; }
60 virtual void SetDataType(Int_t data ) { fDataType = data ; }
62 virtual Int_t GetEventNumber() const {return fEventNumber ; }
63 virtual TString GetCurrentFileName() const {return fCurrentFileName ; }
65 //Minimum pt setters and getters
66 virtual Float_t GetEMCALPtMin() const { return fEMCALPtMin ; }
67 virtual Float_t GetPHOSPtMin() const { return fPHOSPtMin ; }
68 virtual Float_t GetCTSPtMin() const { return fCTSPtMin ; }
70 virtual void SetEMCALPtMin(Float_t pt) { fEMCALPtMin = pt ; }
71 virtual void SetPHOSPtMin(Float_t pt) { fPHOSPtMin = pt ; }
72 virtual void SetCTSPtMin(Float_t pt) { fCTSPtMin = pt ; }
74 //Input setters and getters
76 Bool_t IsCTSSwitchedOn() const { return fFillCTS ; }
77 void SwitchOnCTS() {fFillCTS = kTRUE ; }
78 void SwitchOffCTS() {fFillCTS = kFALSE ; }
80 Bool_t IsEMCALSwitchedOn() const { return fFillEMCAL ; }
81 void SwitchOnEMCAL() {fFillEMCAL = kTRUE ; }
82 void SwitchOffEMCAL() {fFillEMCAL = kFALSE ; }
84 Bool_t IsPHOSSwitchedOn() const { return fFillPHOS ; }
85 void SwitchOnPHOS() {fFillPHOS = kTRUE ; }
86 void SwitchOffPHOS() {fFillPHOS = kFALSE ; }
88 Bool_t IsEMCALCellsSwitchedOn() const { return fFillEMCALCells ; }
89 void SwitchOnEMCALCells() {fFillEMCALCells = kTRUE ; }
90 void SwitchOffEMCALCells() {fFillEMCALCells = kFALSE ; }
92 Bool_t IsPHOSCellsSwitchedOn() const { return fFillPHOSCells ; }
93 void SwitchOnPHOSCells() {fFillPHOSCells = kTRUE ; }
94 void SwitchOffPHOSCells() {fFillPHOSCells = kFALSE ; }
96 virtual Bool_t FillInputEvent(const Int_t iEntry, const char *currentFileName) ;
97 virtual void FillInputCTS() {;}
98 virtual void FillInputEMCAL() {;}
99 virtual void FillInputPHOS() {;}
100 virtual void FillInputEMCALCells() {;}
101 virtual void FillInputPHOSCells() {;}
103 virtual TObjArray* GetAODCTS() const {return fAODCTS ;}
104 virtual TObjArray* GetAODEMCAL() const {return fAODEMCAL ;}
105 virtual TObjArray* GetAODPHOS() const {return fAODPHOS ;}
106 virtual TNamed* GetEMCALCells() const {return fEMCALCells ;}
107 virtual TNamed* GetPHOSCells() const {return fPHOSCells ;}
110 //Kinematics and galice.root available
111 virtual AliStack* GetStack() const ;
112 virtual AliHeader* GetHeader() const ;
113 virtual AliGenEventHeader* GetGenEventHeader() const ;
114 //Filtered kinematics in AOD
115 virtual TClonesArray* GetAODMCParticles(Int_t input = 0) const ;
116 virtual AliAODMCHeader* GetAODMCHeader(Int_t input = 0) const ;
118 virtual AliVEvent* GetInputEvent() const {return fInputEvent;}
119 virtual AliAODEvent* GetOutputEvent() const {return fOutputEvent;}
120 virtual AliMCEvent* GetMC() const {return fMC;}
121 virtual void GetVertex(Double_t *) const {;}
122 virtual void GetSecondInputAODVertex(Double_t *) const {;}
123 virtual Double_t GetBField() const { return 0.;}
127 virtual void SetInputEvent(AliVEvent* const input) {fInputEvent = input;}
128 virtual void SetOutputEvent(AliAODEvent* const aod) {fOutputEvent = aod;}
129 virtual void SetMC(AliMCEvent* const mc) {fMC = mc;}
131 virtual void ResetLists();
133 virtual AliFidutialCut * GetFidutialCut() const {return fFidutialCut ;}
134 virtual void SetFidutialCut(AliFidutialCut * const fc) { fFidutialCut = fc ;}
136 virtual void SetInputOutputMCEvent(AliVEvent* /*esd*/, AliAODEvent* /*aod*/, AliMCEvent* /*mc*/) {;}
138 //Methods for mixing with external input file (AOD)
139 virtual TTree* GetSecondInputAODTree() const {return fSecondInputAODTree ; }
140 //virtual void SetSecondInputAODTree(TTree * tree) {fSecondInputAODTree = tree ;
141 // fSecondInputAODEvent->ReadFromTree(tree);}//Connect tree and AOD event.
143 virtual AliAODEvent* GetSecondInputAODEvent() const { return fSecondInputAODEvent ; }
145 TString GetSecondInputFileName() const {return fSecondInputFileName ; }
146 void SetSecondInputFileName(TString name) { fSecondInputFileName = name ; }
148 Int_t GetSecondInputFirstEvent() const {return fSecondInputFirstEvent ; }
149 void SetSecondInputFirstEvent(Int_t iEvent0) { fSecondInputFirstEvent = iEvent0 ; }
151 Int_t GetAODCTSNormalInputEntries() {if(!fSecondInputAODTree) { fAODCTSNormalInputEntries = fAODCTS->GetEntriesFast() ;}
152 return fAODCTSNormalInputEntries ; }
153 Int_t GetAODEMCALNormalInputEntries() {if(!fSecondInputAODTree) { fAODEMCALNormalInputEntries = fAODEMCAL->GetEntriesFast();}
154 return fAODEMCALNormalInputEntries ; }
155 Int_t GetAODPHOSNormalInputEntries() {if(!fSecondInputAODTree) { fAODPHOSNormalInputEntries = fAODPHOS->GetEntriesFast() ;}
156 return fAODPHOSNormalInputEntries ; }
158 ULong_t GetTrackStatus() const {return fTrackStatus ; }
159 void SetTrackStatus(ULong_t bit) { fTrackStatus = bit ; }
161 void SwitchOnStack() { fReadStack = kTRUE ; }
162 void SwitchOffStack() { fReadStack = kFALSE ; }
163 void SwitchOnAODMCParticles() { fReadAODMCParticles = kTRUE ; }
164 void SwitchOffAODMCParticles() { fReadAODMCParticles = kFALSE ; }
165 Bool_t ReadStack() const { return fReadStack ; }
166 Bool_t ReadAODMCParticles() const { return fReadAODMCParticles ; }
169 Int_t fEventNumber; // Event number
170 TString fCurrentFileName; // Current file name under analysis
171 Int_t fDataType ; // Select MC:Kinematics, Data:ESD/AOD, MCData:Both
172 Int_t fDebug; // Debugging level
173 AliFidutialCut * fFidutialCut; // Acceptance cuts
175 Bool_t fComparePtHardAndJetPt; // In MonteCarlo, jet events, reject fake events with wrong jet energy.
176 Float_t fPtHardAndJetPtFactor; // Factor between ptHard and jet pT to reject/accept event.
178 Float_t fCTSPtMin; // pT Threshold on charged particles
179 Float_t fEMCALPtMin; // pT Threshold on emcal clusters
180 Float_t fPHOSPtMin; // pT Threshold on phos clusters
182 TObjArray * fAODCTS ; //! temporal referenced array with tracks
183 TObjArray * fAODEMCAL ; //! temporal referenced array with EMCAL CaloClusters
184 TObjArray * fAODPHOS ; //! temporal referenced array with PHOS CaloClusters
185 TNamed * fEMCALCells ; //! temporal array with EMCAL CaloCells, ESD or AOD
186 TNamed * fPHOSCells ; //! temporal array with PHOS CaloCells, ESD or AOD
188 AliVEvent * fInputEvent; //! pointer to esd or aod input
189 AliAODEvent * fOutputEvent; //! pointer to aod output
190 AliMCEvent * fMC; //! Monte Carlo Event Handler
192 Bool_t fFillCTS; // use data from CTS
193 Bool_t fFillEMCAL; // use data from EMCAL
194 Bool_t fFillPHOS; // use data from PHOS
195 Bool_t fFillEMCALCells; // use data from EMCAL
196 Bool_t fFillPHOSCells; // use data from PHOS
198 TTree * fSecondInputAODTree; // Tree with second input AOD, for mixing analysis.
199 AliAODEvent* fSecondInputAODEvent; //! pointer to second input AOD event.
200 TString fSecondInputFileName; // File with AOD data to mix with normal stream of data.
201 Int_t fSecondInputFirstEvent; // First event to be considered in the mixing.
203 Int_t fAODCTSNormalInputEntries; // Number of entries in CTS in case of standard input, larger with mixing.
204 Int_t fAODEMCALNormalInputEntries; // Number of entries in EMCAL in case of standard input, larger with mixing.
205 Int_t fAODPHOSNormalInputEntries; // Number of entries in PHOS in case of standard input, larger with mixing.
207 ULong_t fTrackStatus ; // Track selection bit, select tracks refitted in TPC, ITS ...
208 Bool_t fReadStack ; // Access kine information from stack
209 Bool_t fReadAODMCParticles ; // Access kine information from filtered AOD MC particles
211 ClassDef(AliCaloTrackReader,6)
215 #endif //ALICALOTRACKREADER_H