Fork me on GitHub

source: git/modules/EnergyScale.cc@ c614dd7

ImprovedOutputFile Timing
Last change on this file since c614dd7 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 5 years ago

apply .clang-format to all .h, .cc and .cpp files

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