#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())
{
}
AliFastJetFinder::~AliFastJetFinder()
{
delete fInputFJ;
delete fJetBkg;
}
void AliFastJetFinder::FindJets()
{
AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader;
Int_t debug = header->GetDebug();
Bool_t bgMode = header->GetBGMode();
if(debug>0) cout<<"----------in AliFastJetFinder::FindJets() ------------------"<<endl;
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;
}
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);
fastjet::AreaDefinition areaDef;
double ghostEtamax = header->GetGhostEtaMax();
double ghostArea = header->GetGhostArea();
int activeAreaRepeats = header->GetActiveAreaRepeats();
fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
fastjet::AreaType areaType = header->GetAreaType();
areaDef = fastjet::AreaDefinition(areaType,ghostSpec);
fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
vector<fastjet::PseudoJet> jets;
if(bgMode)
{
fastjet::JetAlgorithm algorithmBkg = header->GetBGAlgorithm();
fastjet::JetDefinition jetDefBkg(algorithmBkg, rBkgParam, recombScheme, strategy);
fastjet::ClusterSequenceArea clust_seq_bkg(inputParticles, jetDefBkg, areaDef);
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;
}
double ptmin = header->GetPtMin();
vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets();
double rapMax = header->GetRapMax();
double rapMin = header->GetRapMin();
double phiMax = header->GetPhiMax();
double phiMin = header->GetPhiMin();
fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax);
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);
vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(rho,ptmin);
jets = sorted_by_pt(subJets);
}
else {
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;
}
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++) {
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;
AliAODJet aodjet (jets[j].px(), jets[j].py(), jets[j].pz(), jets[j].E());
aodjet.SetEffArea(area,areaError);
if(debug>0) aodjet.Print("");
for (int i=0; i < nCon; i++)
{
fastjet::PseudoJet mPart=constituents[i];
ind[i]=mPart.user_index();
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);
}
}
}
AddJet(aodjet);
}
}
void AliFastJetFinder::RunTest(const char* datafile)
{
vector<fastjet::PseudoJet> inputParticles;
Float_t px,py,pz,en;
ifstream in;
Int_t nlines = 0;
in.open(datafile);
while (1) {
in >> px >> py >> pz >> en;
if (!in.good()) break;
nlines++;
inputParticles.push_back(fastjet::PseudoJet(px,py,pz,en));
}
in.close();
double rParam = 1.0;
fastjet::Strategy strategy = fastjet::Best;
fastjet::RecombinationScheme recombScheme = fastjet::BIpt_scheme;
fastjet::JetDefinition jetDef(fastjet::kt_algorithm, rParam, recombScheme, strategy);
fastjet::AreaDefinition areaDef;
double ghostEtamax = 7.0;
double ghostArea = 0.05;
int activeAreaRepeats = 1;
fastjet::GhostedAreaSpec ghostSpec(ghostEtamax, activeAreaRepeats, ghostArea);
areaDef = fastjet::AreaDefinition(fastjet::active_area,ghostSpec);
fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef, areaDef);
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;
double ptmin = 5.0;
vector<fastjet::PseudoJet> inclusiveJets = clust_seq.inclusive_jets(ptmin);
cout << "Number of unclustered particles: " << clust_seq.unclustered_particles().size() << endl;
double rapMax = ghostEtamax - rParam;
fastjet::RangeDefinition range(rapMax);
vector<fastjet::PseudoJet> subJets = clust_seq.subtracted_jets(range,ptmin);
cout << "Printing inclusive jets after background subtraction \n";
cout << "------------------------------------------------------\n";
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
{
fHeader->Write();
}
Bool_t AliFastJetFinder::ProcessEvent()
{
fInputFJ->SetHeader(fHeader);
fInputFJ->SetCalTrkEvent(GetCalTrkEvent());
fInputFJ->FillInput();
FindJets();
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;
}