3 // Version 2 of my Cone Jet Finder
4 // improved storage (see also version 1)
6 #ifndef ALITKCONJETFINDERV2_H
7 #define ALITKCONJETFINDERV2_H
10 #include <TParticle.h>
18 tower(const tower& t);
19 tower(Float_t phimin, Float_t phimax, Float_t etamin, Float_t etamax);
22 tower& operator+=(const Float_t E);
23 tower& operator+=(const TParticle *part);
25 Float_t getEtaMin() const {return fEta_min;}
26 Float_t getEtaMax() const {return fEta_max;}
27 Float_t getEta() const {return fEta_center;}
28 Float_t getPhiMin() const {return fPhi_min;}
29 Float_t getPhiMax() const {return fPhi_max;}
30 Float_t getPhi() const {return fPhi_center;}
31 Float_t getEt() const {return fEt;}
32 Float_t getEtaWidth() const {return (fEta_max-fEta_min);}
33 Float_t getPhiWidth() const {return (fPhi_max-fPhi_min);}
35 void setEtaMin(Float_t f) {fEta_min=f;}
36 void setEtaMax(Float_t f) {fEta_max=f;}
37 void setEta(Float_t f) {fEta_center=f;}
38 void setPhiMin(Float_t f) {fPhi_min=f;}
39 void setPhiMax(Float_t f) {fPhi_max=f;}
40 void setPhi(Float_t f) {fPhi_center=f;}
42 void addParticle(const TParticle* part) {fParticles.push_back(part);}
43 void clearParticles() {fParticles.erase(fParticles.begin(),fParticles.end());}
44 void clear() {fEt = 0; clearParticles();}
46 list<const TParticle*> *getParticles();
57 list<const TParticle*> fParticles;
60 ostream& operator<<(ostream& s,const tower& t);
62 //-----------------------------------------------------
64 #include "AliTkEtaPhiVector.h"
70 protojet(const protojet& p);
71 protojet(const protojet *p);
72 ~protojet(){eraseTowers();}
74 Float_t Eta() const {return fCentroid.Eta();}
75 Float_t Phi() const {return fCentroid.Phi();}
76 Float_t getEt() const {
77 if(fUpdate) cerr << "Eta update needed" << endl;
85 Int_t getNTowers() const { return fTowers.size(); }
87 void setCentroidPosition(AliTkEtaPhiVector center) {fCentroid = center;}
89 AliTkEtaPhiVector getCentroidPosition() const {
92 AliTkEtaPhiVector getEtWeightedPosition() const {
93 if(fUpdate) cerr << "Weighted update needed" << endl;
94 return fEtWeightedCentroid;
96 AliTkEtaPhiVector getEtWeightedPosition() {
98 return fEtWeightedCentroid;
101 void addTower(tower *tower) {fTowers.push_back(tower);fUpdate=kTRUE;}
102 void eraseTowers() {fTowers.erase(fTowers.begin(),fTowers.end()); }
103 void eraseTower(tower *pTower) {fTowers.remove(pTower);}
104 void clear(){eraseTowers();}
106 list<tower*> getTowerList() const;
107 list<tower*> getSharedTowerList(const protojet *other) const;
109 bool hasTower(tower *pTower) const;
110 bool shareTowers(const protojet *other) const;
111 Float_t diffToCenter(tower *pTower) const;
112 bool operator<(const protojet &p1) const {return (getEt() < p1.getEt());}
113 bool operator==(const protojet &p1);
117 AliTkEtaPhiVector fCentroid;
118 AliTkEtaPhiVector fEtWeightedCentroid;
121 list<tower*> fTowers;
124 ostream& operator<<(ostream& s,const protojet& p);
126 //-----------------------------------------------------
127 // Real ConeJetFinderClass
128 //-----------------------------------------------------
131 #include <TClonesArray.h>
134 #include "AliTkConeJet.h"
135 #include "AliTkConeJetEvent.h"
136 #ifdef ALICEINTERFACE
137 #include <AliJetEventParticles.h>
140 //if defined produce some debugging histos
143 class AliTkConeJetFinderV2 : public TObject {
146 AliTkConeJetFinderV2();
147 virtual ~AliTkConeJetFinderV2();
150 void defaultSettings();
151 void setSettings(Int_t phibins,Int_t etabins);
152 void setEtMinJet(Float_t et){fEtMinJet=et;} //minimum jet energy required
153 void setEtCut(Float_t et){fEtCut=et;} //min et for seedpoints
154 void setPtCut(Float_t pt){fPtCut=pt;} //set pt cut
155 void setOutput(Bool_t out){fOutput=out;}
156 void setRadius(Float_t r=0.7) {fRadius=r;}
159 void initEvent(const TClonesArray *particles,Int_t type);
160 void initEvent(const TClonesArray *particles,Int_t type,TString desc);
165 // real physics functions...
167 void fillTowersFromTParticles(const TClonesArray *particles);
168 #ifdef ALICEINTERFACE
169 void fillTowersFromAliParticles(const TClonesArray *particles);
170 void initEvent(const AliJetEventParticles *p,TString desc);
173 void createSeedPoints();
174 void findProtojets();
178 // analysis functions...
180 void createEventHistos();
181 void clearEventHistos();
182 void writeEventHistos();
187 void setEvOutFilename(const Char_t *filename);
188 Char_t *getEvOutFilename() {return fEvout_name;}
190 void addProtojet(protojet *pj);
191 void dumpProtojets(Float_t etmin = 5.0);
193 void splitMergeJets(protojet *jet1, protojet *jet2);
194 void splitJets(protojet *jet1,protojet *jet2);
195 void mergeJets(protojet *jet1,protojet *jet2);
196 void addJet(protojet *pj) {fJets.push_back(pj);}
199 Bool_t isJetEnergy(Float_t min,Float_t max);
200 Float_t maxJetEnergy(); //of protojets (before calling finishEvent)
201 AliTkConeJetEvent* getEvent() const {return fEvoutevent;}
204 Bool_t getOutput() const {return fOutput;}
205 Int_t getNTowers() const {return fNTowers;}
206 Int_t getNEtaBins() const {return fEtaBins;}
207 Float_t getEtaMin() const {return fEtaMin;}
208 Float_t getEtaMax() const {return fEtaMax;}
209 Float_t getEtaWidth() const {return fEtaWidth;}
210 Int_t getPhiBins() const {return fPhiBins;}
211 Float_t getPhiMin() const {return fPhiMin;}
212 Float_t getPhiMax() const {return fPhiMax;}
213 Float_t getPhiWidth() const {return fPhiWidth;}
214 Float_t getEtCut() const {return fEtCut;};
215 Float_t getPtCut() const {return fPtCut;};
216 Float_t getEtMinJet() const {return fEtMinJet;}
217 Float_t getRadius() const {return fRadius;}
220 void initEvent_(const TClonesArray *particles,Int_t type);
221 void fillTowerHist();
222 Int_t findTower(Float_t phi, Float_t eta);
223 Bool_t isTParticleAccepted(TParticle *particle);
225 Bool_t fOutput; //outpug some information if true
226 Int_t fNTowers; //# of towers
227 Int_t fEtaBins; //eta-phi grid
235 Float_t fEtCut; // cut for seedpoints
236 Float_t fEtMinJet; // cut for jets
237 Float_t fPtCut; // cut for particles
240 //towers, seeds, protojets and jets
241 vector<tower> *fTowers; //!
242 list<AliTkEtaPhiVector> fSeedPointsNew; //!
243 list<protojet*> fProtojets; //!
244 list<protojet*> fJets; //!
247 TFile *fEvoutfile; //!
248 Char_t *fEvout_name; //!
249 AliTkConeJetEvent *fEvoutevent; //!
250 TTree *fEvouttree; //!
252 #ifdef ALICEINTERFACE
253 TClonesArray *fAliParticles; //!
259 TObjArray fEventHistos;
260 TFile *fHistFile; //!
263 ClassDef(AliTkConeJetFinderV2,1)
266 inline Bool_t AliTkConeJetFinderV2::isTParticleAccepted(TParticle * /* particle */)
268 // check if particle is accepted
269 // makes sense to write this into a own class, but now I'm lazy
271 // check if particle is stable -> a detectable particle
272 #ifndef ALICEINTERFACE
273 if (particle->GetStatusCode()%100 != 1) {
278 // default case: accept