]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliFastJetFinder.cxx
AOD Tag creator task.
[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
19//
a17e6965 20// Author: Rafael.Diaz.Valdes@cern.ch
392c9b47 21//
a17e6965 22//---------------------------------------------------------------------
23
24#include <Riostream.h>
25#include <TLorentzVector.h>
26#include <TFile.h>
27#include <TH1F.h>
28#include <TH2F.h>
29#include <TArrayF.h>
30#include <TRandom.h>
31#include <TClonesArray.h>
32
33#include "AliFastJetFinder.h"
34#include "AliFastJetHeader.h"
35#include "AliJetReaderHeader.h"
36#include "AliJetReader.h"
37#include "AliJet.h"
392c9b47 38
39
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"
46
47#ifdef ENABLE_PLUGIN_SISCONE
48#include "fastjet/SISConePlugin.hh"
49#endif
50
51#include<sstream> // needed for internal io
52#include<vector>
53#include <cmath>
54
55using namespace std;
a17e6965 56
57
58ClassImp(AliFastJetFinder);
59
60
392c9b47 61//____________________________________________________________________________
a17e6965 62
63AliFastJetFinder::AliFastJetFinder():
392c9b47 64 AliJetFinder()
a17e6965 65{
66 // Constructor
67}
68
392c9b47 69//____________________________________________________________________________
a17e6965 70
71AliFastJetFinder::~AliFastJetFinder()
a17e6965 72{
73 // destructor
74}
75
392c9b47 76//______________________________________________________________________________
a17e6965 77void AliFastJetFinder::FindJets()
a17e6965 78{
392c9b47 79
80 Bool_t debug = kFALSE;
81
82 //pick up fastjet header
83 AliFastJetHeader *header = (AliFastJetHeader*)fHeader;
84
85 // check if we are reading AOD jets
86 TRefArray *refs = 0;
87 Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
88 if (fromAod) { refs = fReader->GetReferences(); }
89
90 // RUN ALGORITHM
91 // read input particles -----------------------------
92 vector<fastjet::PseudoJet> input_particles;
a17e6965 93 TClonesArray *lvArray = fReader->GetMomentumArray();
392c9b47 94 if(lvArray == 0) { cout << "Could not get the momentum array" << endl; return; }
a17e6965 95 Int_t nIn = lvArray->GetEntries();
392c9b47 96 if(nIn == 0) { if (debug) cout << "entries = 0 ; Event empty !!!" << endl ; return; }
97 Int_t nJets = 0; // n jets in this event
98 fJets->SetNinput(nIn) ; // number of input objects
99 Float_t px,py,pz,en;
100 // load input vectors
101 for(Int_t i = 0; i < nIn; i++){ // loop for all input particles
102 TLorentzVector *lv = (TLorentzVector*) lvArray->At(i);
103 px = lv->Px();
104 py = lv->Py();
105 pz = lv->Pz();
106 en = lv->Energy();
107
108 fastjet::PseudoJet input_part(px,py,pz,en); // create PseudoJet object
109 input_part.set_user_index(i); //label the particle into Fastjet algortihm
110 input_particles.push_back(input_part); // back of the input_particles vector
111 } // end loop
112
113
114 // create an object that represents your choice of jet algorithm, and
115 // the associated parameters
116 double Rparam = header->GetRparam();
117 fastjet::Strategy strategy = header->GetStrategy();
118 fastjet::RecombinationScheme recomb_scheme = header->GetRecombScheme();
119 fastjet::JetAlgorithm algorithm = header->GetAlgorithm();
120 fastjet::JetDefinition jet_def(algorithm, Rparam, recomb_scheme, strategy);
121
122
123 // create an object that specifies how we to define the area
124 fastjet::AreaDefinition area_def;
125 double ghost_etamax = header->GetGhostEtaMax();
126 double ghost_area = header->GetGhostArea();
127 int active_area_repeats = header->GetActiveAreaRepeats();
128
129 // now create the object that holds info about ghosts
130 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, ghost_area);
131 // and from that get an area definition
132 fastjet::AreaType area_type = header->GetAreaType();
133 area_def = fastjet::AreaDefinition(area_type,ghost_spec);
134
135 // run the jet clustering with the above jet definition
136 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def);
137
138
139 // save a comment in the header
140
141 TString comment = "Running FastJet algorithm with the following setup. ";
142 comment+= "Jet definition: ";
143 comment+= TString(jet_def.description());
144 comment+= ". Area definition: ";
145 comment+= TString(area_def.description());
146 comment+= ". Strategy adopted by FastJet: ";
147 comment+= TString(clust_seq.strategy_string());
148 header->SetComment(comment);
149 if(debug){
150 cout << "--------------------------------------------------------" << endl;
151 cout << comment << endl;
152 cout << "--------------------------------------------------------" << endl;
a17e6965 153 }
392c9b47 154 //header->PrintParameters();
155
156
157 // extract the inclusive jets with pt > ptmin, sorted by pt
158 double ptmin = header->GetPtMin();
159 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
160
161 //cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
162
163
164 //subtract background // ===========================================
165 // set the rapididty , phi range within which to study the background
166 double rap_max = header->GetRapMax();
167 double rap_min = header->GetRapMin();
168 double phi_max = header->GetPhiMax();
169 double phi_min = header->GetPhiMin();
170 fastjet::RangeDefinition range(rap_min, rap_max, phi_min, phi_max);
171 // subtract background
172 vector<fastjet::PseudoJet> sub_jets = clust_seq.subtracted_jets(range,ptmin);
173
174 // print out
175 //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
176 //cout << "---------------------------------------\n";
177 //cout << endl;
178 //printf(" ijet rap phi Pt area +- err\n");
179
180 // sort jets into increasing pt
181 vector<fastjet::PseudoJet> jets = sorted_by_pt(sub_jets);
182 for (size_t j = 0; j < jets.size(); j++) { // loop for jets
183
184 double area = clust_seq.area(jets[j]);
185 double area_error = clust_seq.area_error(jets[j]);
186
187 //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);
188
189 // go to write AOD info
190 AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
191 //cout << "Printing jet " << endl;
192 if(debug) aodjet.Print("");
193 //cout << "Adding jet ... " ;
194 AddJet(aodjet);
195 //cout << "added \n" << endl;
196
197 } // end loop for jets
198
199
a17e6965 200}
201
392c9b47 202//____________________________________________________________________________
203void AliFastJetFinder::RunTest(const char* datafile)
a17e6965 204
a17e6965 205{
a17e6965 206
392c9b47 207 // This simple test run the kt algorithm for an ascii file testdata.dat
208 // read input particles -----------------------------
209 vector<fastjet::PseudoJet> input_particles;
210 Float_t px,py,pz,en;
211 ifstream in;
212 Int_t nlines = 0;
213 // we assume a file basic.dat in the current directory
214 // this file has 3 columns of float data
215 in.open(datafile);
216 while (1) {
217 in >> px >> py >> pz >> en;
218 if (!in.good()) break;
219 //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
220 nlines++;
221 input_particles.push_back(fastjet::PseudoJet(px,py,pz,en));
a17e6965 222 }
392c9b47 223 //printf(" found %d pointsn",nlines);
224 in.close();
225 //////////////////////////////////////////////////
226
227 // create an object that represents your choice of jet algorithm, and
228 // the associated parameters
229 double Rparam = 1.0;
230 fastjet::Strategy strategy = fastjet::Best;
231 fastjet::RecombinationScheme recomb_scheme = fastjet::BIpt_scheme;
232 fastjet::JetDefinition jet_def(fastjet::kt_algorithm, Rparam, recomb_scheme, strategy);
233
234
235
236 // create an object that specifies how we to define the area
237 fastjet::AreaDefinition area_def;
238 double ghost_etamax = 7.0;
239 double ghost_area = 0.05;
240 int active_area_repeats = 1;
241
242
243 // now create the object that holds info about ghosts
244 fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, ghost_area);
245 // and from that get an area definition
246 area_def = fastjet::AreaDefinition(fastjet::active_area,ghost_spec);
247
248
249 // run the jet clustering with the above jet definition
250 fastjet::ClusterSequenceArea clust_seq(input_particles, jet_def, area_def);
251
252
253 // tell the user what was done
254 cout << "--------------------------------------------------------" << endl;
255 cout << "Jet definition was: " << jet_def.description() << endl;
256 cout << "Area definition was: " << area_def.description() << endl;
257 cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
258 cout << "--------------------------------------------------------" << endl;
259
260
261 // extract the inclusive jets with pt > 5 GeV, sorted by pt
262 double ptmin = 5.0;
263 vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
264
265 cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
266
267
268 //subtract background // ===========================================
269 // set the rapididty range within which to study the background
270 double rap_max = ghost_etamax - Rparam;
271 fastjet::RangeDefinition range(rap_max);
272 // subtract background
273 vector<fastjet::PseudoJet> sub_jets = clust_seq.subtracted_jets(range,ptmin);
274
275 // print them out //================================================
276 cout << "Printing inclusive jets after background subtraction \n";
277 cout << "------------------------------------------------------\n";
278 // sort jets into increasing pt
279 vector<fastjet::PseudoJet> jets = sorted_by_pt(sub_jets);
280
281 printf(" ijet rap phi Pt area +- err\n");
282 for (size_t j = 0; j < jets.size(); j++) {
283
284 double area = clust_seq.area(jets[j]);
285 double area_error = clust_seq.area_error(jets[j]);
286
287 printf("%5u %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",j,jets[j].rap(),
288 jets[j].phi(),jets[j].perp(), area, area_error);
289 }
290 cout << endl;
291 // ================================================================
a17e6965 292
392c9b47 293
294
a17e6965 295}
296
392c9b47 297//____________________________________________________________________________
a17e6965 298
299void AliFastJetFinder::WriteJHeaderToFile()
300{
a17e6965 301 fHeader->Write();
302}
303
392c9b47 304//____________________________________________________________________________