Fork me on GitHub

source: git/modules/MomentumSmearing.cc@ 08d82fe

Last change on this file since 08d82fe was fd4b326, checked in by michele <michele.selvaggi@…>, 4 years ago

tracks have now montecarlo mass

  • Property mode set to 100644
File size: 4.1 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[d7d2da3]19/** \class MomentumSmearing
20 *
21 * Performs transverse momentum resolution smearing.
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/MomentumSmearing.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
[341014c]34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
[d7d2da3]36
37#include "TDatabasePDG.h"
[341014c]38#include "TFormula.h"
[d7d2da3]39#include "TLorentzVector.h"
[341014c]40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
[d7d2da3]44
[341014c]45#include <algorithm>
[d7d2da3]46#include <iostream>
47#include <sstream>
[341014c]48#include <stdexcept>
[d7d2da3]49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54MomentumSmearing::MomentumSmearing() :
55 fFormula(0), fItInputArray(0)
56{
57 fFormula = new DelphesFormula;
58}
59
60//------------------------------------------------------------------------------
61
62MomentumSmearing::~MomentumSmearing()
63{
64 if(fFormula) delete fFormula;
65}
66
67//------------------------------------------------------------------------------
68
69void MomentumSmearing::Init()
70{
71 // read resolution formula
72
73 fFormula->Compile(GetString("ResolutionFormula", "0.0"));
74
75 // import input array
76
77 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // create output array
81
82 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
83}
84
85//------------------------------------------------------------------------------
86
87void MomentumSmearing::Finish()
88{
89 if(fItInputArray) delete fItInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void MomentumSmearing::Process()
95{
96 Candidate *candidate, *mother;
[fd4b326]97 Double_t pt, eta, phi, e, m, res;
[d7d2da3]98
99 fItInputArray->Reset();
[341014c]100 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
[d7d2da3]101 {
102 const TLorentzVector &candidatePosition = candidate->Position;
103 const TLorentzVector &candidateMomentum = candidate->Momentum;
104 eta = candidatePosition.Eta();
105 phi = candidatePosition.Phi();
106 pt = candidateMomentum.Pt();
[95aa610]107 e = candidateMomentum.E();
[fd4b326]108 m = candidateMomentum.M();
[f298e48]109 res = fFormula->Eval(pt, eta, phi, e, candidate);
[341014c]110
[d7d2da3]111 // apply smearing formula
[bc58cf5]112 //pt = gRandom->Gaus(pt, fFormula->Eval(pt, eta, phi, e) * pt);
[79d4123]113
[341014c]114 res = (res > 1.0) ? 1.0 : res;
115
116 pt = LogNormal(pt, res * pt);
117
[bc58cf5]118 //if(pt <= 0.0) continue;
[d7d2da3]119
120 mother = candidate;
[341014c]121 candidate = static_cast<Candidate *>(candidate->Clone());
[d7d2da3]122 eta = candidateMomentum.Eta();
123 phi = candidateMomentum.Phi();
[fd4b326]124 candidate->Momentum.SetPtEtaPhiM(pt, eta, phi, m);
[4827699]125 //candidate->TrackResolution = fFormula->Eval(pt, eta, phi, e);
126 candidate->TrackResolution = res;
[d7d2da3]127 candidate->AddCandidate(mother);
[341014c]128
[d7d2da3]129 fOutputArray->Add(candidate);
130 }
131}
[bc58cf5]132//----------------------------------------------------------------
133
134Double_t MomentumSmearing::LogNormal(Double_t mean, Double_t sigma)
135{
136 Double_t a, b;
137
138 if(mean > 0.0)
139 {
[341014c]140 b = TMath::Sqrt(TMath::Log((1.0 + (sigma * sigma) / (mean * mean))));
141 a = TMath::Log(mean) - 0.5 * b * b;
[bc58cf5]142
[341014c]143 return TMath::Exp(a + b * gRandom->Gaus(0.0, 1.0));
[bc58cf5]144 }
145 else
146 {
147 return 0.0;
148 }
149}
150
[d7d2da3]151//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.