Fork me on GitHub

source: git/modules/TimeSmearing.cc@ e70228d

Timing
Last change on this file since e70228d was 7939c6c, checked in by Kaan Yüksel Oyulmaz <kaanyukseloyulmaz@…>, 5 years ago

TimeSmearing Module is updated for neutral particles

  • Property mode set to 100644
File size: 4.2 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() :
[7939c6c]55 fItInputArray(0)
[af88093]56{
57}
58
59//------------------------------------------------------------------------------
60
61TimeSmearing::~TimeSmearing()
[7939c6c]62{
[af88093]63}
64
65//------------------------------------------------------------------------------
66
67void TimeSmearing::Init()
68{
69 // read resolution formula
70
[7939c6c]71 fTimeResolution = GetDouble("TimeResolution", 1.);
[af88093]72
[2b5ff2c]73 // import input array
74 fEtaMax = GetDouble("EtaMax", 6.);
[af88093]75 fInputArray = ImportArray(GetString("InputArray", "MuonMomentumSmearing/muons"));
76 fItInputArray = fInputArray->MakeIterator();
77 // create output array
78
79 fOutputArray = ExportArray(GetString("OutputArray", "muons"));
80}
81
82//------------------------------------------------------------------------------
83
84void TimeSmearing::Finish()
85{
86 if(fItInputArray) delete fItInputArray;
87}
88
89//------------------------------------------------------------------------------
90
91void TimeSmearing::Process()
92{
93 Candidate *candidate, *mother;
[7939c6c]94 Double_t ti, tf_smeared, tf;
95 Double_t pt, eta, phi, e, p, l;
96 Double_t sigma_t, beta_particle;
[2b5ff2c]97
98
[af88093]99 const Double_t c_light = 2.99792458E8;
[6cdc544]100
[7939c6c]101 cout << " STARTINNNGG ---------->" << endl;
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();
[7939c6c]116 p = candidateMomentum.P();
117 beta_particle = p/e;
118 l = candidate->L;
[af88093]119 // apply smearing formula
[7939c6c]120
121 if(candidate->Charge != 0)
[2b5ff2c]122 {
[7939c6c]123 tf_smeared = tf + fTimeResolution*gRandom->Gaus(0, 1);
124 //mother = candidate;
125 //candidate = static_cast<Candidate*>(candidate->Clone()); // I am not sure that we need these lines !!!
126 //candidate->AddCandidate(mother);
127 candidate->InitialPosition.SetT((100+ti)*1.0E3*c_light);
128 candidate->Position.SetT(tf_smeared*1.0E3*c_light);
129 candidate->ErrorT = fTimeResolution*1.0E3*c_light;
130 fOutputArray->Add(candidate);
[2b5ff2c]131 }
[7939c6c]132 else
133 {
134 sigma_t = sqrt(pow(20,2) + pow(150,2)/e);
135 ti = sigma_t - l*1.0E3/(c_light*beta_particle);
136 candidate->InitialPosition.SetT(ti);
137 candidate->ErrorT = sigma_t*1.0E3*c_light; // Do we need to sum with 100 like in upside ?
138 fOutputArray->Add(candidate);
139 }
140
[af88093]141 }
142}
143
144//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.