ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/
 
/* $Id$ */

//---------------------------------------------------------------------
// FastJet v2.3.4 finder algorithm interface
// Last modification: Neutral cell energy included in the jet reconstruction
//
// Authors: Rafael.Diaz.Valdes@cern.ch
//          Magali.estienne@subatech.in2p3.fr (neutral part + bg subtraction option)
//
// ** 2011 magali.estienne@subatech.in2p3.fr &  alexandre.shabetai@cern.ch
// new implementation of background subtraction
// allowing to subtract bkg using a different algo than the one used for signal jets
//---------------------------------------------------------------------

#include <Riostream.h>

#include "AliFastJetFinder.h"
#include "AliFastJetHeaderV1.h"
#include "AliFastJetInput.h"
#include "AliFastJetBkg.h"
#include "AliAODJetEventBackground.h"
#include "AliAODJet.h"

#include "fastjet/PseudoJet.hh"
#include "fastjet/ClusterSequenceArea.hh"
#include "fastjet/AreaDefinition.hh"
#include "fastjet/JetDefinition.hh"

#include<vector> 

using namespace std;

ClassImp(AliFastJetFinder)

////////////////////////////////////////////////////////////////////////

AliFastJetFinder::AliFastJetFinder():
  AliJetFinder(),
  fInputFJ(new AliFastJetInput()),
  fJetBkg(new  AliFastJetBkg())
{
  // Constructor
}

//____________________________________________________________________________
AliFastJetFinder::~AliFastJetFinder()
{
  // destructor
  delete  fInputFJ;
  delete  fJetBkg;

}

