Fork me on GitHub

source: git/modules/EnergyScale.cc@ ad3b7ce

ImprovedOutputFile Timing dual_readout llp
Last change on this file since ad3b7ce was b443089, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix EOL characters in GPLv3 header

  • Property mode set to 100644
File size: 3.0 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 EnergyScale
21 *
22 * Applies energy scale.
23 *
24 * $Date$
25 * $Revision$
26 *
27 *
28 * \author P. Demin - UCL, Louvain-la-Neuve
29 *
30 */
31
32#include "modules/EnergyScale.h"
33
34#include "classes/DelphesClasses.h"
35#include "classes/DelphesFactory.h"
36#include "classes/DelphesFormula.h"
37
38#include "ExRootAnalysis/ExRootResult.h"
39#include "ExRootAnalysis/ExRootFilter.h"
40#include "ExRootAnalysis/ExRootClassifier.h"
41
42#include "TMath.h"
43#include "TString.h"
44#include "TFormula.h"
45#include "TRandom3.h"
46#include "TObjArray.h"
47#include "TDatabasePDG.h"
48#include "TLorentzVector.h"
49
50#include <algorithm>
51#include <stdexcept>
52#include <iostream>
53#include <sstream>
54
55using namespace std;
56
57//------------------------------------------------------------------------------
58
59EnergyScale::EnergyScale() :
60 fFormula(0), fItInputArray(0)
61{
62 fFormula = new DelphesFormula;
63}
64
65//------------------------------------------------------------------------------
66
67EnergyScale::~EnergyScale()
68{
69 if(fFormula) delete fFormula;
70}
71
72//------------------------------------------------------------------------------
73
74void EnergyScale::Init()
75{
76 // read resolution formula
77
78 fFormula->Compile(GetString("ScaleFormula", "0.0"));
79
80 // import input array
81
82 fInputArray = ImportArray(GetString("InputArray", "FastJetFinder/jets"));
83 fItInputArray = fInputArray->MakeIterator();
84
85 // create output array
86
87 fOutputArray = ExportArray(GetString("OutputArray", "jets"));
88}
89
90//------------------------------------------------------------------------------
91
92void EnergyScale::Finish()
93{
94 if(fItInputArray) delete fItInputArray;
95}
96
97//------------------------------------------------------------------------------
98
99void EnergyScale::Process()
100{
101 Candidate *candidate;
102 TLorentzVector momentum;
103 Double_t scale;
104
105 fItInputArray->Reset();
106 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
107 {
108 momentum = candidate->Momentum;
109
110 scale = fFormula->Eval(momentum.Pt(), momentum.Eta());
111
112 if(scale > 0.0) momentum *= scale;
113
114 candidate = static_cast<Candidate*>(candidate->Clone());
115 candidate->Momentum = momentum;
116
117 fOutputArray->Add(candidate);
118 }
119}
120
121//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.