Fork me on GitHub

source: git/modules/TaggingParticlesSkimmer.cc@ 190cfa0

ImprovedOutputFile Timing
Last change on this file since 190cfa0 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 4.6 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 fItPartonInputArray(0), fFilter(0), fClassifier(0)
60{
61}
62
63//------------------------------------------------------------------------------
64
65TaggingParticlesSkimmer::~TaggingParticlesSkimmer()
66{
67}
68
69//------------------------------------------------------------------------------
70
71void TaggingParticlesSkimmer::Init()
72{
73
74 fPTMin = GetDouble("PTMin", 15.0);
75 fEtaMax = GetDouble("EtaMax", 2.5);
76
77 // import input array
78 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
79 fItPartonInputArray = fPartonInputArray->MakeIterator();
80
81 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
82
83 fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
84 fClassifier->fPTMin = GetDouble("PTMin", 15.0);
85 fClassifier->fEtaMax = GetDouble("EtaMax", 2.5);
86
87 fFilter = new ExRootFilter(fPartonInputArray);
88
89 // output array
90 fOutputArray = ExportArray(GetString("OutputArray", "taggingParticles"));
91}
92
93//------------------------------------------------------------------------------
94
95void TaggingParticlesSkimmer::Finish()
96{
97 if(fItPartonInputArray) delete fItPartonInputArray;
98 if(fFilter) delete fFilter;
99 if(fClassifier) delete fClassifier;
100}
101
102//------------------------------------------------------------------------------
103
104void TaggingParticlesSkimmer::Process()
105{
106 Candidate *candidate, *tau, *daughter;
107 TLorentzVector tauMomentum;
108 Double_t pt, eta;
109 TObjArray *tauArray;
110 Int_t pdgCode, i;
111
112 // first select hadronic taus and replace them by visible part
113 fFilter->Reset();
114 tauArray = fFilter->GetSubArray(fClassifier, 0);
115
116 if(tauArray == 0) return;
117
118 TIter itTauArray(tauArray);
119
120 // loop over all input taus
121 itTauArray.Reset();
122 while((tau = static_cast<Candidate *>(itTauArray.Next())))
123 {
124 if(tau->D1 < 0) continue;
125
126 if(tau->D1 >= fParticleInputArray->GetEntriesFast() || tau->D2 >= fParticleInputArray->GetEntriesFast())
127 {
128 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
129 }
130
131 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
132
133 for(i = tau->D1; i <= tau->D2; ++i)
134 {
135 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
136 if(TMath::Abs(daughter->PID) == 16) continue;
137 tauMomentum += daughter->Momentum;
138 }
139
140 candidate = static_cast<Candidate *>(tau->Clone());
141 candidate->Momentum = tauMomentum;
142
143 fOutputArray->Add(candidate);
144 }
145
146 // then add all other partons (except tau's to avoid double counting)
147
148 fItPartonInputArray->Reset();
149 while((candidate = static_cast<Candidate *>(fItPartonInputArray->Next())))
150 {
151 pdgCode = TMath::Abs(candidate->PID);
152 if(pdgCode == 15) continue;
153
154 pt = candidate->Momentum.Pt();
155 if(pt < fPTMin) continue;
156
157 eta = TMath::Abs(candidate->Momentum.Eta());
158 if(eta > fEtaMax) continue;
159
160 fOutputArray->Add(candidate);
161 }
162}
Note: See TracBrowser for help on using the repository browser.