//______________________________________________________________________________
void AliFastJetFinder::FindJets()
{
  // runs a FASTJET based algo

  //pick up fastjet header
  AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
  Int_t  debug  = header->GetDebug();     // debug option
  Bool_t bgMode = header->GetBGMode();    // choose to subtract BG or not
  if(debug>0) cout<<"----------in AliFastJetFinder::FindJets() ------------------"<<endl;

  // RUN ALGORITHM  
  // read input particles -----------------------------

  vector<fastjet::PseudoJet> inputParticles=fInputFJ->GetInputParticles();
  if(inputParticles.size()==0){
    if(debug>0) Printf("%s:%d No input particles found, skipping event",(char*)__FILE__,__LINE__);
    return;
  }

  // create an object that represents your choice of jet algorithm, and 
  // the associated parameters
  double rParam = header->GetRparam();
  double rBkgParam = header->GetRparamBkg();
  fastjet::Strategy strategy = header->GetStrategy();
  fastjet::RecombinationScheme recombScheme = header->GetRecombScheme();
  fastjet::JetAlgorithm algorithm = header->GetAlgorithm(); 
  fastjet::JetDefinition jetDef(algorithm, rParam, recombScheme, strategy);

  // create an object that specifies how we to define the area
  fastjet::AreaDefinition areaDef;
  double ghostEtamax       = header->GetGhostEtaMax(); 
  double ghostArea         = header->GetGhostArea(); 
  int    activeAreaRepeats = header->GetActiveAreaRepeats(); 
  
  // now create the object that holds info about ghosts
  fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
  // and from that get an area definition
  fastjet::AreaType areaType = header->GetAreaType();
  areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
  
  //***************************** JETS FINDING
  // run the jet clustering with the above jet definition
  fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);

  vector<fastjet::PseudoJet> jets;

  if(bgMode) // Do BG subtraction directly with the same algorithm (cambridge or kt) for jet signal and background
    {
      //***************************** JETS FINDING FOR RHO ESTIMATION
      // run the jet clustering with the above jet definition
      fastjet::JetAlgorithm algorithmBkg = header->GetBGAlgorithm();
      fastjet::JetDefinition jetDefBkg(algorithmBkg, rBkgParam, recombScheme, strategy);
      fastjet::ClusterSequenceArea clust_seq_bkg(inputParticles, jetDefBkg, areaDef);

      // save a comment in the header
      TString comment = "Running FastJet algorithm with the following setup. ";
      comment+= "Jet definition: ";
      comment+= TString(jetDef.description());
      comment+= "Jet bckg definition: ";
      comment+= TString(jetDefBkg.description());
      comment+= ". Area definition: ";
      comment+= TString(areaDef.description());
      comment+= ". Strategy adopted by FastJet: ";
      comment+= TString(clust_seq.strategy_string());
      header->SetComment(comment);
      if(debug>0){
	cout << "--------------------------------------------------------" << endl;
	cout << comment << endl;
	cout << "--------------------------------------------------------" << endl;
      }
      
      // extract the inclusive jets sorted by pt
      double ptmin = header->GetPtMin(); 
      vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets();
      
      //subtract background // ===========================================
      // set the rapididty , phi range within which to study the background 
      double rapMax = header->GetRapMax(); 
      double rapMin = header->GetRapMin();
      double phiMax = header->GetPhiMax();
      double phiMin = header->GetPhiMin();
      fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
     
      // Extract rho and sigma
      Double_t rho = 0.;
      Double_t sigma = 0.;
      Double_t meanarea = 0.;
      Bool_t kUse4VectorArea = header->Use4VectorArea();
      vector<fastjet::PseudoJet> bkgJets = clust_seq_bkg.inclusive_jets();
      clust_seq_bkg.get_median_rho_and_sigma(bkgJets,range, kUse4VectorArea, rho, sigma, meanarea, false);

      // subtract background and extract jets bkg subtracted
      vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(rho,ptmin);
      
      // sort jets into increasing pt
      jets = sorted_by_pt(subJets);  

    }

  else { // No BG subtraction!!!!!!!! Default header is bgmode=0.
    
    // save a comment in the header
    TString comment = "Running FastJet algorithm with the following setup. ";
    comment+= "Jet definition: ";
    comment+= TString(jetDef.description());
    comment+= ". Strategy adopted by FastJet: ";
    comment+= TString(clust_seq.strategy_string());
    header->SetComment(comment);
    if(debug>0){
      cout << "--------------------------------------------------------" << endl;
      cout << comment << endl;
      cout << "--------------------------------------------------------" << endl;
    }
  
    // extract the inclusive jets with pt > ptmin, sorted by pt
    double ptmin = header->GetPtMin(); 
    vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
    
    jets = sorted_by_pt(inclusiveJets); 
  
  }
   
  for (size_t j = 0; j < jets.size(); j++) { // loop for jets     

    double area      = clust_seq.area(jets[j]);
    double areaError = clust_seq.area_error(jets[j]);

    if(debug>0) printf("Jet found %5d %9.5f %8.5f %10.3f \n",(Int_t)j,jets[j].rap(),jets[j].phi(),jets[j].perp());

    vector<fastjet::PseudoJet> constituents = clust_seq.constituents(jets[j]);
    int nCon= constituents.size();
    TArrayI ind(nCon);
      
    if ((jets[j].eta() > (header->GetJetEtaMax())) ||
	(jets[j].eta() < (header->GetJetEtaMin())) ||
	(jets[j].phi() > (header->GetJetPhiMax())) ||
	(jets[j].phi() < (header->GetJetPhiMin())) ||
	(jets[j].perp() < header->GetPtMin())) continue; // acceptance eta range and etmin

    // go to write AOD  info
    AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
    aodjet.SetEffArea(area,areaError);
    //cout << "Printing jet " << endl;
    if(debug>0) aodjet.Print("");
      
    for (int i=0; i < nCon; i++)
      {
	fastjet::PseudoJet mPart=constituents[i];
	ind[i]=mPart.user_index();

	// Jet constituents (charged tracks) added to the AliAODJet
	AliJetCalTrkEvent* calEvt  = GetCalTrkEvent();
	for(Int_t itrack=0; itrack<calEvt->GetNCalTrkTracks(); itrack++)
	  {
	    if(itrack==ind[i])
	      {
		TObject *track = calEvt->GetCalTrkTrack(itrack)->GetTrackObject();
		aodjet.AddTrack(track);
	      }
	  }
      } // End loop on Constituents

    AddJet(aodjet);
           
  } // end loop for jets

}

