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
19 // Last modification: Neutral cell energy included in the jet reconstruction
21 // Authors: Rafael.Diaz.Valdes@cern.ch
22 // Magali.estienne@subatech.in2p3.fr (neutral part + bg subtraction option)
24 //---------------------------------------------------------------------
27 #include <Riostream.h>
28 #include <TLorentzVector.h>
33 #include <TClonesArray.h>
35 #include "AliFastJetFinder.h"
36 #include "AliFastJetHeaderV1.h"
37 #include "AliJetReaderHeader.h"
38 #include "AliJetReader.h"
39 #include "AliJetUnitArray.h"
40 #include "AliFastJetInput.h"
41 #include "AliJetBkg.h"
42 #include "AliAODJetEventBackground.h"
43 #include "AliAODTrack.h"
45 #include "fastjet/PseudoJet.hh"
46 #include "fastjet/ClusterSequenceArea.hh"
47 #include "fastjet/AreaDefinition.hh"
48 #include "fastjet/JetDefinition.hh"
49 // get info on how fastjet was configured
50 #include "fastjet/config.h"
52 #ifdef ENABLE_PLUGIN_SISCONE
53 #include "fastjet/SISConePlugin.hh"
56 #include<sstream> // needed for internal io
63 ClassImp(AliFastJetFinder)
66 //____________________________________________________________________________
68 AliFastJetFinder::AliFastJetFinder():
70 fInputFJ(new AliFastJetInput()),
71 fJetBkg(new AliJetBkg())
76 //____________________________________________________________________________
78 AliFastJetFinder::~AliFastJetFinder()
81 delete fInputFJ; fInputFJ = 0;
82 delete fJetBkg; fJetBkg = 0;
85 //______________________________________________________________________________
86 void AliFastJetFinder::FindJets()
89 //pick up fastjet header
90 AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
91 Int_t debug = header->GetDebug(); // debug option
92 Bool_t bgMode = header->GetBGMode(); // choose to subtract BG or not
93 if(debug)cout<<"----------in AliFastJetFinder::FindJets() ------------------"<<endl;
95 // check if we are reading AOD jets
98 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
99 if (fromAod) { refs = fReader->GetReferences(); }
102 // read input particles -----------------------------
104 vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
105 if(inputParticles.size()==0){
106 if(debug)Printf("%s:%d No input particles found, skipping event",(char*)__FILE__,__LINE__);
110 // create an object that represents your choice of jet algorithm, and
111 // the associated parameters
112 double rParam = header->GetRparam();
113 fastjet::Strategy strategy = header->GetStrategy();
114 fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
115 fastjet::JetAlgorithm algorithm = header->GetAlgorithm();
116 fastjet::JetDefinition jetDef(algorithm, rParam, recombScheme, strategy);
118 // create an object that specifies how we to define the area
119 fastjet::AreaDefinition areaDef;
120 double ghostEtamax = header->GetGhostEtaMax();
121 double ghostArea = header->GetGhostArea();
122 int activeAreaRepeats = header->GetActiveAreaRepeats();
124 // now create the object that holds info about ghosts
125 fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
126 // and from that get an area definition
127 fastjet::AreaType areaType = header->GetAreaType();
128 areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
134 if(bgMode) // BG subtraction!!!!!!!!!! Not entering here
136 //***************************** JETS FINDING AND EXTRACTION
137 // run the jet clustering with the above jet definition
138 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
140 // save a comment in the header
142 TString comment = "Running FastJet algorithm with the following setup. ";
143 comment+= "Jet definition: ";
144 comment+= TString(jetDef.description());
145 comment+= ". Area definition: ";
146 comment+= TString(areaDef.description());
147 comment+= ". Strategy adopted by FastJet: ";
148 comment+= TString(clust_seq.strategy_string());
149 header->SetComment(comment);
151 cout << "--------------------------------------------------------" << endl;
152 cout << comment << endl;
153 cout << "--------------------------------------------------------" << endl;
155 //header->PrintParameters();
158 // extract the inclusive jets with pt > ptmin, sorted by pt
159 double ptmin = header->GetPtMin();
160 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
162 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
165 //subtract background // ===========================================
166 // set the rapididty , phi range within which to study the background
167 double rapMax = header->GetRapMax();
168 double rapMin = header->GetRapMin();
169 double phiMax = header->GetPhiMax();
170 double phiMin = header->GetPhiMin();
171 fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
173 // subtract background
174 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
177 //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
178 //cout << "---------------------------------------\n";
180 //printf(" ijet rap phi Pt area +- err\n");
182 // sort jets into increasing pt
183 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
184 for (size_t j = 0; j < jets.size(); j++) { // loop for jets
186 double area = clust_seq.area(jets[j]);
187 double areaError = clust_seq.area_error(jets[j]);
189 if(debug) 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);
191 // go to write AOD info
192 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
193 //cout << "Printing jet " << endl;
194 if(debug) aodjet.Print("");
195 //cout << "Adding jet ... " ;
197 //cout << "added \n" << endl;
207 else { // No BG subtraction!!!!!!!! Default header is bgmode=0.
209 TClonesArray* fUnit = fReader->GetUnitArray(); //Big mmentum array
210 if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return; }
211 Int_t nIn = fUnit->GetEntries();
214 //fastjet::ClusterSequence clust_seq(inputParticles, jetDef);
215 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
217 // save a comment in the header
219 TString comment = "Running FastJet algorithm with the following setup. ";
220 comment+= "Jet definition: ";
221 comment+= TString(jetDef.description());
222 comment+= ". Strategy adopted by FastJet: ";
223 comment+= TString(clust_seq.strategy_string());
224 header->SetComment(comment);
226 cout << "--------------------------------------------------------" << endl;
227 cout << comment << endl;
228 cout << "--------------------------------------------------------" << endl;
230 //header->PrintParameters();
232 // extract the inclusive jets with pt > ptmin, sorted by pt
233 double ptmin = header->GetPtMin();
234 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
236 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
238 vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); // Added by me
239 for (size_t j = 0; j < jets.size(); j++) { // loop for jets // Added by me
241 if(debug) printf("Jet found %5d %9.5f %8.5f %10.3f \n",(Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp());
243 vector<fastjet::PseudoJet> constituents = clust_seq.constituents(jets[j]);
244 int nCon= constituents.size();
246 Double_t area=clust_seq.area(jets[j]);
247 // go to write AOD info
248 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
249 aodjet.SetEffArea(area,0);
250 //cout << "Printing jet " << endl;
251 if(debug) aodjet.Print("");
253 for (int i=0; i < nCon; i++)
255 fastjet::PseudoJet mPart=constituents[i];
256 ind[i]=mPart.user_index();
257 // cout<<i<<" index="<<ind[i]<<endl;
259 //internal loop over all the unit cells
262 for(Int_t ii=0; ii<nIn; ii++)
264 AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(ii);
265 if(uArray->GetUnitEnergy()>0.){
266 uArray->SetUnitTrackID(ipart);//used to have the index available in AliJEtBkg
268 TRefArray* trackArray = (TRefArray*)uArray->GetUnitTrackRef();
269 Int_t tracksInCell = trackArray->GetEntries();
270 for(int ji = 0; ji < tracksInCell; ji++){
271 AliAODTrack * track = (AliAODTrack*)trackArray->At(ji);
272 aodjet.AddTrack(track);
286 } // end loop for jets
291 //____________________________________________________________________________
292 void AliFastJetFinder::RunTest(const char* datafile)
296 // This simple test run the kt algorithm for an ascii file testdata.dat
297 // read input particles -----------------------------
298 vector<fastjet::PseudoJet> inputParticles;
302 // we assume a file basic.dat in the current directory
303 // this file has 3 columns of float data
306 in >> px >> py >> pz >> en;
307 if (!in.good()) break;
308 //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
310 inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en));
312 //printf(" found %d pointsn",nlines);
314 //////////////////////////////////////////////////
316 // create an object that represents your choice of jet algorithm, and
317 // the associated parameters
319 fastjet::Strategy strategy = fastjet::Best;
320 fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
321 fastjet::JetDefinition jetDef(fastjet::kt_algorithm, rParam, recombScheme, strategy);
325 // create an object that specifies how we to define the area
326 fastjet::AreaDefinition areaDef;
327 double ghostEtamax = 7.0;
328 double ghostArea = 0.05;
329 int activeAreaRepeats = 1;
332 // now create the object that holds info about ghosts
333 fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
334 // and from that get an area definition
335 areaDef = fastjet::AreaDefinition(fastjet::active_area,ghostSpec);
338 // run the jet clustering with the above jet definition
339 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
342 // tell the user what was done
343 cout << "--------------------------------------------------------" << endl;
344 cout << "Jet definition was: " << jetDef.description() << endl;
345 cout << "Area definition was: " << areaDef.description() << endl;
346 cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
347 cout << "--------------------------------------------------------" << endl;
350 // extract the inclusive jets with pt > 5 GeV, sorted by pt
352 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
354 cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
357 //subtract background // ===========================================
358 // set the rapididty range within which to study the background
359 double rapMax = ghostEtamax - rParam;
360 fastjet::RangeDefinition range(rapMax);
361 // subtract background
362 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
364 // print them out //================================================
365 cout << "Printing inclusive jets after background subtraction \n";
366 cout << "------------------------------------------------------\n";
367 // sort jets into increasing pt
368 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
370 printf(" ijet rap phi Pt area +- err\n");
371 for (size_t j = 0; j < jets.size(); j++) {
373 double area = clust_seq.area(jets[j]);
374 double areaError = clust_seq.area_error(jets[j]);
376 printf("%5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",(Int_t)j,jets[j].rap(),
377 jets[j].phi(),jets[j].perp(), area, areaError);
380 // ================================================================
386 //____________________________________________________________________________
388 void AliFastJetFinder::WriteJHeaderToFile() const
393 //____________________________________________________________________________
395 Float_t AliFastJetFinder::EtaToTheta(Float_t arg)
397 // return (180./TMath::Pi())*2.*atan(exp(-arg));
398 return 2.*atan(exp(-arg));
403 //____________________________________________________________________________
405 void AliFastJetFinder::InitTask(TChain *tree)
408 printf("Fast jet finder initialization ******************");
409 fReader->CreateTasks(tree);
414 Bool_t AliFastJetFinder::ProcessEvent()
418 // from meomntum array
420 Bool_t ok = fReader->FillMomentumArray();
422 if (!ok) return kFALSE;
423 fInputFJ->SetHeader(fHeader);
424 fInputFJ->SetReader(fReader);
425 fInputFJ->FillInput();
431 fJetBkg->SetHeader(fHeader);
432 fJetBkg->SetReader(fReader);
433 Double_t sigma1 = 0,meanarea1= 0,sigma2 = 0,meanarea2 = 0;
434 Double_t bkg1 = 0,bkg2 = 0;
436 fJetBkg->SetFastJetInput(fInputFJ);
437 fJetBkg->BkgFastJetb(bkg1,sigma1,meanarea1);
438 fJetBkg->BkgFastJetWoHardest(bkg2,sigma2,meanarea2);
439 fAODEvBkg->SetBackground(0,bkg1,sigma1,meanarea1);
440 fAODEvBkg->SetBackground(1,bkg2,sigma2,meanarea2);
447 fJetBkg->SetFastJetInput(fInputFJ);
448 Double_t bkg1=fJetBkg->BkgFastJet();
449 Double_t bkg2=fJetBkg->BkgChargedFastJet();
450 Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
451 Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
453 fAODEvBkg->SetBackground(0,bkg1);
454 fAODEvBkg->SetBackground(1,bkg2);
455 fAODEvBkg->SetBackground(2,bkg3);
456 fAODEvBkg->SetBackground(3,bkg4);
463 Bool_t AliFastJetFinder::ProcessEvent2()
467 // Charged only or charged+neutral jets
470 TRefArray* ref = new TRefArray();
471 Bool_t procid = kFALSE;
472 Bool_t ok = fReader->ExecTasks(procid,ref);
474 // Delete reference pointer
475 if (!ok) {delete ref; return kFALSE;}
478 fInputFJ->SetHeader(fHeader);
479 fInputFJ->SetReader(fReader);
480 fInputFJ->FillInput();
486 fJetBkg->SetHeader(fHeader);
487 fJetBkg->SetReader(fReader);
488 fJetBkg->SetFastJetInput(fInputFJ);
489 Double_t sigma1,meanarea1,sigma2,meanarea2;
491 fJetBkg->BkgFastJetb(bkg1,sigma1,meanarea1);
492 fJetBkg->BkgFastJetWoHardest(bkg2,sigma2,meanarea2);
493 fAODEvBkg->SetBackground(0,bkg1,sigma1,meanarea1);
494 fAODEvBkg->SetBackground(1,bkg2,sigma2,meanarea2);
498 // Double_t bkg1=fJetBkg->BkgFastJet();
499 // Double_t bkg2=fJetBkg->BkgChargedFastJet();
500 // Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
501 // Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
503 // fAODEvBkg->SetBackground(0,bkg1);
504 // fAODEvBkg->SetBackground(1,bkg2);
505 // fAODEvBkg->SetBackground(2,bkg3);
506 // fAODEvBkg->SetBackground(3,bkg4);
508 Int_t nEntRef = ref->GetEntries();
510 for(Int_t i=0; i<nEntRef; i++)
512 // Reset the UnitArray content which were referenced
513 ((AliJetUnitArray*)ref->At(i))->SetUnitTrackID(0);
514 ((AliJetUnitArray*)ref->At(i))->SetUnitEnergy(0.);
515 ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag(kPtSmaller);
516 ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag2(kPtSmaller);
517 ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlag(kBad);
518 ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlagC(kTRUE,kBad);
519 ((AliJetUnitArray*)ref->At(i))->SetUnitDetectorFlag(kTpc);
520 ((AliJetUnitArray*)ref->At(i))->SetUnitFlag(kOutJet);
521 ((AliJetUnitArray*)ref->At(i))->ClearUnitTrackRef();
524 AliJetUnitArray* uA = (AliJetUnitArray*)ref->At(i);
525 uA->ResetBit(kIsReferenced);
529 // Delete the reference pointer