Fork me on GitHub

source: git/modules/TimeSmearing.cc@ d612dec

Last change on this file since d612dec was 1ad8eca, checked in by michele <michele.selvaggi@…>, 3 years ago

fixed momentum cut and generalised to neutrals

  • Property mode set to 100644
File size: 3.5 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 *
[c27c038]21 * Performs time smearing.
[af88093]22 *
[c27c038]23 * \author M. Selvaggi - CERN
[af88093]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
53TimeSmearing::TimeSmearing() :
[1ad8eca]54 fItInputArray(0), fResolutionFormula(0)
[af88093]55{
[c27c038]56 fResolutionFormula = new DelphesFormula;
[af88093]57}
58
59//------------------------------------------------------------------------------
60
61TimeSmearing::~TimeSmearing()
62{
[c27c038]63 if(fResolutionFormula) delete fResolutionFormula;
[af88093]64}
65
66//------------------------------------------------------------------------------
67
68void TimeSmearing::Init()
69{
70 // read resolution formula
71
[c27c038]72 // read time resolution formula in seconds
73 fResolutionFormula->Compile(GetString("TimeResolution", "30e-12"));
[af88093]74
[c27c038]75 // import track input array
[1ad8eca]76 fInputArray = ImportArray(GetString("InputArray", "MuonMomentumSmearing/muons"));
77 fItInputArray = fInputArray->MakeIterator();
[af88093]78
79
[c27c038]80 // create output array
81 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
[af88093]82}
83
84//------------------------------------------------------------------------------
85
86void TimeSmearing::Finish()
87{
[1ad8eca]88 if(fItInputArray) delete fItInputArray;
[af88093]89}
90
91//------------------------------------------------------------------------------
92
93void TimeSmearing::Process()
94{
95 Candidate *candidate, *mother;
[c27c038]96 Double_t tf_smeared, tf;
97 Double_t eta, energy;
98 Double_t timeResolution;
99
[af88093]100 const Double_t c_light = 2.99792458E8;
[6cdc544]101
[1ad8eca]102 fItInputArray->Reset();
103 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
[04e2040]104 {
[c438758]105
[efc1546]106 const TLorentzVector &candidateFinalPosition = candidate->Position;
[c27c038]107 const TLorentzVector &candidateMomentum = candidate->Momentum;
[efc1546]108
[341014c]109 tf = candidateFinalPosition.T() * 1.0E-3 / c_light;
[04e2040]110
[c27c038]111 eta = candidateMomentum.Eta();
112 energy = candidateMomentum.E();
[6cdc544]113
[af88093]114 // apply smearing formula
[c27c038]115 timeResolution = fResolutionFormula->Eval(0.0, eta, 0.0, energy);
116 tf_smeared = gRandom->Gaus(tf, timeResolution);
[341014c]117
[af88093]118 mother = candidate;
[341014c]119 candidate = static_cast<Candidate *>(candidate->Clone());
[efc1546]120
[c27c038]121 candidate->Position.SetT(tf_smeared * 1.0E3 * c_light);
122 candidate->ErrorT = timeResolution * 1.0E3 * c_light;
[6cdc544]123
[af88093]124 candidate->AddCandidate(mother);
125 fOutputArray->Add(candidate);
126 }
127}
Note: See TracBrowser for help on using the repository browser.