1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 //---------------------------------------------------------------------
18 // FastJet v2.3.4 finder algorithm interface
20 // Author: Rafael.Diaz.Valdes@cern.ch
22 //---------------------------------------------------------------------
24 #include <Riostream.h>
25 #include <TLorentzVector.h>
31 #include <TClonesArray.h>
33 #include "AliFastJetFinder.h"
34 #include "AliFastJetHeader.h"
35 #include "AliJetReaderHeader.h"
36 #include "AliJetReader.h"
40 #include "fastjet/PseudoJet.hh"
41 #include "fastjet/ClusterSequenceArea.hh"
42 #include "fastjet/AreaDefinition.hh"
43 #include "fastjet/JetDefinition.hh"
44 // get info on how fastjet was configured
45 #include "fastjet/config.h"
47 #ifdef ENABLE_PLUGIN_SISCONE
48 #include "fastjet/SISConePlugin.hh"
51 #include<sstream> // needed for internal io
58 ClassImp(AliFastJetFinder)
59 //____________________________________________________________________________
61 AliFastJetFinder::AliFastJetFinder():
67 //____________________________________________________________________________
69 AliFastJetFinder::~AliFastJetFinder()
74 //______________________________________________________________________________
75 void AliFastJetFinder::FindJets()
78 Bool_t debug = kFALSE;
80 //pick up fastjet header
81 AliFastJetHeader *header = (AliFastJetHeader*)fHeader;
83 // check if we are reading AOD jets
85 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
86 if (fromAod) { refs = fReader->GetReferences(); }
89 // read input particles -----------------------------
90 vector<fastjet::PseudoJet> input_particles;
91 TClonesArray *lvArray = fReader->GetMomentumArray();
92 if(lvArray == 0) { cout << "Could not get the momentum array" << endl; return; }
93 Int_t nIn = lvArray->GetEntries();
94 if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
95 //Int_t nJets = 0; // n jets in this event
96 fJets->SetNinput(nIn) ; // number of input objects
99 for(Int_t i = 0; i < nIn; i++){ // loop for all input particles
100 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
106 fastjet::PseudoJet input_part(px,py,pz,en); // create PseudoJet object
107 input_part.set_user_index(i); //label the particle into Fastjet algortihm
108 input_particles.push_back(input_part); // back of the input_particles vector
112 // create an object that represents your choice of jet algorithm, and
113 // the associated parameters
114 double Rparam = header->GetRparam();
115 fastjet::Strategy strategy = header->GetStrategy();
116 fastjet::RecombinationScheme recomb_scheme = header->GetRecombScheme();
117 fastjet::JetAlgorithm algorithm = header->GetAlgorithm();
118 fastjet::JetDefinition jet_def(algorithm, Rparam, recomb_scheme, strategy);
121 // create an object that specifies how we to define the area
122 fastjet::AreaDefinition area_def;
123 double ghost_etamax = header->GetGhostEtaMax();
124 double ghost_area = header->GetGhostArea();
125 int active_area_repeats = header->GetActiveAreaRepeats();
127 // now create the object that holds info about ghosts
128 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, ghost_area);
129 // and from that get an area definition
130 fastjet::AreaType area_type = header->GetAreaType();
131 area_def = fastjet::AreaDefinition(area_type,ghost_spec);
133 // run the jet clustering with the above jet definition
134 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def);
137 // save a comment in the header
139 TString comment = "Running FastJet algorithm with the following setup. ";
140 comment+= "Jet definition: ";
141 comment+= TString(jet_def.description());
142 comment+= ". Area definition: ";
143 comment+= TString(area_def.description());
144 comment+= ". Strategy adopted by FastJet: ";
145 comment+= TString(clust_seq.strategy_string());
146 header->SetComment(comment);
148 cout << "--------------------------------------------------------" << endl;
149 cout << comment << endl;
150 cout << "--------------------------------------------------------" << endl;
152 //header->PrintParameters();
155 // extract the inclusive jets with pt > ptmin, sorted by pt
156 double ptmin = header->GetPtMin();
157 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
159 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
162 //subtract background // ===========================================
163 // set the rapididty , phi range within which to study the background
164 double rap_max = header->GetRapMax();
165 double rap_min = header->GetRapMin();
166 double phi_max = header->GetPhiMax();
167 double phi_min = header->GetPhiMin();
168 fastjet::RangeDefinition range(rap_min, rap_max, phi_min, phi_max);
169 // subtract background
170 vector<fastjet::PseudoJet> sub_jets = clust_seq.subtracted_jets(range,ptmin);
173 //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
174 //cout << "---------------------------------------\n";
176 //printf(" ijet rap phi Pt area +- err\n");
178 // sort jets into increasing pt
179 vector<fastjet::PseudoJet> jets = sorted_by_pt(sub_jets);
180 for (size_t j = 0; j < jets.size(); j++) { // loop for jets
182 // double area = clust_seq.area(jets[j]);
183 // double area_error = clust_seq.area_error(jets[j]);
185 //printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),jets[j].phi(),jets[j].perp(), area, area_error);
187 // go to write AOD info
188 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
189 //cout << "Printing jet " << endl;
190 if(debug) aodjet.Print("");
191 //cout << "Adding jet ... " ;
193 //cout << "added \n" << endl;
195 } // end loop for jets
200 //____________________________________________________________________________
201 void AliFastJetFinder::RunTest(const char* datafile)
205 // This simple test run the kt algorithm for an ascii file testdata.dat
206 // read input particles -----------------------------
207 vector<fastjet::PseudoJet> input_particles;
211 // we assume a file basic.dat in the current directory
212 // this file has 3 columns of float data
215 in >> px >> py >> pz >> en;
216 if (!in.good()) break;
217 //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
219 input_particles.push_back(fastjet::PseudoJet(px,py,pz,en));
221 //printf(" found %d pointsn",nlines);
223 //////////////////////////////////////////////////
225 // create an object that represents your choice of jet algorithm, and
226 // the associated parameters
228 fastjet::Strategy strategy = fastjet::Best;
229 fastjet::RecombinationScheme recomb_scheme = fastjet::BIpt_scheme;
230 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
234 // create an object that specifies how we to define the area
235 fastjet::AreaDefinition area_def;
236 double ghost_etamax = 7.0;
237 double ghost_area = 0.05;
238 int active_area_repeats = 1;
241 // now create the object that holds info about ghosts
242 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, ghost_area);
243 // and from that get an area definition
244 area_def = fastjet::AreaDefinition(fastjet::active_area,ghost_spec);
247 // run the jet clustering with the above jet definition
248 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def);
251 // tell the user what was done
252 cout << "--------------------------------------------------------" << endl;
253 cout << "Jet definition was: " << jet_def.description() << endl;
254 cout << "Area definition was: " << area_def.description() << endl;
255 cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
256 cout << "--------------------------------------------------------" << endl;
259 // extract the inclusive jets with pt > 5 GeV, sorted by pt
261 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
263 cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
266 //subtract background // ===========================================
267 // set the rapididty range within which to study the background
268 double rap_max = ghost_etamax - Rparam;
269 fastjet::RangeDefinition range(rap_max);
270 // subtract background
271 vector<fastjet::PseudoJet> sub_jets = clust_seq.subtracted_jets(range,ptmin);
273 // print them out //================================================
274 cout << "Printing inclusive jets after background subtraction \n";
275 cout << "------------------------------------------------------\n";
276 // sort jets into increasing pt
277 vector<fastjet::PseudoJet> jets = sorted_by_pt(sub_jets);
279 printf(" ijet rap phi Pt area +- err\n");
280 for (size_t j = 0; j < jets.size(); j++) {
282 double area = clust_seq.area(jets[j]);
283 double area_error = clust_seq.area_error(jets[j]);
285 printf("%5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),
286 jets[j].phi(),jets[j].perp(), area, area_error);
289 // ================================================================
295 //____________________________________________________________________________
297 void AliFastJetFinder::WriteJHeaderToFile()
302 //____________________________________________________________________________