Fork me on GitHub

source: git/modules/TimeSmearing.cc@ 03b9c0f

Timing
Last change on this file since 03b9c0f was 2b5ff2c, checked in by Michele Selvaggi <michele.selvaggi@…>, 5 years ago

take Olmos Timesmearing and use TFormula as time resolution parameter

  • 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
[af88093]19/** \class TimeSmearing
20 *
21 * Performs transverse momentum resolution smearing.
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/TimeSmearing.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"
[af88093]36
37#include "TDatabasePDG.h"
[341014c]38#include "TFormula.h"
[af88093]39#include "TLorentzVector.h"
[341014c]40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
[af88093]44
[6cdc544]45#include <algorithm>
[af88093]46#include <iostream>
47#include <sstream>
[341014c]48#include <stdexcept>
[af88093]49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54TimeSmearing::TimeSmearing() :
[2b5ff2c]55 fFormula(0), fItInputArray(0)
[af88093]56{
[2b5ff2c]57 fFormula = new DelphesFormula;
[af88093]58}
59
60//------------------------------------------------------------------------------
61
62TimeSmearing::~TimeSmearing()
63{
[2b5ff2c]64 if(fFormula) delete fFormula;
[af88093]65}
66
67//------------------------------------------------------------------------------
68
69void TimeSmearing::Init()
70{
71 // read resolution formula
72
[2b5ff2c]73 fFormula->Compile(GetString("TimeResolution", "1.0"));
[af88093]74
[2b5ff2c]75 // import input array
76 fEtaMax = GetDouble("EtaMax", 6.);
[af88093]77 fInputArray = ImportArray(GetString("InputArray", "MuonMomentumSmearing/muons"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // create output array
81
82 fOutputArray = ExportArray(GetString("OutputArray", "muons"));
83}
84
85//------------------------------------------------------------------------------
86
87void TimeSmearing::Finish()
88{
89 if(fItInputArray) delete fItInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void TimeSmearing::Process()
95{
96 Candidate *candidate, *mother;
[2b5ff2c]97 Double_t ti, tf_smeared, tf, timeResolution;
98 Double_t pt, eta, phi, e, d0, dz, ctgTheta;
99
100
[af88093]101 const Double_t c_light = 2.99792458E8;
[6cdc544]102
[af88093]103 fItInputArray->Reset();
[2b5ff2c]104 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
[af88093]105 {
[2b5ff2c]106 ti = candidate->InitialPosition.T()*1.0E-3/c_light;
107 tf = candidate->Position.T()*1.0E-3/c_light;
108
109 // dummy, only need to properly call TFormula
110 const TLorentzVector &candidatePosition = candidate->Position;
111 const TLorentzVector &candidateMomentum = candidate->Momentum;
112 eta = candidatePosition.Eta();
113 phi = candidatePosition.Phi();
114 pt = candidateMomentum.Pt();
115 e = candidateMomentum.E();
116 d0 = candidate->D0;
117 dz = candidate->DZ;
118 ctgTheta = candidate->CtgTheta;
[6cdc544]119
[af88093]120 // apply smearing formula
[2b5ff2c]121
122 timeResolution = fFormula->Eval(pt, eta, phi, e, d0, dz, ctgTheta);
123 if(fabs(candidate->Position.Eta())<fEtaMax)
124 {
125 tf_smeared = tf + timeResolution*gRandom->Gaus(0, 1);
126 }
127 else continue;
[efc1546]128
[2b5ff2c]129 // double beta_particle = candidate->Momentum.P()/candidate->Momentum.E();
130 // ti = tf_smeared - candidate->Ld*1.0E-3/(c_light*beta_particle);
[6cdc544]131
[2b5ff2c]132 mother = candidate;
133 candidate = static_cast<Candidate*>(candidate->Clone());
[af88093]134 candidate->AddCandidate(mother);
[2b5ff2c]135 candidate->InitialPosition.SetT((100+ti)*1.0E3*c_light);
136 candidate->Position.SetT(tf_smeared*1.0E3*c_light);
137 candidate->ErrorT = timeResolution*1.0E3*c_light;
[6cdc544]138
[af88093]139 fOutputArray->Add(candidate);
140 }
141}
142
143//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.