Fork me on GitHub

source: git/modules/AngularSmearing.cc@ cd75093

ImprovedOutputFile Timing dual_readout llp
Last change on this file since cd75093 was edcf81a, checked in by Michele Selvaggi <michele.selvaggi@…>, 10 years ago

Added AngularSmearing module, and modified FCC card

  • Property mode set to 100644
File size: 2.9 KB
Line 
1
2/** \class AngularSmearing
3 *
4 * Performs transverse angular resolution smearing.
5 *
6 * $Date: 2014-06-17 16:58:53 +0100 $
7 *
8 *
9 *
10 * \author M. Selvaggi - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include "modules/AngularSmearing.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
41AngularSmearing::AngularSmearing() :
42 fFormulaEta(0), fFormulaPhi(0), fItInputArray(0)
43{
44 fFormulaEta = new DelphesFormula;
45 fFormulaPhi = new DelphesFormula;
46}
47
48//------------------------------------------------------------------------------
49
50AngularSmearing::~AngularSmearing()
51{
52 if(fFormulaEta) delete fFormulaEta;
53 if(fFormulaPhi) delete fFormulaPhi;
54}
55
56//------------------------------------------------------------------------------
57
58void AngularSmearing::Init()
59{
60 // read resolution formula
61
62 fFormulaEta->Compile(GetString("EtaResolutionFormula", "0.0"));
63 fFormulaPhi->Compile(GetString("PhiResolutionFormula", "0.0"));
64
65
66 // import input array
67
68 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
69 fItInputArray = fInputArray->MakeIterator();
70
71 // create output array
72
73 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
74}
75
76//------------------------------------------------------------------------------
77
78void AngularSmearing::Finish()
79{
80 if(fItInputArray) delete fItInputArray;
81}
82
83//------------------------------------------------------------------------------
84
85void AngularSmearing::Process()
86{
87 Candidate *candidate, *mother;
88 Double_t pt, eta, phi;
89
90 fItInputArray->Reset();
91 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
92 {
93 const TLorentzVector &candidatePosition = candidate->Position;
94 const TLorentzVector &candidateMomentum = candidate->Momentum;
95 eta = candidatePosition.Eta();
96 phi = candidatePosition.Phi();
97 pt = candidateMomentum.Pt();
98
99 // apply smearing formula for eta,phi
100
101 eta = gRandom->Gaus(eta, fFormulaEta->Eval(pt, eta));
102 phi = gRandom->Gaus(phi, fFormulaPhi->Eval(pt, eta));
103
104 if(pt <= 0.0) continue;
105
106 mother = candidate;
107 candidate = static_cast<Candidate*>(candidate->Clone());
108 eta = candidateMomentum.Eta();
109 phi = candidateMomentum.Phi();
110 candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt*TMath::CosH(eta));
111 candidate->AddCandidate(mother);
112
113 fOutputArray->Add(candidate);
114 }
115}
116
117//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.