small corrections
[u/mrichter/AliRoot.git] / JETAN / AliFastJetFinder.cxx
CommitLineData
a17e6965 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
392c9b47 16
a17e6965 17//---------------------------------------------------------------------
392c9b47 18// FastJet v2.3.4 finder algorithm interface
8838ab7a 19// Last modification: Neutral cell energy included in the jet reconstruction
20//
21// Authors: Rafael.Diaz.Valdes@cern.ch
22// Magali.estienne@subatech.in2p3.fr (neutral part + bg subtraction option)
392c9b47 23//
a17e6965 24//---------------------------------------------------------------------
25
8838ab7a 26
a17e6965 27#include <Riostream.h>
28#include <TLorentzVector.h>
29#include <TFile.h>
30#include <TH1F.h>
31#include <TH2F.h>
32#include <TArrayF.h>
33#include <TRandom.h>
34#include <TClonesArray.h>
35
36#include "AliFastJetFinder.h"
8838ab7a 37#include "AliFastJetHeaderV1.h"
a17e6965 38#include "AliJetReaderHeader.h"
39#include "AliJetReader.h"
8838ab7a 40#include "AliJetUnitArray.h"
392c9b47 41
42#include "fastjet/PseudoJet.hh"
43#include "fastjet/ClusterSequenceArea.hh"
44#include "fastjet/AreaDefinition.hh"
45#include "fastjet/JetDefinition.hh"
46// get info on how fastjet was configured
47#include "fastjet/config.h"
48
49#ifdef ENABLE_PLUGIN_SISCONE
50#include "fastjet/SISConePlugin.hh"
51#endif
52
53#include<sstream> // needed for internal io
54#include<vector>
55#include <cmath>
56
57using namespace std;
a17e6965 58
59
9157b9c9 60ClassImp(AliFastJetFinder)
8838ab7a 61
62
392c9b47 63//____________________________________________________________________________
a17e6965 64
65AliFastJetFinder::AliFastJetFinder():
392c9b47 66 AliJetFinder()
a17e6965 67{
68 // Constructor
69}
70
392c9b47 71//____________________________________________________________________________
a17e6965 72
73AliFastJetFinder::~AliFastJetFinder()
a17e6965 74{
75 // destructor
76}
77
392c9b47 78//______________________________________________________________________________
a17e6965 79void AliFastJetFinder::FindJets()
a17e6965 80{
392c9b47 81
392c9b47 82 //pick up fastjet header
8838ab7a 83 AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
84 Bool_t debug = header->GetDebug(); // debug option
85 Bool_t bgMode = header->GetBGMode(); // choose to subtract BG or not
86 Int_t fOpt = fReader->GetReaderHeader()->GetDetector();
392c9b47 87
88 // check if we are reading AOD jets
89 TRefArray *refs = 0;
90 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
91 if (fromAod) { refs = fReader->GetReferences(); }
92
93 // RUN ALGORITHM
94 // read input particles -----------------------------
3dda898d 95 vector<fastjet::PseudoJet> inputParticles;
8838ab7a 96 if(fOpt==0)
97 {
98 TClonesArray *lvArray = fReader->GetMomentumArray();
99 if(lvArray == 0) { cout << "Could not get the momentum array" << endl; return; }
100 Int_t nIn = lvArray->GetEntries();
101 if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
8838ab7a 102 Float_t px,py,pz,en;
103 // load input vectors
104 for(Int_t i = 0; i < nIn; i++){ // loop for all input particles
105 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
106 px = lv->Px();
107 py = lv->Py();
108 pz = lv->Pz();
109 en = lv->Energy();
110
3dda898d 111 fastjet::PseudoJet inputPart(px,py,pz,en); // create PseudoJet object
112 inputPart.set_user_index(i); //label the particle into Fastjet algortihm
113 inputParticles.push_back(inputPart); // back of the inputParticles vector
8838ab7a 114 } // end loop
115 }
116 else {
117 TClonesArray* fUnit = fReader->GetUnitArray();
118 if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return; }
8838ab7a 119 Int_t nIn = fUnit->GetEntries();
120 if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
8838ab7a 121 // Information extracted from fUnitArray
122 // load input vectors and calculate total energy in array
123 Float_t pt,eta,phi,theta,px,py,pz,en;
124 Int_t ipart = 0;
125 for(Int_t i=0; i<nIn; i++)
126 {
127 AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(i);
128
129 if(uArray->GetUnitEnergy()>0.){
130
131 // It is not necessary anymore to cut on particle pt
132 pt = uArray->GetUnitEnergy();
133 eta = uArray->GetUnitEta();
134 phi = uArray->GetUnitPhi();
135 theta = EtaToTheta(eta);
136 en = (TMath::Abs(TMath::Sin(theta)) == 0) ? pt : pt/TMath::Abs(TMath::Sin(theta));
137 px = TMath::Cos(phi)*pt;
138 py = TMath::Sin(phi)*pt;
139 pz = en*TMath::TanH(eta);
140 if(debug) cout << "pt: " << pt << ", eta: " << eta << ", phi: " << phi << ", en: " << en << ", px: " << px << ", py: " << py << ", pz: " << pz << endl;
141
3dda898d 142 fastjet::PseudoJet inputPart(px,py,pz,en); // create PseudoJet object
143 inputPart.set_user_index(ipart); //label the particle into Fastjet algortihm
144 inputParticles.push_back(inputPart); // back of the inputParticles vector
8838ab7a 145 ipart++;
146 }
147 } // End loop on UnitArray
148 }
392c9b47 149
150 // create an object that represents your choice of jet algorithm, and
151 // the associated parameters
3dda898d 152 double rParam = header->GetRparam();
392c9b47 153 fastjet::Strategy strategy = header->GetStrategy();
3dda898d 154 fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
392c9b47 155 fastjet::JetAlgorithm algorithm = header->GetAlgorithm();
3dda898d 156 fastjet::JetDefinition jet_def(algorithm, rParam, recombScheme, strategy);
8838ab7a 157
392c9b47 158 // create an object that specifies how we to define the area
3dda898d 159 fastjet::AreaDefinition areaDef;
160 double ghostEtamax = header->GetGhostEtaMax();
161 double ghostArea = header->GetGhostArea();
162 int activeAreaRepeats = header->GetActiveAreaRepeats();
392c9b47 163
164 // now create the object that holds info about ghosts
3dda898d 165 fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea);
392c9b47 166 // and from that get an area definition
3dda898d 167 fastjet::AreaType areaType = header->GetAreaType();
168 areaDef = fastjet::AreaDefinition(areaType,ghost_spec);
392c9b47 169
8838ab7a 170 if(bgMode) // BG subtraction
171 {
172 //***************************** JETS FINDING AND EXTRACTION
173 // run the jet clustering with the above jet definition
3dda898d 174 fastjet::ClusterSequenceArea clust_seq(inputParticles, jet_def, areaDef);
392c9b47 175
8838ab7a 176 // save a comment in the header
177
178 TString comment = "Running FastJet algorithm with the following setup. ";
179 comment+= "Jet definition: ";
180 comment+= TString(jet_def.description());
181 comment+= ". Area definition: ";
3dda898d 182 comment+= TString(areaDef.description());
8838ab7a 183 comment+= ". Strategy adopted by FastJet: ";
184 comment+= TString(clust_seq.strategy_string());
185 header->SetComment(comment);
186 if(debug){
187 cout << "--------------------------------------------------------" << endl;
188 cout << comment << endl;
189 cout << "--------------------------------------------------------" << endl;
190 }
191 //header->PrintParameters();
192
193
194 // extract the inclusive jets with pt > ptmin, sorted by pt
195 double ptmin = header->GetPtMin();
3dda898d 196 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
8838ab7a 197
198 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
199
200
201 //subtract background // ===========================================
202 // set the rapididty , phi range within which to study the background
3dda898d 203 double rapMax = header->GetRapMax();
204 double rapMin = header->GetRapMin();
205 double phiMax = header->GetPhiMax();
206 double phiMin = header->GetPhiMin();
207 fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
8838ab7a 208
209 // subtract background
3dda898d 210 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
8838ab7a 211
212 // print out
213 //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
214 //cout << "---------------------------------------\n";
215 //cout << endl;
216 //printf(" ijet rap phi Pt area +- err\n");
217
218 // sort jets into increasing pt
3dda898d 219 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
8838ab7a 220 for (size_t j = 0; j < jets.size(); j++) { // loop for jets
221
222 double area = clust_seq.area(jets[j]);
3dda898d 223 double areaError = clust_seq.area_error(jets[j]);
8838ab7a 224
3dda898d 225 printf("Jet found %5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n", (Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp(), area, areaError);
392c9b47 226
227 // go to write AOD info
8838ab7a 228 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
229 //cout << "Printing jet " << endl;
230 if(debug) aodjet.Print("");
231 //cout << "Adding jet ... " ;
232 AddJet(aodjet);
233 //cout << "added \n" << endl;
234
235 }
236 }
237 else { // No BG subtraction
392c9b47 238
3dda898d 239 fastjet::ClusterSequence clust_seq(inputParticles, jet_def);
8838ab7a 240
241 // save a comment in the header
242
243 TString comment = "Running FastJet algorithm with the following setup. ";
244 comment+= "Jet definition: ";
245 comment+= TString(jet_def.description());
246 comment+= ". Strategy adopted by FastJet: ";
247 comment+= TString(clust_seq.strategy_string());
248 header->SetComment(comment);
249 if(debug){
250 cout << "--------------------------------------------------------" << endl;
251 cout << comment << endl;
252 cout << "--------------------------------------------------------" << endl;
253 }
254 //header->PrintParameters();
255
256 // extract the inclusive jets with pt > ptmin, sorted by pt
257 double ptmin = header->GetPtMin();
3dda898d 258 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
8838ab7a 259
260 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
392c9b47 261
3dda898d 262 vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); // Added by me
8838ab7a 263 for (size_t j = 0; j < jets.size(); j++) { // loop for jets // Added by me
264
265 printf("Jet found %5d %9.5f %8.5f %10.3f \n",(Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp());
266
267 // go to write AOD info
268 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
269 //cout << "Printing jet " << endl;
270 if(debug) aodjet.Print("");
271 //cout << "Adding jet ... " ;
272 AddJet(aodjet);
273 //cout << "added \n" << endl;
274
275 } // end loop for jets
276 }
277
a17e6965 278}
279
392c9b47 280//____________________________________________________________________________
281void AliFastJetFinder::RunTest(const char* datafile)
a17e6965 282
a17e6965 283{
a17e6965 284
392c9b47 285 // This simple test run the kt algorithm for an ascii file testdata.dat
286 // read input particles -----------------------------
3dda898d 287 vector<fastjet::PseudoJet> inputParticles;
392c9b47 288 Float_t px,py,pz,en;
289 ifstream in;
290 Int_t nlines = 0;
291 // we assume a file basic.dat in the current directory
292 // this file has 3 columns of float data
293 in.open(datafile);
294 while (1) {
295 in >> px >> py >> pz >> en;
296 if (!in.good()) break;
297 //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
298 nlines++;
3dda898d 299 inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en));
a17e6965 300 }
392c9b47 301 //printf(" found %d pointsn",nlines);
302 in.close();
303 //////////////////////////////////////////////////
304
305 // create an object that represents your choice of jet algorithm, and
306 // the associated parameters
3dda898d 307 double rParam = 1.0;
392c9b47 308 fastjet::Strategy strategy = fastjet::Best;
3dda898d 309 fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
310 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, rParam, recombScheme, strategy);
392c9b47 311
312
313
314 // create an object that specifies how we to define the area
3dda898d 315 fastjet::AreaDefinition areaDef;
316 double ghostEtamax = 7.0;
317 double ghostArea = 0.05;
318 int activeAreaRepeats = 1;
392c9b47 319
320
321 // now create the object that holds info about ghosts
3dda898d 322 fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea);
392c9b47 323 // and from that get an area definition
3dda898d 324 areaDef = fastjet::AreaDefinition(fastjet::active_area,ghost_spec);
392c9b47 325
326
327 // run the jet clustering with the above jet definition
3dda898d 328 fastjet::ClusterSequenceArea clust_seq(inputParticles, jet_def, areaDef);
392c9b47 329
330
331 // tell the user what was done
332 cout << "--------------------------------------------------------" << endl;
333 cout << "Jet definition was: " << jet_def.description() << endl;
3dda898d 334 cout << "Area definition was: " << areaDef.description() << endl;
392c9b47 335 cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
336 cout << "--------------------------------------------------------" << endl;
337
338
339 // extract the inclusive jets with pt > 5 GeV, sorted by pt
340 double ptmin = 5.0;
3dda898d 341 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
392c9b47 342
343 cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
344
345
346 //subtract background // ===========================================
347 // set the rapididty range within which to study the background
3dda898d 348 double rapMax = ghostEtamax - rParam;
349 fastjet::RangeDefinition range(rapMax);
392c9b47 350 // subtract background
3dda898d 351 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
392c9b47 352
353 // print them out //================================================
354 cout << "Printing inclusive jets after background subtraction \n";
355 cout << "------------------------------------------------------\n";
356 // sort jets into increasing pt
3dda898d 357 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
392c9b47 358
359 printf(" ijet rap phi Pt area +- err\n");
360 for (size_t j = 0; j < jets.size(); j++) {
361
362 double area = clust_seq.area(jets[j]);
3dda898d 363 double areaError = clust_seq.area_error(jets[j]);
392c9b47 364
8838ab7a 365 printf("%5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",(Int_t)j,jets[j].rap(),
3dda898d 366 jets[j].phi(),jets[j].perp(), area, areaError);
392c9b47 367 }
368 cout << endl;
369 // ================================================================
a17e6965 370
392c9b47 371
372
a17e6965 373}
374
392c9b47 375//____________________________________________________________________________
a17e6965 376
377void AliFastJetFinder::WriteJHeaderToFile()
378{
a17e6965 379 fHeader->Write();
380}
381
9157b9c9 382//____________________________________________________________________________
8838ab7a 383
384Float_t AliFastJetFinder::EtaToTheta(Float_t arg)
385{
386 // return (180./TMath::Pi())*2.*atan(exp(-arg));
387 return 2.*atan(exp(-arg));
388
389
390}
391
392//____________________________________________________________________________
393
394void AliFastJetFinder::InitTask(TChain *tree)
395{
396
397 printf("Fast jet finder initialization ******************");
398 fReader->CreateTasks(tree);
399
400}