Fork me on GitHub

source: git/modules/TaggingParticlesSkimmer.cc@ 7993cad

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 7993cad was 4a93785, checked in by Michele Selvaggi <michele.selvaggi@…>, 10 years ago

Inverted eta check (#603)

  • 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/ExRootResult.h"
38#include "ExRootAnalysis/ExRootFilter.h"
39#include "ExRootAnalysis/ExRootClassifier.h"
40
41#include "TMath.h"
42#include "TString.h"
43#include "TFormula.h"
44#include "TRandom3.h"
45#include "TObjArray.h"
46#include "TDatabasePDG.h"
47#include "TLorentzVector.h"
48
49#include <algorithm>
50#include <stdexcept>
51#include <iostream>
52#include <sstream>
53
54using namespace std;
55
56
57//------------------------------------------------------------------------------
58
59TaggingParticlesSkimmer::TaggingParticlesSkimmer() :
60 fItPartonInputArray(0), fFilter(0), fClassifier(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
85 fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
86 fClassifier->fPTMin = GetDouble("PTMin", 15.0);
87 fClassifier->fEtaMax = GetDouble("EtaMax", 2.5);
88
89
90 fFilter = new ExRootFilter(fPartonInputArray);
91
92 // output array
93 fOutputArray = ExportArray(GetString("OutputArray", "taggingParticles"));
94}
95
96//------------------------------------------------------------------------------
97
98void TaggingParticlesSkimmer::Finish()
99{
100 if(fItPartonInputArray) delete fItPartonInputArray;
101 if(fFilter) delete fFilter;
102 if(fClassifier) delete fClassifier;
103}
104
105//------------------------------------------------------------------------------
106
107void TaggingParticlesSkimmer::Process()
108{
109 Candidate *candidate, *tau, *daughter;
110 TLorentzVector tauMomentum;
111 Double_t pt, eta;
112 TObjArray *tauArray;
113 Int_t pdgCode, i;
114
115 // first select hadronic taus and replace them by visible part
116 fFilter->Reset();
117 tauArray = fFilter->GetSubArray(fClassifier, 0);
118
119 if(tauArray == 0) return;
120
121 TIter itTauArray(tauArray);
122
123 // loop over all input taus
124 itTauArray.Reset();
125 while((tau = static_cast<Candidate *>(itTauArray.Next())))
126 {
127 if(tau->D1 < 0) continue;
128
129 if(tau->D1 >= fParticleInputArray->GetEntriesFast() ||
130 tau->D2 >= fParticleInputArray->GetEntriesFast())
131 {
132 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
133 }
134
135 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
136
137 for(i = tau->D1; i <= tau->D2; ++i)
138 {
139 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
140 if(TMath::Abs(daughter->PID) == 16) continue;
141 tauMomentum += daughter->Momentum;
142 }
143
144 candidate = static_cast<Candidate*>(tau->Clone());
145 candidate->Momentum = tauMomentum;
146
147
148 fOutputArray->Add(candidate);
149
150 }
151
152 // then add all other partons (except tau's to avoid double counting)
153
154 fItPartonInputArray->Reset();
155 while((candidate = static_cast<Candidate*>(fItPartonInputArray->Next())))
156 {
157 pdgCode = TMath::Abs(candidate->PID);
158 if(pdgCode == 15) continue;
159
160 pt = candidate->Momentum.Pt();
161 if(pt < fPTMin) continue;
162
163 eta = TMath::Abs(candidate->Momentum.Eta());
164 if(eta > fEtaMax) continue;
165
166 fOutputArray->Add(candidate);
167 }
168
169
170}
171
Note: See TracBrowser for help on using the repository browser.