/************************************************************************** * 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$ */ //--------------------------------------------------------------------- // Class to calculate the background per unit area // manages the search for jets // Authors: Elena Bruna elena.bruna@yale.edu // Sevil Salur ssalur@lbl.gov // // 2011 : // renamed from AliJetBkg to AliFastJetBkg as this class uses only FASTJET based algos. //--------------------------------------------------------------------- #include #include #include #include #include "AliJetHeader.h" #include "AliFastJetHeaderV1.h" #include "AliAODJet.h" #include "AliFastJetInput.h" #include "AliFastJetBkg.h" #include "fastjet/PseudoJet.hh" #include "fastjet/ClusterSequenceArea.hh" #include "fastjet/AreaDefinition.hh" #include "fastjet/JetDefinition.hh" #include using namespace std; ClassImp(AliFastJetBkg) //////////////////////////////////////////////////////////////////////// AliFastJetBkg::AliFastJetBkg(): TObject(), fHeader(0), fInputFJ(0) { // Default constructor } //______________________________________________________________________ AliFastJetBkg::AliFastJetBkg(const AliFastJetBkg& input): TObject(input), fHeader(input.fHeader), fInputFJ(input.fInputFJ) { // copy constructor } //______________________________________________________________________ AliFastJetBkg& AliFastJetBkg::operator=(const AliFastJetBkg& source) { // Assignment operator. if(this!=&source){ TObject::operator=(source); fHeader = source.fHeader; fInputFJ = source.fInputFJ; } return *this; } //___________________________________________________________________ void AliFastJetBkg::BkgFastJetb(Double_t& rho,Double_t& sigma, Double_t& meanarea) { // Bkg estimation AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option if(debug>0) cout<<"=============== AliFastJetBkg::BkgFastJetb() =========== "< inputParticles=fInputFJ->GetInputParticles(); double rParamBkg = header->GetRparamBkg(); //Radius for background calculation Double_t medianb,sigmab,meanareab; CalcRhob(medianb,sigmab,meanareab,inputParticles,rParamBkg,"All"); rho=medianb; sigma=sigmab; meanarea=meanareab; } //_________________________________________________________________ void AliFastJetBkg::BkgFastJetWoHardest(Double_t& rho,Double_t& sigma, Double_t& meanarea) { // Bkg estimation without hardest jet AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option if(debug) cout<<"=============== AliFastJetBkg::BkgWoHardest() =========== "< inputParticles=fInputFJ->GetInputParticles(); double rParamBkg = header->GetRparamBkg(); //Radius for background calculation Double_t medianb,sigmab,meanareab; CalcRhoWoHardest(medianb,sigmab,meanareab,inputParticles,rParamBkg,"All"); rho=medianb; sigma=sigmab; meanarea=meanareab; } //____________________________________________________________________ void AliFastJetBkg::CalcRhob(Double_t& median,Double_t& sigma,Double_t& meanarea,vector inputParticles,Double_t rParamBkg,TString method) { // calculate rho using the fastjet method AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option fastjet::Strategy strategy = header->GetStrategy(); fastjet::RecombinationScheme recombScheme = header->GetRecombScheme(); fastjet::JetAlgorithm algorithm = header->GetBGAlgorithm(); fastjet::JetDefinition jetDef(algorithm, rParamBkg, 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 ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea); // and from that get an area definition fastjet::AreaType areaType = header->GetAreaType(); areaDef = fastjet::AreaDefinition(areaType,ghost_spec); //fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef); fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef,areaDef); TString comment = "Running FastJet algorithm for BKG calculation with the following setup. "; comment+= "Jet definition: "; comment+= TString(jetDef.description()); // comment+= ". Area definition: "; // comment+= TString(areaDef.description()); comment+= ". Strategy adopted by FastJet: "; comment+= TString(clust_seq.strategy_string()); comment+= Form("Method: %s",method.Data()); header->SetComment(comment); if(debug>0){ cout << "--------------------------------------------------------" << endl; cout << comment << endl; cout << "--------------------------------------------------------" << endl; } vector inclusiveJets = clust_seq.inclusive_jets(0.); vector jets = sorted_by_pt(inclusiveJets); double phiMin = 0, phiMax = 0, rapMin = 0, rapMax = 0; phiMin = 0; phiMax = 2*TMath::Pi(); rapMax = ghostEtamax - rParamBkg; rapMin = - ghostEtamax + rParamBkg; fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax); double medianb, sigmab, meanareab; clust_seq.get_median_rho_and_sigma(inclusiveJets, range, false, medianb, sigmab, meanareab, false); median=medianb; sigma=sigmab; meanarea=meanareab; } //____________________________________________________________________ void AliFastJetBkg::CalcRhoWoHardest(Double_t& median,Double_t& sigma,Double_t& meanarea,vector inputParticles,Double_t rParamBkg,TString method) { // calculate rho (without the hardest jet) using the fastjet method AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option fastjet::Strategy strategy = header->GetStrategy(); fastjet::RecombinationScheme recombScheme = header->GetRecombScheme(); fastjet::JetAlgorithm algorithm = header->GetBGAlgorithm(); fastjet::JetDefinition jetDef(algorithm, rParamBkg, 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 ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea); // and from that get an area definition fastjet::AreaType areaType = header->GetAreaType(); areaDef = fastjet::AreaDefinition(areaType,ghost_spec); //fastjet::ClusterSequenceArea clust_seq(inputParticles, jetDef); fastjet::ClusterSequenceArea clust_seq(inputParticles,jetDef,areaDef); TString comment = "Running FastJet algorithm for BKG calculation with the following setup. "; comment+= "Jet definition: "; comment+= TString(jetDef.description()); // comment+= ". Area definition: "; // comment+= TString(areaDef.description()); comment+= ". Strategy adopted by FastJet: "; comment+= TString(clust_seq.strategy_string()); comment+= Form("Method: %s",method.Data()); header->SetComment(comment); if(debug>0){ cout << "--------------------------------------------------------" << endl; cout << comment << endl; cout << "--------------------------------------------------------" << endl; } vector inclusiveJets = clust_seq.inclusive_jets(0.); vector jets = sorted_by_pt(inclusiveJets); vector jets2=sorted_by_pt(inclusiveJets); if(jets2.size()>=2) jets2.erase(jets2.begin(),jets2.begin()+1); double phiMin = 0, phiMax = 0, rapMin = 0, rapMax = 0; phiMin = 0; phiMax = 2*TMath::Pi(); rapMax = ghostEtamax - rParamBkg; rapMin = - ghostEtamax + rParamBkg; fastjet::RangeDefinition range(rapMin, rapMax, phiMin, phiMax); double medianb, sigmab, meanareab; clust_seq.get_median_rho_and_sigma(jets2, range, false, medianb, sigmab, meanareab, false); median=medianb; sigma=sigmab; meanarea=meanareab; } //___________________________________________________________________ Float_t AliFastJetBkg::BkgFastJet() { // Return background AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option if(debug>0) cout<<"=============== AliFastJetBkg::BkgFastJet() =========== "< inputParticles=fInputFJ->GetInputParticles(); if(debug>0) cout<<"printing inputParticles for BKG "<GetRparamBkg(); //Radius for background calculation Double_t rho=CalcRho(inputParticles,rParamBkg,"All"); if(debug) cout<<"-------- rho (from all part)="<GetDebug(); // debug option if(debug>0) cout<<"=============== AliFastJetBkg::BkgChargedFastJet() =========== "< inputParticlesCharged=fInputFJ->GetInputParticlesCh(); if(debug>0) cout<<"printing CHARGED inputParticles for BKG "<GetRparam(); Double_t rho=CalcRho(inputParticlesCharged,rParam,"Charg"); if(debug>0) cout<<"-------- rho (from CHARGED part)="<GetDebug(); // debug option if(debug>0) cout<<"==============AliFastJetBkg::BkgStat()============="<GetDebug(); // debug option if(debug>0) cout<<"==============AliFastJetBkg::SubtractFastJetBackgCone()============="<GetRparam(); //Hard wired Calorimeter area (get it later from the AliJetReaderHeader.h) Double_t accEMCal=2*0.7*110./180*TMath::Pi();//2.68 area of EMCal Int_t nJ=fAODJets->GetEntries(); //this must be the # of jets... if(debug>0) cout<<"nJets: "<At(0)); jeteta=jettmp->Eta(); jetphi=jettmp->Phi(); acc=EmcalAcceptance(jeteta,jetphi,rCone); if(acc==1)restarea= accEMCal-TMath::Pi()*rc*rc; if(debug) cout<<" acc "<=2) { AliAODJet *jettmp = (AliAODJet*)(fAODJets->At(0)); AliAODJet *jettmp1 = (AliAODJet*)(fAODJets->At(1)); jeteta=jettmp->Eta(); jetphi=jettmp->Phi(); jeteta1=jettmp1->Eta(); jetphi1=jettmp1->Phi(); acc=EmcalAcceptance(jeteta,jetphi,rCone); acc1=EmcalAcceptance(jeteta1,jetphi1,rCone); if(acc1==1 && acc==1)restarea= accEMCal-2*TMath::Pi()*rc*rc; if(acc1==1 && acc==0)restarea= accEMCal-TMath::Pi()*rc*rc; if(acc1==0 && acc==1)restarea= accEMCal-TMath::Pi()*rc*rc; if(debug) cout<<" acc1="<1) cout<<" out of the cone "<0) cout<<"total area left "<0) cout<<"sumpt="< inputParticles,Double_t rParamBkg,TString method) { // calculate rho using the fastjet method AliFastJetHeaderV1 *header = (AliFastJetHeaderV1*)fHeader; Int_t debug = header->GetDebug(); // debug option fastjet::Strategy strategy = header->GetStrategy(); fastjet::RecombinationScheme recombScheme = header->GetRecombScheme(); fastjet::JetAlgorithm algorithm = header->GetBGAlgorithm(); fastjet::JetDefinition jetDef(algorithm, rParamBkg, 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 if (method.Contains("Charg"))ghostEtamax=0.9; fastjet::GhostedAreaSpec ghost_spec(ghostEtamax, activeAreaRepeats, ghostArea); // and from that get an area definition fastjet::AreaType areaType = header->GetAreaType(); areaDef = fastjet::AreaDefinition(areaType,ghost_spec); if(debug>0) cout<<"rParamBkg="<SetComment(comment); if(debug>0){ cout << "--------------------------------------------------------" << endl; cout << comment << endl; cout << "--------------------------------------------------------" << endl; } double ptmin = header->GetPtMin(); vector inclusiveJets = clust_seq.inclusive_jets(ptmin); vector jets = sorted_by_pt(inclusiveJets); if (debug>0) { cout<<"# of BKG jets = "<0) cout<<"bkg in R="<