Fork me on GitHub

source: git/modules/JetPileUpSubtractor.cc@ 239e1d0

ImprovedOutputFile Timing dual_readout llp 3.3.3pre11
Last change on this file since 239e1d0 was d4b9697, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix rapidity ranges for GridMedianBackgroundEstimator

  • Property mode set to 100644
File size: 3.7 KB
Line 
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
19
20/** \class JetPileUpSubtractor
21 *
22 * Subtract pile-up contribution from jets using the fastjet area method
23 *
24 * \author M. Selvaggi - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/JetPileUpSubtractor.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootResult.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootClassifier.h"
37
38#include "TMath.h"
39#include "TString.h"
40#include "TFormula.h"
41#include "TRandom3.h"
42#include "TObjArray.h"
43#include "TDatabasePDG.h"
44#include "TLorentzVector.h"
45
46#include <algorithm>
47#include <stdexcept>
48#include <iostream>
49#include <sstream>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55JetPileUpSubtractor::JetPileUpSubtractor() :
56 fItJetInputArray(0), fItRhoInputArray(0)
57{
58
59}
60
61//------------------------------------------------------------------------------
62
63JetPileUpSubtractor::~JetPileUpSubtractor()
64{
65
66}
67
68//------------------------------------------------------------------------------
69
70void JetPileUpSubtractor::Init()
71{
72 fJetPTMin = GetDouble("JetPTMin", 20.0);
73
74 // import input array(s)
75
76 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
77 fItJetInputArray = fJetInputArray->MakeIterator();
78
79 fRhoInputArray = ImportArray(GetString("RhoInputArray", "Rho/rho"));
80 fItRhoInputArray = fRhoInputArray->MakeIterator();
81
82 // create output array(s)
83
84 fOutputArray = ExportArray(GetString("OutputArray", "jets"));
85
86}
87
88//------------------------------------------------------------------------------
89
90void JetPileUpSubtractor::Finish()
91{
92 if(fItRhoInputArray) delete fItRhoInputArray;
93 if(fItJetInputArray) delete fItJetInputArray;
94}
95
96//------------------------------------------------------------------------------
97
98void JetPileUpSubtractor::Process()
99{
100 Candidate *candidate, *object;
101 TLorentzVector momentum, area;
102 Double_t eta = 0.0;
103 Double_t rho = 0.0;
104
105 // loop over all input candidates
106 fItJetInputArray->Reset();
107 while((candidate = static_cast<Candidate*>(fItJetInputArray->Next())))
108 {
109 momentum = candidate->Momentum;
110 area = candidate->Area;
111 eta = momentum.Eta();
112
113 // find rho
114 rho = 0.0;
115 if(fRhoInputArray)
116 {
117 fItRhoInputArray->Reset();
118 while((object = static_cast<Candidate*>(fItRhoInputArray->Next())))
119 {
120 if(eta >= object->Edges[0] && eta < object->Edges[1])
121 {
122 rho = object->Momentum.Pt();
123 }
124 }
125 }
126
127 // apply pile-up correction
128 if(momentum.Pt() <= rho * area.Pt()) continue;
129
130 momentum -= rho * area;
131
132 if(momentum.Pt() <= fJetPTMin) continue;
133
134 candidate = static_cast<Candidate*>(candidate->Clone());
135 candidate->Momentum = momentum;
136
137 fOutputArray->Add(candidate);
138 }
139}
140
141//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.