[7c0fcd5] | 1 | /*
|
---|
| 2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
| 3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
| 4 | *
|
---|
| 5 | * This program is free software: you can redistribute it and/or modify
|
---|
| 6 | * it under the terms of the GNU General Public License as published by
|
---|
| 7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
| 8 | * (at your option) any later version.
|
---|
| 9 | *
|
---|
| 10 | * This program is distributed in the hope that it will be useful,
|
---|
| 11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 13 | * GNU General Public License for more details.
|
---|
| 14 | *
|
---|
| 15 | * You should have received a copy of the GNU General Public License
|
---|
| 16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
| 17 | */
|
---|
| 18 |
|
---|
[d7d2da3] | 19 |
|
---|
| 20 | /** \class JetPileUpSubtractor
|
---|
| 21 | *
|
---|
| 22 | * Subtract pile-up contribution from jets using the fastjet area method
|
---|
| 23 | *
|
---|
| 24 | * $Date: 2012-11-18 15:57:08 +0100 (Sun, 18 Nov 2012) $
|
---|
| 25 | * $Revision: 814 $
|
---|
| 26 | *
|
---|
| 27 | * \author M. Selvaggi - UCL, Louvain-la-Neuve
|
---|
| 28 | *
|
---|
| 29 | */
|
---|
| 30 |
|
---|
| 31 | #include "modules/JetPileUpSubtractor.h"
|
---|
| 32 |
|
---|
| 33 | #include "classes/DelphesClasses.h"
|
---|
| 34 | #include "classes/DelphesFactory.h"
|
---|
| 35 | #include "classes/DelphesFormula.h"
|
---|
| 36 |
|
---|
| 37 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 38 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 39 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 40 |
|
---|
| 41 | #include "TMath.h"
|
---|
| 42 | #include "TString.h"
|
---|
| 43 | #include "TFormula.h"
|
---|
| 44 | #include "TRandom3.h"
|
---|
| 45 | #include "TObjArray.h"
|
---|
| 46 | #include "TDatabasePDG.h"
|
---|
| 47 | #include "TLorentzVector.h"
|
---|
| 48 |
|
---|
| 49 | #include <algorithm>
|
---|
| 50 | #include <stdexcept>
|
---|
| 51 | #include <iostream>
|
---|
| 52 | #include <sstream>
|
---|
| 53 |
|
---|
| 54 | using namespace std;
|
---|
| 55 |
|
---|
| 56 | //------------------------------------------------------------------------------
|
---|
| 57 |
|
---|
| 58 | JetPileUpSubtractor::JetPileUpSubtractor() :
|
---|
[8336b6e] | 59 | fItJetInputArray(0), fItRhoInputArray(0)
|
---|
[d7d2da3] | 60 | {
|
---|
| 61 |
|
---|
| 62 | }
|
---|
| 63 |
|
---|
| 64 | //------------------------------------------------------------------------------
|
---|
| 65 |
|
---|
| 66 | JetPileUpSubtractor::~JetPileUpSubtractor()
|
---|
| 67 | {
|
---|
| 68 |
|
---|
| 69 | }
|
---|
| 70 |
|
---|
| 71 | //------------------------------------------------------------------------------
|
---|
| 72 |
|
---|
| 73 | void JetPileUpSubtractor::Init()
|
---|
| 74 | {
|
---|
| 75 | fJetPTMin = GetDouble("JetPTMin", 20.0);
|
---|
| 76 |
|
---|
| 77 | // import input array(s)
|
---|
| 78 |
|
---|
| 79 | fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
|
---|
| 80 | fItJetInputArray = fJetInputArray->MakeIterator();
|
---|
| 81 |
|
---|
| 82 | fRhoInputArray = ImportArray(GetString("RhoInputArray", "Rho/rho"));
|
---|
[8336b6e] | 83 | fItRhoInputArray = fRhoInputArray->MakeIterator();
|
---|
[d7d2da3] | 84 |
|
---|
| 85 | // create output array(s)
|
---|
| 86 |
|
---|
| 87 | fOutputArray = ExportArray(GetString("OutputArray", "jets"));
|
---|
| 88 |
|
---|
| 89 | }
|
---|
| 90 |
|
---|
| 91 | //------------------------------------------------------------------------------
|
---|
| 92 |
|
---|
| 93 | void JetPileUpSubtractor::Finish()
|
---|
| 94 | {
|
---|
[8336b6e] | 95 | if(fItRhoInputArray) delete fItRhoInputArray;
|
---|
[d7d2da3] | 96 | if(fItJetInputArray) delete fItJetInputArray;
|
---|
| 97 | }
|
---|
| 98 |
|
---|
| 99 | //------------------------------------------------------------------------------
|
---|
| 100 |
|
---|
| 101 | void JetPileUpSubtractor::Process()
|
---|
| 102 | {
|
---|
[e9023b9] | 103 | Candidate *candidate, *object;
|
---|
[d7d2da3] | 104 | TLorentzVector momentum, area;
|
---|
[8336b6e] | 105 | Double_t eta = 0.0;
|
---|
[d7d2da3] | 106 | Double_t rho = 0.0;
|
---|
| 107 |
|
---|
| 108 | // loop over all input candidates
|
---|
| 109 | fItJetInputArray->Reset();
|
---|
| 110 | while((candidate = static_cast<Candidate*>(fItJetInputArray->Next())))
|
---|
| 111 | {
|
---|
| 112 | momentum = candidate->Momentum;
|
---|
| 113 | area = candidate->Area;
|
---|
[8336b6e] | 114 | eta = TMath::Abs(momentum.Eta());
|
---|
| 115 |
|
---|
| 116 | // find rho
|
---|
| 117 | rho = 0.0;
|
---|
[e9023b9] | 118 | if(fRhoInputArray)
|
---|
[8336b6e] | 119 | {
|
---|
[e9023b9] | 120 | fItRhoInputArray->Reset();
|
---|
| 121 | while((object = static_cast<Candidate*>(fItRhoInputArray->Next())))
|
---|
[8336b6e] | 122 | {
|
---|
[e9023b9] | 123 | if(eta >= object->Edges[0] && eta < object->Edges[1])
|
---|
| 124 | {
|
---|
| 125 | rho = object->Momentum.Pt();
|
---|
| 126 | }
|
---|
[8336b6e] | 127 | }
|
---|
[e9023b9] | 128 | }
|
---|
[d7d2da3] | 129 |
|
---|
| 130 | // apply pile-up correction
|
---|
| 131 | if(momentum.Pt() <= rho * area.Pt()) continue;
|
---|
| 132 |
|
---|
| 133 | momentum -= rho * area;
|
---|
[e9023b9] | 134 |
|
---|
[d7d2da3] | 135 | if(momentum.Pt() <= fJetPTMin) continue;
|
---|
| 136 |
|
---|
| 137 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
| 138 | candidate->Momentum = momentum;
|
---|
| 139 |
|
---|
| 140 | fOutputArray->Add(candidate);
|
---|
| 141 | }
|
---|
| 142 | }
|
---|
| 143 |
|
---|
| 144 | //------------------------------------------------------------------------------
|
---|