Fork me on GitHub

source: git/modules/AngularSmearing.cc@ 7143a0e

Last change on this file since 7143a0e was f298e48, checked in by Roberto Preghenella <preghenella@…>, 5 years ago

Add radius (distance from beam line) to formula parameterisations

also redefines the DelphesFormula::Eval function to allow one to
include future parameters for parameterisations that can be
taken from Candidates in a generalised way.

  • 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/** \class AngularSmearing
20 *
21 * Performs transverse angular resolution smearing.
22 *
23 * \author M. Selvaggi - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/AngularSmearing.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
54AngularSmearing::AngularSmearing() :
55 fFormulaEta(0), fFormulaPhi(0), fItInputArray(0)
56{
57 fFormulaEta = new DelphesFormula;
58 fFormulaPhi = new DelphesFormula;
59}
60
61//------------------------------------------------------------------------------
62
63AngularSmearing::~AngularSmearing()
64{
65 if(fFormulaEta) delete fFormulaEta;
66 if(fFormulaPhi) delete fFormulaPhi;
67}
68
69//------------------------------------------------------------------------------
70
71void AngularSmearing::Init()
72{
73 // read resolution formula
74
75 fFormulaEta->Compile(GetString("EtaResolutionFormula", "0.0"));
76 fFormulaPhi->Compile(GetString("PhiResolutionFormula", "0.0"));
77
78 // import input array
79
80 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
81 fItInputArray = fInputArray->MakeIterator();
82
83 // create output array
84
85 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
86}
87
88//------------------------------------------------------------------------------
89
90void AngularSmearing::Finish()
91{
92 if(fItInputArray) delete fItInputArray;
93}
94
95//------------------------------------------------------------------------------
96
97void AngularSmearing::Process()
98{
99 Candidate *candidate, *mother;
100 Double_t pt, eta, phi, e;
101
102 fItInputArray->Reset();
103 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
104 {
105 const TLorentzVector &candidatePosition = candidate->Position;
106 const TLorentzVector &candidateMomentum = candidate->Momentum;
107 eta = candidatePosition.Eta();
108 phi = candidatePosition.Phi();
109 pt = candidateMomentum.Pt();
110 e = candidateMomentum.E();
111
112 // apply smearing formula for eta,phi
113
114 eta = gRandom->Gaus(eta, fFormulaEta->Eval(pt, eta, phi, e, candidate));
115 phi = gRandom->Gaus(phi, fFormulaPhi->Eval(pt, eta, phi, e, candidate));
116
117 if(pt <= 0.0) continue;
118
119 mother = candidate;
120 candidate = static_cast<Candidate *>(candidate->Clone());
121 eta = candidateMomentum.Eta();
122 phi = candidateMomentum.Phi();
123 candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt * TMath::CosH(eta));
124 candidate->AddCandidate(mother);
125
126 fOutputArray->Add(candidate);
127 }
128}
129
130//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.