Fork me on GitHub

source: git/modules/DecayFilter.cc@ d4fb268

Last change on this file since d4fb268 was 4d3fb73, checked in by Roberto Preghenella <preghenella@…>, 5 years ago

Implementation of a DecayFilter module

  • Property mode set to 100644
File size: 3.9 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
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.
9 *
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.
14 *
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
19/** \class DecayFilter
20 *
21 * This module randomly generates decays along the particle trajectory length
22 * according to actual particle decay length, taking into account for the boost
23 * and using ROOT TDatabasePDG as a source for the particle lifetime.
24 *
25 * This module is to to be used after a PropagateParticle step or a similar module
26 * that calculates and store a trajectory length.
27 *
28 * Particles that decay are not added to the OutputArray.
29 *
30 * \author R. Preghenella - INFN, Bologna
31 *
32 */
33
34#include "modules/DecayFilter.h"
35
36#include "classes/DelphesClasses.h"
37#include "classes/DelphesFactory.h"
38#include "classes/DelphesFormula.h"
39
40#include "ExRootAnalysis/ExRootClassifier.h"
41#include "ExRootAnalysis/ExRootFilter.h"
42#include "ExRootAnalysis/ExRootResult.h"
43
44#include "TDatabasePDG.h"
45#include "TParticlePDG.h"
46#include "TFormula.h"
47#include "TLorentzVector.h"
48#include "TMath.h"
49#include "TObjArray.h"
50#include "TRandom3.h"
51#include "TString.h"
52
53#include <algorithm>
54#include <iostream>
55#include <sstream>
56#include <stdexcept>
57
58using namespace std;
59
60//------------------------------------------------------------------------------
61
62DecayFilter::DecayFilter() :
63 fItInputArray(0)
64{}
65
66//------------------------------------------------------------------------------
67
68DecayFilter::~DecayFilter()
69{}
70
71//------------------------------------------------------------------------------
72
73void DecayFilter::Init()
74{
75
76 // import input array(s)
77
78 fInputArray = ImportArray(GetString("InputArray", "FastJetFinder/jets"));
79 fItInputArray = fInputArray->MakeIterator();
80
81 // create output array(s)
82
83 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
84}
85
86//------------------------------------------------------------------------------
87
88void DecayFilter::Finish()
89{
90 if(fItInputArray) delete fItInputArray;
91}
92
93//------------------------------------------------------------------------------
94
95void DecayFilter::Process()
96{
97 Candidate *candidate;
98 TDatabasePDG *pdgdb = TDatabasePDG::Instance();
99 const Double_t c = TMath::C(); // [m/s]
100 Double_t m, t, p, bgct, L, l;
101 Bool_t hasDecayed = kFALSE;
102
103 // loop over all input candidates
104 fItInputArray->Reset();
105 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
106 {
107 // get particle information from PDG
108 TParticlePDG *pdg = pdgdb->GetParticle(candidate->PID);
109 if (!pdg) { // don't know this particle
110 fOutputArray->Add(candidate);
111 continue;
112 }
113 m = pdg->Mass();
114 t = pdg->Lifetime(); // [s]
115 if (t == 0.) { // does not decay
116 fOutputArray->Add(candidate);
117 continue;
118 }
119
120 // compute boosted decay length (beta gamma c tau)
121 p = candidate->P;
122 bgct = p / m * c * t; // [m]
123
124 // get full trajectory length and generate random decay length
125 L = candidate->L * 1.0E-3; // [m]
126 l = gRandom->Exp(bgct);
127
128 // if random decay happens before end of trajectory, reject track
129 if (l < L) continue;
130
131 // else particle did not decay within the trajectory
132 fOutputArray->Add(candidate);
133 }
134}
135
136//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.