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
97 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
98 if (fromAod) { refs = fReader->GetReferences(); }
101 // read input particles -----------------------------
103 vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
104 if(inputParticles.size()==0){
105 if(debug)Printf("%s:%d No input particles found, skipping event",(char*)__FILE__,__LINE__);
109 // create an object that represents your choice of jet algorithm, and
110 // the associated parameters
111 double rParam = header->GetRparam();
112 fastjet::Strategy strategy = header->GetStrategy();
113 fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
114 fastjet::JetAlgorithm algorithm = header->GetAlgorithm();
115 fastjet::JetDefinition jetDef(algorithm, rParam, recombScheme, strategy);
117 // create an object that specifies how we to define the area
118 fastjet::AreaDefinition areaDef;
119 double ghostEtamax = header->GetGhostEtaMax();
120 double ghostArea = header->GetGhostArea();
121 int activeAreaRepeats = header->GetActiveAreaRepeats();
123 // now create the object that holds info about ghosts
124 fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
125 // and from that get an area definition
126 fastjet::AreaType areaType = header->GetAreaType();
127 areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
133 if(bgMode) // BG subtraction!!!!!!!!!! Not entering here
135 //***************************** JETS FINDING AND EXTRACTION
136 // run the jet clustering with the above jet definition
137 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
139 // save a comment in the header
141 TString comment = "Running FastJet algorithm with the following setup. ";
142 comment+= "Jet definition: ";
143 comment+= TString(jetDef.description());
144 comment+= ". Area definition: ";
145 comment+= TString(areaDef.description());
146 comment+= ". Strategy adopted by FastJet: ";
147 comment+= TString(clust_seq.strategy_string());
148 header->SetComment(comment);
150 cout << "--------------------------------------------------------" << endl;
151 cout << comment << endl;
152 cout << "--------------------------------------------------------" << endl;
154 //header->PrintParameters();
157 // extract the inclusive jets with pt > ptmin, sorted by pt
158 double ptmin = header->GetPtMin();
159 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
161 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
164 //subtract background // ===========================================
165 // set the rapididty , phi range within which to study the background
166 double rapMax = header->GetRapMax();
167 double rapMin = header->GetRapMin();
168 double phiMax = header->GetPhiMax();
169 double phiMin = header->GetPhiMin();
170 fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
172 // subtract background
173 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
176 //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
177 //cout << "---------------------------------------\n";
179 //printf(" ijet rap phi Pt area +- err\n");
181 // sort jets into increasing pt
182 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
183 for (size_t j = 0; j < jets.size(); j++) { // loop for jets
185 double area = clust_seq.area(jets[j]);
186 double areaError = clust_seq.area_error(jets[j]);
188 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);
190 // go to write AOD info
191 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
192 //cout << "Printing jet " << endl;
193 if(debug) aodjet.Print("");
194 //cout << "Adding jet ... " ;
196 //cout << "added \n" << endl;
206 else { // No BG subtraction!!!!!!!! Default header is bgmode=0.
208 TClonesArray* fUnit = fReader->GetUnitArray(); //Big mmentum array
209 if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return; }
210 Int_t nIn = fUnit->GetEntries();
213 //fastjet::ClusterSequence clust_seq(inputParticles, jetDef);
214 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
216 // save a comment in the header
218 TString comment = "Running FastJet algorithm with the following setup. ";
219 comment+= "Jet definition: ";
220 comment+= TString(jetDef.description());
221 comment+= ". Strategy adopted by FastJet: ";
222 comment+= TString(clust_seq.strategy_string());
223 header->SetComment(comment);
225 cout << "--------------------------------------------------------" << endl;
226 cout << comment << endl;
227 cout << "--------------------------------------------------------" << endl;
229 //header->PrintParameters();
231 // extract the inclusive jets with pt > ptmin, sorted by pt
232 double ptmin = header->GetPtMin();
233 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
235 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
237 vector<fastjet::PseudoJet> jets = sorted_by_pt(inclusiveJets); // Added by me
238 for (size_t j = 0; j < jets.size(); j++) { // loop for jets // Added by me
240 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());
242 vector<fastjet::PseudoJet> constituents = clust_seq.constituents(jets[j]);
243 int nCon= constituents.size();
245 Double_t area=clust_seq.area(jets[j]);
246 // go to write AOD info
247 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
248 aodjet.SetEffArea(area,0);
249 //cout << "Printing jet " << endl;
250 if(debug) aodjet.Print("");
252 for (int i=0; i < nCon; i++)
254 fastjet::PseudoJet mPart=constituents[i];
255 ind[i]=mPart.user_index();
256 // cout<<i<<" index="<<ind[i]<<endl;
258 //internal loop over all the unit cells
261 for(Int_t ii=0; ii<nIn; ii++)
263 AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(ii);
264 if(uArray->GetUnitEnergy()>0.){
265 uArray->SetUnitTrackID(ipart);//used to have the index available in AliJEtBkg
267 TRefArray* trackArray = (TRefArray*)uArray->GetUnitTrackRef();
268 Int_t tracksInCell = trackArray->GetEntries();
269 for(int ji = 0; ji < tracksInCell; ji++){
270 AliAODTrack * track = (AliAODTrack*)trackArray->At(ji);
271 aodjet.AddTrack(track);
285 } // end loop for jets
290 //____________________________________________________________________________
291 void AliFastJetFinder::RunTest(const char* datafile)
295 // This simple test run the kt algorithm for an ascii file testdata.dat
296 // read input particles -----------------------------
297 vector<fastjet::PseudoJet> inputParticles;
301 // we assume a file basic.dat in the current directory
302 // this file has 3 columns of float data
305 in >> px >> py >> pz >> en;
306 if (!in.good()) break;
307 //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
309 inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en));
311 //printf(" found %d pointsn",nlines);
313 //////////////////////////////////////////////////
315 // create an object that represents your choice of jet algorithm, and
316 // the associated parameters
318 fastjet::Strategy strategy = fastjet::Best;
319 fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
320 fastjet::JetDefinition jetDef(fastjet::kt_algorithm, rParam, recombScheme, strategy);
324 // create an object that specifies how we to define the area
325 fastjet::AreaDefinition areaDef;
326 double ghostEtamax = 7.0;
327 double ghostArea = 0.05;
328 int activeAreaRepeats = 1;
331 // now create the object that holds info about ghosts
332 fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
333 // and from that get an area definition
334 areaDef = fastjet::AreaDefinition(fastjet::active_area,ghostSpec);
337 // run the jet clustering with the above jet definition
338 fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
341 // tell the user what was done
342 cout << "--------------------------------------------------------" << endl;
343 cout << "Jet definition was: " << jetDef.description() << endl;
344 cout << "Area definition was: " << areaDef.description() << endl;
345 cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
346 cout << "--------------------------------------------------------" << endl;
349 // extract the inclusive jets with pt > 5 GeV, sorted by pt
351 vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
353 cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
356 //subtract background // ===========================================
357 // set the rapididty range within which to study the background
358 double rapMax = ghostEtamax - rParam;
359 fastjet::RangeDefinition range(rapMax);
360 // subtract background
361 vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
363 // print them out //================================================
364 cout << "Printing inclusive jets after background subtraction \n";
365 cout << "------------------------------------------------------\n";
366 // sort jets into increasing pt
367 vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);
369 printf(" ijet rap phi Pt area +- err\n");
370 for (size_t j = 0; j < jets.size(); j++) {
372 double area = clust_seq.area(jets[j]);
373 double areaError = clust_seq.area_error(jets[j]);
375 printf("%5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",(Int_t)j,jets[j].rap(),
376 jets[j].phi(),jets[j].perp(), area, areaError);
379 // ================================================================
385 //____________________________________________________________________________
387 void AliFastJetFinder::WriteJHeaderToFile() const
392 //____________________________________________________________________________
394 Float_t AliFastJetFinder::EtaToTheta(Float_t arg)
396 // return (180./TMath::Pi())*2.*atan(exp(-arg));
397 return 2.*atan(exp(-arg));
402 //____________________________________________________________________________
404 void AliFastJetFinder::InitTask(TChain *tree)
407 printf("Fast jet finder initialization ******************");
408 fReader->CreateTasks(tree);
413 Bool_t AliFastJetFinder::ProcessEvent()
417 // from meomntum array
419 Bool_t ok = fReader->FillMomentumArray();
421 if (!ok) return kFALSE;
422 fInputFJ->SetHeader(fHeader);
423 fInputFJ->SetReader(fReader);
424 fInputFJ->FillInput();
428 fJetBkg->SetHeader(fHeader);
429 fJetBkg->SetReader(fReader);
431 fJetBkg->SetFastJetInput(fInputFJ);
432 Double_t bkg1=fJetBkg->BkgFastJet();
433 Double_t bkg2=fJetBkg->BkgChargedFastJet();
434 Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
435 Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
437 fAODEvBkg->SetBackground(0,bkg1);
438 fAODEvBkg->SetBackground(1,bkg2);
439 fAODEvBkg->SetBackground(2,bkg3);
440 fAODEvBkg->SetBackground(3,bkg4);
447 Bool_t AliFastJetFinder::ProcessEvent2()
451 // Charged only or charged+neutral jets
454 TRefArray* ref = new TRefArray();
455 Bool_t procid = kFALSE;
456 Bool_t ok = fReader->ExecTasks(procid,ref);
458 // Delete reference pointer
459 if (!ok) {delete ref; return kFALSE;}
462 fInputFJ->SetHeader(fHeader);
463 fInputFJ->SetReader(fReader);
464 fInputFJ->FillInput();
469 fJetBkg->SetHeader(fHeader);
470 fJetBkg->SetReader(fReader);
471 fJetBkg->SetFastJetInput(fInputFJ);
472 // Double_t bkg1=fJetBkg->BkgFastJet();
473 // Double_t bkg2=fJetBkg->BkgChargedFastJet();
474 // Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
475 // Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
477 // fAODEvBkg->SetBackground(0,bkg1);
478 // fAODEvBkg->SetBackground(1,bkg2);
479 // fAODEvBkg->SetBackground(2,bkg3);
480 // fAODEvBkg->SetBackground(3,bkg4);
482 Int_t nEntRef = ref->GetEntries();
484 for(Int_t i=0; i<nEntRef; i++)
486 // Reset the UnitArray content which were referenced
487 ((AliJetUnitArray*)ref->At(i))->SetUnitTrackID(0);
488 ((AliJetUnitArray*)ref->At(i))->SetUnitEnergy(0.);
489 ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag(kPtSmaller);
490 ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag2(kPtSmaller);
491 ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlag(kBad);
492 ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlagC(kTRUE,kBad);
493 ((AliJetUnitArray*)ref->At(i))->SetUnitDetectorFlag(kTpc);
494 ((AliJetUnitArray*)ref->At(i))->SetUnitFlag(kOutJet);
495 ((AliJetUnitArray*)ref->At(i))->ClearUnitTrackRef();
498 AliJetUnitArray* uA = (AliJetUnitArray*)ref->At(i);
499 uA->ResetBit(kIsReferenced);
503 // Delete the reference pointer