Fork me on GitHub

source: svn/trunk/modules/JetPileUpSubtractor.cc@ 1316

Last change on this file since 1316 was 1315, checked in by Pavel Demin, 11 years ago

add eta ranges for rho

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