//____________________________________________________________________________
void AliFastJetFinder::RunTest(const char* datafile)
{
  // This simple test run the kt algorithm for an ascii file testdata.dat

  // read input particles -----------------------------
  vector<fastjet::PseudoJet> inputParticles;
  Float_t px,py,pz,en;
  ifstream in;
  Int_t nlines = 0;
  // we assume a file basic.dat in the current directory
  // this file has 3 columns of float data
  in.open(datafile);
  while (1) {
    in >> px >> py >> pz >> en;
    if (!in.good()) break;
    //printf("px=%8f, py=%8f, pz=%8fn",px,py,pz);
    nlines++;
    inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en)); 
  }
  //printf(" found %d pointsn",nlines);
  in.close();
  //////////////////////////////////////////////////
 
  // create an object that represents your choice of jet algorithm, and 
  // the associated parameters
  double rParam = 1.0;
  fastjet::Strategy strategy = fastjet::Best;
  fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
  fastjet::JetDefinition jetDef(fastjet::kt_algorithm, rParam, recombScheme, strategy);
  
  // create an object that specifies how we to define the area
  fastjet::AreaDefinition areaDef;
  double ghostEtamax = 7.0;
  double ghostArea    = 0.05;
  int    activeAreaRepeats = 1;
  
  // now create the object that holds info about ghosts
  fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
  // and from that get an area definition
  areaDef = fastjet::AreaDefinition(fastjet::active_area,ghostSpec);

  // run the jet clustering with the above jet definition
  fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
  
  // tell the user what was done
  cout << "--------------------------------------------------------" << endl;
  cout << "Jet definition was: " << jetDef.description() << endl;
  cout << "Area definition was: " << areaDef.description() << endl;
  cout << "Strategy adopted by FastJet was "<< clust_seq.strategy_string()<<endl<<endl;
  cout << "--------------------------------------------------------" << endl;
  
  // extract the inclusive jets with pt > 5 GeV, sorted by pt
  double ptmin = 5.0;
  vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
  
  cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
 
  //subtract background // ===========================================
  // set the rapididty range within which to study the background 
  double rapMax = ghostEtamax - rParam;
  fastjet::RangeDefinition range(rapMax);
  // subtract background
  vector<fastjet::PseudoJet> subJets =  clust_seq.subtracted_jets(range,ptmin);  
  
  // print them out //================================================
  cout << "Printing inclusive jets  after background subtraction \n";
  cout << "------------------------------------------------------\n";
  // sort jets into increasing pt
  vector<fastjet::PseudoJet> jets = sorted_by_pt(subJets);  

  printf(" ijet   rap      phi        Pt         area  +-   err\n");
  for (size_t j = 0; j < jets.size(); j++) {
    double area      = clust_seq.area(jets[j]);
    double areaError = clust_seq.area_error(jets[j]);
    printf("%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);
  }
  cout << endl;
  // ================================================================

}

//____________________________________________________________________________
void AliFastJetFinder::WriteJHeaderToFile() const
{
  // Write Jet Header to file
  fHeader->Write();

}

