Fork me on GitHub

source: svn/trunk/modules/EnergyScale.cc@ 1235

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

add EnergyScale

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id Revision Date
File size: 2.2 KB
Line 
1
2/** \class EnergyScale
3 *
4 * Performs transverse momentum resolution smearing.
5 *
6 * $Date: 2013-06-25 12:42:33 +0000 (Tue, 25 Jun 2013) $
7 * $Revision: 1141 $
8 *
9 *
10 * \author P. Demin - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include "modules/EnergyScale.h"
15
16#include "classes/DelphesClasses.h"
17#include "classes/DelphesFactory.h"
18#include "classes/DelphesFormula.h"
19
20#include "ExRootAnalysis/ExRootResult.h"
21#include "ExRootAnalysis/ExRootFilter.h"
22#include "ExRootAnalysis/ExRootClassifier.h"
23
24#include "TMath.h"
25#include "TString.h"
26#include "TFormula.h"
27#include "TRandom3.h"
28#include "TObjArray.h"
29#include "TDatabasePDG.h"
30#include "TLorentzVector.h"
31
32#include <algorithm>
33#include <stdexcept>
34#include <iostream>
35#include <sstream>
36
37using namespace std;
38
39//------------------------------------------------------------------------------
40
41EnergyScale::EnergyScale() :
42 fFormula(0), fItInputArray(0)
43{
44 fFormula = new DelphesFormula;
45}
46
47//------------------------------------------------------------------------------
48
49EnergyScale::~EnergyScale()
50{
51 if(fFormula) delete fFormula;
52}
53
54//------------------------------------------------------------------------------
55
56void EnergyScale::Init()
57{
58 // read resolution formula
59
60 fFormula->Compile(GetString("ScaleFormula", "0.0"));
61
62 // import input array
63
64 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
65 fItInputArray = fInputArray->MakeIterator();
66
67 // create output array
68
69 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
70}
71
72//------------------------------------------------------------------------------
73
74void EnergyScale::Finish()
75{
76 if(fItInputArray) delete fItInputArray;
77}
78
79//------------------------------------------------------------------------------
80
81void EnergyScale::Process()
82{
83 Candidate *candidate;
84 TLorentzVector momentum;
85 Double_t scale;
86
87 fItInputArray->Reset();
88 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
89 {
90 momentum = candidate->Momentum;
91
92 scale = fFormula->Eval(momentum.Pt(),momentum.Eta());
93
94 if(scale>0)momentum *= scale;
95
96 candidate = static_cast<Candidate*>(candidate->Clone());
97 candidate->Momentum = momentum;
98
99 fOutputArray->Add(candidate);
100
101 }
102}
103
104//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.