]> git.uio.no Git - u/mrichter/AliRoot.git/blob - JETAN/AliFastJetFinder.cxx
changed order of adding histograms to the TList
[u/mrichter/AliRoot.git] / JETAN / AliFastJetFinder.cxx
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  
16
17 //---------------------------------------------------------------------
18 // FastJet v2.3.4 finder algorithm interface
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)
23 //
24 //---------------------------------------------------------------------
25
26
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 <TClonesArray.h>
34
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"
44
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"
51
52 #ifdef ENABLE_PLUGIN_SISCONE
53 #include "fastjet/SISConePlugin.hh"
54 #endif
55
56 #include<sstream>  // needed for internal io
57 #include<vector> 
58 #include <cmath> 
59
60 using namespace std;
61
62
63 ClassImp(AliFastJetFinder)
64
65
66 //____________________________________________________________________________
67
68 AliFastJetFinder::AliFastJetFinder():
69   AliJetFinder(),
70   fInputFJ(new AliFastJetInput()),
71   fJetBkg(new  AliJetBkg())
72 {
73   // Constructor
74 }
75
76 //____________________________________________________________________________
77
78 AliFastJetFinder::~AliFastJetFinder()
79 {
80   // destructor
81   delete  fInputFJ; fInputFJ = 0;
82   delete  fJetBkg; fJetBkg = 0;
83 }
84
85 //______________________________________________________________________________
86 void AliFastJetFinder::FindJets()
87 {
88
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;
94
95   // check if we are reading AOD jets
96   TRefArray *refs = 0;
97   Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
98   if (fromAod) { refs = fReader->GetReferences(); }
99
100   // RUN ALGORITHM  
101   // read input particles -----------------------------
102
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__);
106     return;
107   }
108
109   // create an object that represents your choice of jet algorithm, and 
110   // the associated parameters
111   double rParam = header->GetRparam();
112   double rBkgParam = header->GetRparamBkg();
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);
117
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(); 
123   
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);
129   
130   //***************************** JETS FINDING
131   // run the jet clustering with the above jet definition
132   fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
133
134   vector<fastjet::PseudoJet> jets;
135
136   if(bgMode) // Do BG subtraction directly with the same algorithm (cambridge or kt) for jet signal and background
137     {
138       //***************************** CLUSTER JETS FINDING FOR RHO ESTIMATION
139       // run the jet clustering with the above jet definition
140           fastjet::JetAlgorithm algorithmBkg = header->GetBGAlgorithm();
141           fastjet::JetDefinition jetDefBkg(algorithmBkg, rBkgParam, recombScheme, strategy);
142       fastjet::ClusterSequenceArea clust_seq_bkg(inputParticles, jetDefBkg, areaDef);
143
144       // save a comment in the header
145       TString comment = "Running FastJet algorithm with the following setup. ";
146       comment+= "Jet definition: ";
147       comment+= TString(jetDef.description());
148       comment+= "Jet bckg definition: ";
149       comment+= TString(jetDefBkg.description());
150       comment+= ". Area definition: ";
151       comment+= TString(areaDef.description());
152       comment+= ". Strategy adopted by FastJet and bkg: ";
153       comment+= TString(clust_seq.strategy_string());
154       header->SetComment(comment);
155       if(debug){
156           cout << "--------------------------------------------------------" << endl;
157           cout << comment << endl;
158           cout << "--------------------------------------------------------" << endl;
159       }
160       //header->PrintParameters();
161       
162       // extract the inclusive jets with pt > ptmin, sorted by pt
163       double ptmin = header->GetPtMin(); 
164       vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets();
165
166       //subtract background // ===========================================
167       // set the rapididty , phi range within which to study the background 
168       double rapMax = header->GetRapMax(); 
169       double rapMin = header->GetRapMin();
170       double phiMax = header->GetPhiMax();
171       double phiMin = header->GetPhiMin();
172       fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
173       
174       // Extract rho and sigma
175       Double_t rho = 0.;
176       Double_t sigma = 0.;
177       Double_t meanarea = 0.;
178       Bool_t Use4VectorArea = header->Use4VectorArea();
179       vector<fastjet::PseudoJet> bkgJets = clust_seq_bkg.inclusive_jets();
180       clust_seq_bkg.get_median_rho_and_sigma(bkgJets,range, Use4VectorArea, rho, sigma, meanarea, false);
181
182           // subtract background and extract jets bkg subtracted
183       vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(rho,ptmin);
184
185           // print out
186       //cout << "Printing inclusive sub jets with pt > "<< ptmin<<" GeV\n";
187       //cout << "---------------------------------------\n";
188       //cout << endl;
189       //printf(" ijet   rap      phi        Pt         area  +-   err\n");
190       
191       // sort jets into increasing pt
192       jets = sorted_by_pt(subJets);  
193
194     }
195   else { // No BG subtraction!!!!!!!! Default header is bgmode=0.
196
197           // save a comment in the header
198           TString comment = "Running FastJet algorithm with the following setup. ";
199           comment+= "Jet definition: ";
200           comment+= TString(jetDef.description());
201           comment+= ". Strategy adopted by FastJet: ";
202           comment+= TString(clust_seq.strategy_string());
203           header->SetComment(comment);
204           if(debug){
205                   cout << "--------------------------------------------------------" << endl;
206                   cout << comment << endl;
207                   cout << "--------------------------------------------------------" << endl;
208           }
209           //header->PrintParameters();
210
211           // extract the inclusive jets with pt > ptmin, sorted by pt
212           double ptmin = header->GetPtMin();
213           vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
214
215           jets = sorted_by_pt(inclusiveJets); // Added by me
216
217   }
218
219   for (size_t j = 0; j < jets.size(); j++) { // loop for jets
220
221           double area      = clust_seq.area(jets[j]);
222           double areaError = clust_seq.area_error(jets[j]);
223
224           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);
225
226       vector<fastjet::PseudoJet> constituents = clust_seq.constituents(jets[j]);
227       int nCon= constituents.size();
228       TArrayI ind(nCon);
229       
230       if ((jets[j].eta() > (header->GetJetEtaMax())) ||
231           (jets[j].eta() < (header->GetJetEtaMin())) ||
232           (jets[j].phi() > (header->GetJetPhiMax())) ||
233           (jets[j].phi() < (header->GetJetPhiMin())) ||
234           (jets[j].perp() < header->GetPtMin())) continue; // acceptance eta range and etmin
235
236           // go to write AOD  info
237           AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
238       aodjet.SetEffArea(area,areaError);
239           //cout << "Printing jet " << endl;
240           if(debug) aodjet.Print("");
241
242           Int_t count=0;
243           for (int i=0; i < nCon; i++)
244           {
245                   fastjet::PseudoJet mPart=constituents[i];
246                   ind[i]=mPart.user_index();
247                   //    cout<<i<<"  index="<<ind[i]<<endl;
248
249                   if(fromAod)
250                   {
251                           if(fReader->GetReaderHeader()->GetDetector()==0)
252                           {
253                                   for (Int_t iref = 0; iref < refs->GetEntries(); iref++)
254                                   {
255                                           if(iref==ind[i]){
256                                                   AliAODTrack * track = (AliAODTrack*)refs->At(iref);
257                                                   aodjet.AddTrack(track);
258                                           }
259                                   }
260
261                           } else {
262
263                                   TClonesArray* fUnit = fReader->GetUnitArray(); //Big mmentum array
264                                   if(fUnit == 0) { cout << "Could not get the momentum array" << endl; return; }
265                                   Int_t         nIn = fUnit->GetEntries();
266
267                                   //internal loop over all the unit cells
268                                   Int_t ipart = 0;
269
270                                   for(Int_t ii=0; ii<nIn; ii++)
271                                   {
272                                           AliJetUnitArray *uArray = (AliJetUnitArray*)fUnit->At(ii);
273                                           if(uArray->GetUnitEnergy()>0.){
274                                                   uArray->SetUnitTrackID(ipart);//used to have the index available in AliJEtBkg
275                                                   if(ipart==ind[i]){
276                                                           TRefArray* trackArray = (TRefArray*)uArray->GetUnitTrackRef();
277                                                           Int_t tracksInCell = trackArray->GetEntries();
278                                                           for(int ji = 0; ji < tracksInCell; ji++){
279                                                                   AliAODTrack * track = (AliAODTrack*)trackArray->At(ji);
280                                                                   aodjet.AddTrack(track);
281                                                           }
282
283                                                           count++;
284                                                   }
285                                                   ipart++;
286                                           }
287                                   }
288                           }
289                   }
290          } 
291
292  AddJet(aodjet);
293
294
295   } // End loop on jets
296
297 }
298
299 //____________________________________________________________________________
300 void AliFastJetFinder::RunTest(const char* datafile)
301
302 {
303
304    // This simple test run the kt algorithm for an ascii file testdata.dat
305    // read input particles -----------------------------
306   vector<fastjet::PseudoJet> inputParticles;
307   Float_t px,py,pz,en;
308   ifstream in;
309   Int_t nlines = 0;
310   // we assume a file basic.dat in the current directory
311   // this file has 3 columns of float data
312   in.open(datafile);
313   while (1) {
314       in >> px >> py >> pz >> en;
315       if (!in.good()) break;
316       //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
317       nlines++;
318       inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en)); 
319    }
320    //printf(" found %d pointsn",nlines);
321    in.close();
322    //////////////////////////////////////////////////
323  
324   // create an object that represents your choice of jet algorithm, and 
325   // the associated parameters
326   double rParam = 1.0;
327   fastjet::Strategy strategy = fastjet::Best;
328   fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
329   fastjet::JetDefinition jetDef(fastjet::kt_algorithm, rParam, recombScheme, strategy);
330   
331   
332  
333   // create an object that specifies how we to define the area
334   fastjet::AreaDefinition areaDef;
335   double ghostEtamax = 7.0;
336   double ghostArea    = 0.05;
337   int    activeAreaRepeats = 1;
338   
339
340   // now create the object that holds info about ghosts
341   fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
342   // and from that get an area definition
343   areaDef = fastjet::AreaDefinition(fastjet::active_area,ghostSpec);
344   
345
346   // run the jet clustering with the above jet definition
347   fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
348   
349   
350   // tell the user what was done
351   cout << "--------------------------------------------------------" << endl;
352   cout << "Jet definition was: " << jetDef.description() << endl;
353   cout << "Area definition was: " << areaDef.description() << endl;
354   cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
355   cout << "--------------------------------------------------------" << endl;
356  
357   
358   // extract the inclusive jets with pt > 5 GeV, sorted by pt
359   double ptmin = 5.0;
360   vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
361   
362   cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
363  
364  
365   //subtract background // ===========================================
366   // set the rapididty range within which to study the background 
367   double rapMax = ghostEtamax - rParam;
368   fastjet::RangeDefinition range(rapMax);
369   // subtract background
370   vector<fastjet::PseudoJet> subJets =  clust_seq.subtracted_jets(range,ptmin);  
371   
372   // print them out //================================================
373   cout << "Printing inclusive jets  after background subtraction \n";
374   cout << "------------------------------------------------------\n";
375   // sort jets into increasing pt
376   vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);  
377
378   printf(" ijet   rap      phi        Pt         area  +-   err\n");
379   for (size_t j = 0; j < jets.size(); j++) {
380
381     double area     = clust_seq.area(jets[j]);
382     double areaError = clust_seq.area_error(jets[j]);
383
384     printf("%5d %9.5f %8.5f %10.3f %8.3f +- %6.3f\n",(Int_t)j,jets[j].rap(),
385            jets[j].phi(),jets[j].perp(), area, areaError);
386   }
387   cout << endl;
388   // ================================================================
389
390   
391  
392 }
393
394 //____________________________________________________________________________
395
396 void AliFastJetFinder::WriteJHeaderToFile() const
397 {
398   fHeader->Write();
399 }
400
401 //____________________________________________________________________________
402
403 Float_t  AliFastJetFinder::EtaToTheta(Float_t arg)
404 {
405   //  return (180./TMath::Pi())*2.*atan(exp(-arg));
406   return 2.*atan(exp(-arg));
407
408
409 }
410
411 //____________________________________________________________________________
412
413 void AliFastJetFinder::InitTask(TChain *tree)
414 {
415
416   printf("Fast jet finder initialization ******************");
417   fReader->CreateTasks(tree);
418
419 }
420
421
422 Bool_t AliFastJetFinder::ProcessEvent()
423 {
424   //
425   // Process one event
426   // from meomntum array
427
428   Bool_t ok = fReader->FillMomentumArray();
429
430   if (!ok) return kFALSE;
431   fInputFJ->SetHeader(fHeader);
432   fInputFJ->SetReader(fReader);
433   fInputFJ->FillInput();
434   // Jets
435   FindJets(); 
436
437
438   if( fAODEvBkg){
439     fJetBkg->SetHeader(fHeader);
440     fJetBkg->SetReader(fReader);
441     Double_t sigma1 = 0,meanarea1= 0,sigma2 = 0,meanarea2 = 0;
442     Double_t bkg1 = 0,bkg2 = 0;
443     
444     fJetBkg->SetFastJetInput(fInputFJ);
445     fJetBkg->BkgFastJetb(bkg1,sigma1,meanarea1);
446     fJetBkg->BkgFastJetWoHardest(bkg2,sigma2,meanarea2);
447     fAODEvBkg->SetBackground(0,bkg1,sigma1,meanarea1);
448     fAODEvBkg->SetBackground(1,bkg2,sigma2,meanarea2);
449   }
450
451
452
453
454   /*
455   fJetBkg->SetFastJetInput(fInputFJ);
456   Double_t bkg1=fJetBkg->BkgFastJet();
457   Double_t bkg2=fJetBkg->BkgChargedFastJet();
458   Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
459   Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
460   
461   fAODEvBkg->SetBackground(0,bkg1);
462   fAODEvBkg->SetBackground(1,bkg2);
463   fAODEvBkg->SetBackground(2,bkg3);
464   fAODEvBkg->SetBackground(3,bkg4);
465   */
466   Reset();  
467   return kTRUE;
468
469 }
470
471 Bool_t AliFastJetFinder::ProcessEvent2()
472 {
473   //
474   // Process one event
475   // Charged only or charged+neutral jets
476   //
477
478   TRefArray* ref = new TRefArray();
479   Bool_t procid = kFALSE;
480   Bool_t ok = fReader->ExecTasks(procid,ref);
481
482   // Delete reference pointer  
483   if (!ok) {delete ref; return kFALSE;}
484   
485   // Leading particles
486   fInputFJ->SetHeader(fHeader);
487   fInputFJ->SetReader(fReader);
488   fInputFJ->FillInput();
489   
490   // Jets
491   FindJets();
492   
493   if( fAODEvBkg){
494     fJetBkg->SetHeader(fHeader);
495     fJetBkg->SetReader(fReader);
496     fJetBkg->SetFastJetInput(fInputFJ);
497     Double_t sigma1,meanarea1,sigma2,meanarea2;
498     Double_t bkg1,bkg2;
499     fJetBkg->BkgFastJetb(bkg1,sigma1,meanarea1);
500     fJetBkg->BkgFastJetWoHardest(bkg2,sigma2,meanarea2);
501     fAODEvBkg->SetBackground(0,bkg1,sigma1,meanarea1);
502     fAODEvBkg->SetBackground(1,bkg2,sigma2,meanarea2);
503   }
504
505
506 //  Double_t bkg1=fJetBkg->BkgFastJet();
507 //  Double_t bkg2=fJetBkg->BkgChargedFastJet();
508 //  Double_t bkg3=fJetBkg->BkgFastJetCone(fAODjets);
509 //  Double_t bkg4=fJetBkg->BkgRemoveJetLeading(fAODjets);
510   
511 //  fAODEvBkg->SetBackground(0,bkg1);
512 //  fAODEvBkg->SetBackground(1,bkg2);
513 //  fAODEvBkg->SetBackground(2,bkg3);
514 //  fAODEvBkg->SetBackground(3,bkg4);
515   
516   Int_t nEntRef    = ref->GetEntries();
517
518   for(Int_t i=0; i<nEntRef; i++)
519     { 
520       // Reset the UnitArray content which were referenced
521       ((AliJetUnitArray*)ref->At(i))->SetUnitTrackID(0);
522       ((AliJetUnitArray*)ref->At(i))->SetUnitEnergy(0.);
523       ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag(kPtSmaller);
524       ((AliJetUnitArray*)ref->At(i))->SetUnitCutFlag2(kPtSmaller);
525       ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlag(kBad);
526       ((AliJetUnitArray*)ref->At(i))->SetUnitSignalFlagC(kTRUE,kBad);
527       ((AliJetUnitArray*)ref->At(i))->SetUnitDetectorFlag(kTpc);
528       ((AliJetUnitArray*)ref->At(i))->SetUnitFlag(kOutJet);
529       ((AliJetUnitArray*)ref->At(i))->ClearUnitTrackRef();
530
531       // Reset process ID
532       AliJetUnitArray* uA = (AliJetUnitArray*)ref->At(i);
533       uA->ResetBit(kIsReferenced);
534       uA->SetUniqueID(0);     
535     }
536
537   // Delete the reference pointer
538   ref->Delete();
539   delete ref;
540
541   Reset();
542
543   return kTRUE;
544 }