//____________________________________________________________________________
Bool_t AliFastJetFinder::ProcessEvent()
{
  // Process one event
  // Charged only or charged+neutral jets

  fInputFJ->SetHeader(fHeader);
  fInputFJ->SetCalTrkEvent(GetCalTrkEvent());
  fInputFJ->FillInput();
  
  // Find Jets
  FindJets(); 

  // Background
  if(fAODEvBkg){
    AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;

    fJetBkg->SetHeader(fHeader);
    fJetBkg->SetFastJetInput(fInputFJ);
    
     Int_t count = 0;  
     if(header->GetBkgFastJetb()){
        Double_t sigma1 = 0 , meanarea1= 0, bkg1 = 0;
        fJetBkg->BkgFastJetb(bkg1,sigma1,meanarea1);
        fAODEvBkg->SetBackground(count,bkg1,sigma1,meanarea1);
        count++;
     }

     if(header->GetBkgFastJetWoHardest()){
        Double_t sigma2 = 0 , meanarea2 = 0, bkg2 = 0;
        fJetBkg->BkgFastJetWoHardest(bkg2,sigma2,meanarea2);
        fAODEvBkg->SetBackground(count,bkg2,sigma2,meanarea2);
        count++;
     }

  }  

  Reset();  
  return kTRUE;

}
 AliFastJetFinder.cxx:1
 AliFastJetFinder.cxx:2
 AliFastJetFinder.cxx:3
 AliFastJetFinder.cxx:4
 AliFastJetFinder.cxx:5
 AliFastJetFinder.cxx:6
 AliFastJetFinder.cxx:7
 AliFastJetFinder.cxx:8
 AliFastJetFinder.cxx:9
 AliFastJetFinder.cxx:10
 AliFastJetFinder.cxx:11
 AliFastJetFinder.cxx:12
 AliFastJetFinder.cxx:13
 AliFastJetFinder.cxx:14
 AliFastJetFinder.cxx:15
 AliFastJetFinder.cxx:16
 AliFastJetFinder.cxx:17
 AliFastJetFinder.cxx:18
 AliFastJetFinder.cxx:19
 AliFastJetFinder.cxx:20
 AliFastJetFinder.cxx:21
 AliFastJetFinder.cxx:22
 AliFastJetFinder.cxx:23
 AliFastJetFinder.cxx:24
 AliFastJetFinder.cxx:25
 AliFastJetFinder.cxx:26
 AliFastJetFinder.cxx:27
 AliFastJetFinder.cxx:28
 AliFastJetFinder.cxx:29
 AliFastJetFinder.cxx:30
 AliFastJetFinder.cxx:31
 AliFastJetFinder.cxx:32
 AliFastJetFinder.cxx:33
 AliFastJetFinder.cxx:34
 AliFastJetFinder.cxx:35
 AliFastJetFinder.cxx:36
 AliFastJetFinder.cxx:37
 AliFastJetFinder.cxx:38
 AliFastJetFinder.cxx:39
 AliFastJetFinder.cxx:40
 AliFastJetFinder.cxx:41
 AliFastJetFinder.cxx:42
 AliFastJetFinder.cxx:43
 AliFastJetFinder.cxx:44
 AliFastJetFinder.cxx:45
 AliFastJetFinder.cxx:46
 AliFastJetFinder.cxx:47
 AliFastJetFinder.cxx:48
 AliFastJetFinder.cxx:49
 AliFastJetFinder.cxx:50
 AliFastJetFinder.cxx:51
 AliFastJetFinder.cxx:52
 AliFastJetFinder.cxx:53
 AliFastJetFinder.cxx:54
 AliFastJetFinder.cxx:55
 AliFastJetFinder.cxx:56
 AliFastJetFinder.cxx:57
 AliFastJetFinder.cxx:58
 AliFastJetFinder.cxx:59
 AliFastJetFinder.cxx:60
 AliFastJetFinder.cxx:61
 AliFastJetFinder.cxx:62
 AliFastJetFinder.cxx:63
 AliFastJetFinder.cxx:64
 AliFastJetFinder.cxx:65
 AliFastJetFinder.cxx:66
 AliFastJetFinder.cxx:67
 AliFastJetFinder.cxx:68
 AliFastJetFinder.cxx:69
 AliFastJetFinder.cxx:70
 AliFastJetFinder.cxx:71
 AliFastJetFinder.cxx:72
 AliFastJetFinder.cxx:73
 AliFastJetFinder.cxx:74
 AliFastJetFinder.cxx:75
 AliFastJetFinder.cxx:76
 AliFastJetFinder.cxx:77
 AliFastJetFinder.cxx:78
 AliFastJetFinder.cxx:79
 AliFastJetFinder.cxx:80
 AliFastJetFinder.cxx:81
 AliFastJetFinder.cxx:82
 AliFastJetFinder.cxx:83
 AliFastJetFinder.cxx:84
 AliFastJetFinder.cxx:85
 AliFastJetFinder.cxx:86
 AliFastJetFinder.cxx:87
 AliFastJetFinder.cxx:88
 AliFastJetFinder.cxx:89
 AliFastJetFinder.cxx:90
 AliFastJetFinder.cxx:91
 AliFastJetFinder.cxx:92
 AliFastJetFinder.cxx:93
 AliFastJetFinder.cxx:94
 AliFastJetFinder.cxx:95
 AliFastJetFinder.cxx:96
 AliFastJetFinder.cxx:97
 AliFastJetFinder.cxx:98
 AliFastJetFinder.cxx:99
 AliFastJetFinder.cxx:100
 AliFastJetFinder.cxx:101
 AliFastJetFinder.cxx:102
 AliFastJetFinder.cxx:103
 AliFastJetFinder.cxx:104
 AliFastJetFinder.cxx:105
 AliFastJetFinder.cxx:106
 AliFastJetFinder.cxx:107
 AliFastJetFinder.cxx:108
 AliFastJetFinder.cxx:109
 AliFastJetFinder.cxx:110
 AliFastJetFinder.cxx:111
 AliFastJetFinder.cxx:112
 AliFastJetFinder.cxx:113
 AliFastJetFinder.cxx:114
 AliFastJetFinder.cxx:115
 AliFastJetFinder.cxx:116
 AliFastJetFinder.cxx:117
 AliFastJetFinder.cxx:118
 AliFastJetFinder.cxx:119
 AliFastJetFinder.cxx:120
 AliFastJetFinder.cxx:121
 AliFastJetFinder.cxx:122
 AliFastJetFinder.cxx:123
 AliFastJetFinder.cxx:124
 AliFastJetFinder.cxx:125
 AliFastJetFinder.cxx:126
 AliFastJetFinder.cxx:127
 AliFastJetFinder.cxx:128
 AliFastJetFinder.cxx:129
 AliFastJetFinder.cxx:130
 AliFastJetFinder.cxx:131
 AliFastJetFinder.cxx:132
 AliFastJetFinder.cxx:133
 AliFastJetFinder.cxx:134
 AliFastJetFinder.cxx:135
 AliFastJetFinder.cxx:136
 AliFastJetFinder.cxx:137
 AliFastJetFinder.cxx:138
 AliFastJetFinder.cxx:139
 AliFastJetFinder.cxx:140
 AliFastJetFinder.cxx:141
 AliFastJetFinder.cxx:142
 AliFastJetFinder.cxx:143
 AliFastJetFinder.cxx:144
 AliFastJetFinder.cxx:145
 AliFastJetFinder.cxx:146
 AliFastJetFinder.cxx:147
 AliFastJetFinder.cxx:148
 AliFastJetFinder.cxx:149
 AliFastJetFinder.cxx:150
 AliFastJetFinder.cxx:151
 AliFastJetFinder.cxx:152
 AliFastJetFinder.cxx:153
 AliFastJetFinder.cxx:154
 AliFastJetFinder.cxx:155
 AliFastJetFinder.cxx:156
 AliFastJetFinder.cxx:157
 AliFastJetFinder.cxx:158
 AliFastJetFinder.cxx:159
 AliFastJetFinder.cxx:160
 AliFastJetFinder.cxx:161
 AliFastJetFinder.cxx:162
 AliFastJetFinder.cxx:163
 AliFastJetFinder.cxx:164
 AliFastJetFinder.cxx:165
 AliFastJetFinder.cxx:166
 AliFastJetFinder.cxx:167
 AliFastJetFinder.cxx:168
 AliFastJetFinder.cxx:169
 AliFastJetFinder.cxx:170
 AliFastJetFinder.cxx:171
 AliFastJetFinder.cxx:172
 AliFastJetFinder.cxx:173
 AliFastJetFinder.cxx:174
 AliFastJetFinder.cxx:175
 AliFastJetFinder.cxx:176
 AliFastJetFinder.cxx:177
 AliFastJetFinder.cxx:178
 AliFastJetFinder.cxx:179
 AliFastJetFinder.cxx:180
 AliFastJetFinder.cxx:181
 AliFastJetFinder.cxx:182
 AliFastJetFinder.cxx:183
 AliFastJetFinder.cxx:184
 AliFastJetFinder.cxx:185
 AliFastJetFinder.cxx:186
 AliFastJetFinder.cxx:187
 AliFastJetFinder.cxx:188
 AliFastJetFinder.cxx:189
 AliFastJetFinder.cxx:190
 AliFastJetFinder.cxx:191
 AliFastJetFinder.cxx:192
 AliFastJetFinder.cxx:193
 AliFastJetFinder.cxx:194
 AliFastJetFinder.cxx:195
 AliFastJetFinder.cxx:196
 AliFastJetFinder.cxx:197
 AliFastJetFinder.cxx:198
 AliFastJetFinder.cxx:199
 AliFastJetFinder.cxx:200
 AliFastJetFinder.cxx:201
 AliFastJetFinder.cxx:202
 AliFastJetFinder.cxx:203
 AliFastJetFinder.cxx:204
 AliFastJetFinder.cxx:205
 AliFastJetFinder.cxx:206
 AliFastJetFinder.cxx:207
 AliFastJetFinder.cxx:208
 AliFastJetFinder.cxx:209
 AliFastJetFinder.cxx:210
 AliFastJetFinder.cxx:211
 AliFastJetFinder.cxx:212
 AliFastJetFinder.cxx:213
 AliFastJetFinder.cxx:214
 AliFastJetFinder.cxx:215
 AliFastJetFinder.cxx:216
 AliFastJetFinder.cxx:217
 AliFastJetFinder.cxx:218
 AliFastJetFinder.cxx:219
 AliFastJetFinder.cxx:220
 AliFastJetFinder.cxx:221
 AliFastJetFinder.cxx:222
 AliFastJetFinder.cxx:223
 AliFastJetFinder.cxx:224
 AliFastJetFinder.cxx:225
 AliFastJetFinder.cxx:226
 AliFastJetFinder.cxx:227
 AliFastJetFinder.cxx:228
 AliFastJetFinder.cxx:229
 AliFastJetFinder.cxx:230
 AliFastJetFinder.cxx:231
 AliFastJetFinder.cxx:232
 AliFastJetFinder.cxx:233
 AliFastJetFinder.cxx:234
 AliFastJetFinder.cxx:235
 AliFastJetFinder.cxx:236
 AliFastJetFinder.cxx:237
 AliFastJetFinder.cxx:238
 AliFastJetFinder.cxx:239
 AliFastJetFinder.cxx:240
 AliFastJetFinder.cxx:241
 AliFastJetFinder.cxx:242
 AliFastJetFinder.cxx:243
 AliFastJetFinder.cxx:244
 AliFastJetFinder.cxx:245
 AliFastJetFinder.cxx:246
 AliFastJetFinder.cxx:247
 AliFastJetFinder.cxx:248
 AliFastJetFinder.cxx:249
 AliFastJetFinder.cxx:250
 AliFastJetFinder.cxx:251
 AliFastJetFinder.cxx:252
 AliFastJetFinder.cxx:253
 AliFastJetFinder.cxx:254
 AliFastJetFinder.cxx:255
 AliFastJetFinder.cxx:256
 AliFastJetFinder.cxx:257
 AliFastJetFinder.cxx:258
 AliFastJetFinder.cxx:259
 AliFastJetFinder.cxx:260
 AliFastJetFinder.cxx:261
 AliFastJetFinder.cxx:262
 AliFastJetFinder.cxx:263
 AliFastJetFinder.cxx:264
 AliFastJetFinder.cxx:265
 AliFastJetFinder.cxx:266
 AliFastJetFinder.cxx:267
 AliFastJetFinder.cxx:268
 AliFastJetFinder.cxx:269
 AliFastJetFinder.cxx:270
 AliFastJetFinder.cxx:271
 AliFastJetFinder.cxx:272
 AliFastJetFinder.cxx:273
 AliFastJetFinder.cxx:274
 AliFastJetFinder.cxx:275
 AliFastJetFinder.cxx:276
 AliFastJetFinder.cxx:277
 AliFastJetFinder.cxx:278
 AliFastJetFinder.cxx:279
 AliFastJetFinder.cxx:280
 AliFastJetFinder.cxx:281
 AliFastJetFinder.cxx:282
 AliFastJetFinder.cxx:283
 AliFastJetFinder.cxx:284
 AliFastJetFinder.cxx:285
 AliFastJetFinder.cxx:286
 AliFastJetFinder.cxx:287
 AliFastJetFinder.cxx:288
 AliFastJetFinder.cxx:289
 AliFastJetFinder.cxx:290
 AliFastJetFinder.cxx:291
 AliFastJetFinder.cxx:292
 AliFastJetFinder.cxx:293
 AliFastJetFinder.cxx:294
 AliFastJetFinder.cxx:295
 AliFastJetFinder.cxx:296
 AliFastJetFinder.cxx:297
 AliFastJetFinder.cxx:298
 AliFastJetFinder.cxx:299
 AliFastJetFinder.cxx:300
 AliFastJetFinder.cxx:301
 AliFastJetFinder.cxx:302
 AliFastJetFinder.cxx:303
 AliFastJetFinder.cxx:304
 AliFastJetFinder.cxx:305
 AliFastJetFinder.cxx:306
 AliFastJetFinder.cxx:307
 AliFastJetFinder.cxx:308
 AliFastJetFinder.cxx:309
 AliFastJetFinder.cxx:310
 AliFastJetFinder.cxx:311
 AliFastJetFinder.cxx:312
 AliFastJetFinder.cxx:313
 AliFastJetFinder.cxx:314
 AliFastJetFinder.cxx:315
 AliFastJetFinder.cxx:316
 AliFastJetFinder.cxx:317
 AliFastJetFinder.cxx:318
 AliFastJetFinder.cxx:319
 AliFastJetFinder.cxx:320
 AliFastJetFinder.cxx:321
 AliFastJetFinder.cxx:322
 AliFastJetFinder.cxx:323
 AliFastJetFinder.cxx:324
 AliFastJetFinder.cxx:325
 AliFastJetFinder.cxx:326
 AliFastJetFinder.cxx:327
 AliFastJetFinder.cxx:328
 AliFastJetFinder.cxx:329
 AliFastJetFinder.cxx:330
 AliFastJetFinder.cxx:331
 AliFastJetFinder.cxx:332
 AliFastJetFinder.cxx:333
 AliFastJetFinder.cxx:334
 AliFastJetFinder.cxx:335
 AliFastJetFinder.cxx:336
 AliFastJetFinder.cxx:337
 AliFastJetFinder.cxx:338
 AliFastJetFinder.cxx:339
 AliFastJetFinder.cxx:340
 AliFastJetFinder.cxx:341
 AliFastJetFinder.cxx:342
 AliFastJetFinder.cxx:343
 AliFastJetFinder.cxx:344
 AliFastJetFinder.cxx:345
 AliFastJetFinder.cxx:346
 AliFastJetFinder.cxx:347
 AliFastJetFinder.cxx:348
 AliFastJetFinder.cxx:349
 AliFastJetFinder.cxx:350
 AliFastJetFinder.cxx:351
 AliFastJetFinder.cxx:352
 AliFastJetFinder.cxx:353
 AliFastJetFinder.cxx:354
 AliFastJetFinder.cxx:355
 AliFastJetFinder.cxx:356
 AliFastJetFinder.cxx:357
 AliFastJetFinder.cxx:358
 AliFastJetFinder.cxx:359
 AliFastJetFinder.cxx:360
 AliFastJetFinder.cxx:361
 AliFastJetFinder.cxx:362
 AliFastJetFinder.cxx:363
 AliFastJetFinder.cxx:364
 AliFastJetFinder.cxx:365
 AliFastJetFinder.cxx:366
 AliFastJetFinder.cxx:367
 AliFastJetFinder.cxx:368
 AliFastJetFinder.cxx:369