Fork me on GitHub

source: git/modules/TaggingParticlesSkimmer.cc@ 952bbbc

Last change on this file since 952bbbc was 1a4956e, checked in by Pavel Demin <pavel-demin@…>, 5 years ago

reorder member initializer lists

  • Property mode set to 100644
File size: 4.7 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//------------------------------------------------------------------------------
20
21/** \class TaggingParticlesSkimmer
22 *
23 * Filters particle collection by only keeping gen particles useful for b and tau tagging.
24 * tau particles are replaced by their "visible decay".
25 *
26 * \author M. Selvaggi
27 *
28 */
29
30#include "modules/TaggingParticlesSkimmer.h"
31#include "modules/TauTagging.h"
32
33#include "classes/DelphesClasses.h"
34#include "classes/DelphesFactory.h"
35#include "classes/DelphesFormula.h"
36
37#include "ExRootAnalysis/ExRootClassifier.h"
38#include "ExRootAnalysis/ExRootFilter.h"
39#include "ExRootAnalysis/ExRootResult.h"
40
41#include "TDatabasePDG.h"
42#include "TFormula.h"
43#include "TLorentzVector.h"
44#include "TMath.h"
45#include "TObjArray.h"
46#include "TRandom3.h"
47#include "TString.h"
48
49#include <algorithm>
50#include <iostream>
51#include <sstream>
52#include <stdexcept>
53
54using namespace std;
55
56//------------------------------------------------------------------------------
57
58TaggingParticlesSkimmer::TaggingParticlesSkimmer() :
59 fClassifier(0), fFilter(0), fItPartonInputArray(0),
60 fPartonInputArray(0), fParticleInputArray(0), fOutputArray(0)
61{
62}
63
64//------------------------------------------------------------------------------
65
66TaggingParticlesSkimmer::~TaggingParticlesSkimmer()
67{
68}
69
70//------------------------------------------------------------------------------
71
72void TaggingParticlesSkimmer::Init()
73{
74
75 fPTMin = GetDouble("PTMin", 15.0);
76 fEtaMax = GetDouble("EtaMax", 2.5);
77
78 // import input array
79 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
80 fItPartonInputArray = fPartonInputArray->MakeIterator();
81
82 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
83
84 fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
85 fClassifier->fPTMin = GetDouble("PTMin", 15.0);
86 fClassifier->fEtaMax = GetDouble("EtaMax", 2.5);
87
88 fFilter = new ExRootFilter(fPartonInputArray);
89
90 // output array
91 fOutputArray = ExportArray(GetString("OutputArray", "taggingParticles"));
92}
93
94//------------------------------------------------------------------------------
95
96void TaggingParticlesSkimmer::Finish()
97{
98 if(fItPartonInputArray) delete fItPartonInputArray;
99 if(fFilter) delete fFilter;
100 if(fClassifier) delete fClassifier;
101}
102
103//------------------------------------------------------------------------------
104
105void TaggingParticlesSkimmer::Process()
106{
107 Candidate *candidate, *tau, *daughter;
108 TLorentzVector tauMomentum;
109 Double_t pt, eta;
110 TObjArray *tauArray;
111 Int_t pdgCode, i;
112
113 // first select hadronic taus and replace them by visible part
114 fFilter->Reset();
115 tauArray = fFilter->GetSubArray(fClassifier, 0);
116
117 if(tauArray == 0) return;
118
119 TIter itTauArray(tauArray);
120
121 // loop over all input taus
122 itTauArray.Reset();
123 while((tau = static_cast<Candidate *>(itTauArray.Next())))
124 {
125 if(tau->D1 < 0) continue;
126
127 if(tau->D1 >= fParticleInputArray->GetEntriesFast() || tau->D2 >= fParticleInputArray->GetEntriesFast())
128 {
129 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
130 }
131
132 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
133
134 for(i = tau->D1; i <= tau->D2; ++i)
135 {
136 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
137 if(TMath::Abs(daughter->PID) == 16) continue;
138 tauMomentum += daughter->Momentum;
139 }
140
141 candidate = static_cast<Candidate *>(tau->Clone());
142 candidate->Momentum = tauMomentum;
143
144 fOutputArray->Add(candidate);
145 }
146
147 // then add all other partons (except tau's to avoid double counting)
148
149 fItPartonInputArray->Reset();
150 while((candidate = static_cast<Candidate *>(fItPartonInputArray->Next())))
151 {
152 pdgCode = TMath::Abs(candidate->PID);
153 if(pdgCode == 15) continue;
154
155 pt = candidate->Momentum.Pt();
156 if(pt < fPTMin) continue;
157
158 eta = TMath::Abs(candidate->Momentum.Eta());
159 if(eta > fEtaMax) continue;
160
161 fOutputArray->Add(candidate);
162 }
163}
Note: See TracBrowser for help on using the repository browser.