/*
* Delphes: a framework for fast simulation of a generic collider experiment
* Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see .
*/
/** \class JetPileUpSubtractor
*
* Subtract pile-up contribution from jets using the fastjet area method
*
* \author M. Selvaggi - UCL, Louvain-la-Neuve
*
*/
#include "modules/JetPileUpSubtractor.h"
#include "classes/DelphesClasses.h"
#include "classes/DelphesFactory.h"
#include "classes/DelphesFormula.h"
#include "ExRootAnalysis/ExRootClassifier.h"
#include "ExRootAnalysis/ExRootFilter.h"
#include "ExRootAnalysis/ExRootResult.h"
#include "TDatabasePDG.h"
#include "TFormula.h"
#include "TLorentzVector.h"
#include "TMath.h"
#include "TObjArray.h"
#include "TRandom3.h"
#include "TString.h"
#include
#include
#include
#include
using namespace std;
//------------------------------------------------------------------------------
JetPileUpSubtractor::JetPileUpSubtractor() :
fItJetInputArray(0), fItRhoInputArray(0)
{
}
//------------------------------------------------------------------------------
JetPileUpSubtractor::~JetPileUpSubtractor()
{
}
//------------------------------------------------------------------------------
void JetPileUpSubtractor::Init()
{
fJetPTMin = GetDouble("JetPTMin", 20.0);
// import input array(s)
fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
fItJetInputArray = fJetInputArray->MakeIterator();
fRhoInputArray = ImportArray(GetString("RhoInputArray", "Rho/rho"));
fItRhoInputArray = fRhoInputArray->MakeIterator();
// create output array(s)
fOutputArray = ExportArray(GetString("OutputArray", "jets"));
}
//------------------------------------------------------------------------------
void JetPileUpSubtractor::Finish()
{
if(fItRhoInputArray) delete fItRhoInputArray;
if(fItJetInputArray) delete fItJetInputArray;
}
//------------------------------------------------------------------------------
void JetPileUpSubtractor::Process()
{
Candidate *candidate, *object;
TLorentzVector momentum, area;
Double_t eta = 0.0;
Double_t rho = 0.0;
// loop over all input candidates
fItJetInputArray->Reset();
while((candidate = static_cast(fItJetInputArray->Next())))
{
momentum = candidate->Momentum;
area = candidate->Area;
eta = momentum.Eta();
// find rho
rho = 0.0;
if(fRhoInputArray)
{
fItRhoInputArray->Reset();
while((object = static_cast(fItRhoInputArray->Next())))
{
if(eta >= object->Edges[0] && eta < object->Edges[1])
{
rho = object->Momentum.Pt();
}
}
}
// apply pile-up correction
if(momentum.Pt() <= rho * area.Pt()) continue;
momentum -= rho * area;
if(momentum.Pt() <= fJetPTMin) continue;
candidate = static_cast(candidate->Clone());
candidate->Momentum = momentum;
fOutputArray->Add(candidate);
}
}
//------------------------------------------------------------------------------