Fork me on GitHub

source: git/modules/DenseTrackFilter.cc@ 4692fd9

Last change on this file since 4692fd9 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: 7.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/** \class DenseTrackFilter
20 *
21 * Applies reconstruction inefficiency on tracks in dense environment
22 *
23 * \author M. Selvaggi - CERN
24 *
25 */
26
27#include "modules/DenseTrackFilter.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
36
37#include "TDatabasePDG.h"
38#include "TFormula.h"
39#include "TLorentzVector.h"
40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
44
45#include <algorithm>
46#include <iostream>
47#include <sstream>
48#include <stdexcept>
49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54DenseTrackFilter::DenseTrackFilter() :
55 fItTrackInputArray(0)
56{
57}
58
59//------------------------------------------------------------------------------
60
61DenseTrackFilter::~DenseTrackFilter()
62{
63}
64
65//------------------------------------------------------------------------------
66
67void DenseTrackFilter::Init()
68{
69 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
70 Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
71 TBinMap::iterator itEtaBin;
72 set<Double_t>::iterator itPhiBin;
73 vector<Double_t> *phiBins;
74
75 // read eta and phi bins
76 param = GetParam("EtaPhiBins");
77 size = param.GetSize();
78 fBinMap.clear();
79 fEtaBins.clear();
80 fPhiBins.clear();
81 for(i = 0; i < size / 2; ++i)
82 {
83 paramEtaBins = param[i * 2];
84 sizeEtaBins = paramEtaBins.GetSize();
85 paramPhiBins = param[i * 2 + 1];
86 sizePhiBins = paramPhiBins.GetSize();
87
88 for(j = 0; j < sizeEtaBins; ++j)
89 {
90 for(k = 0; k < sizePhiBins; ++k)
91 {
92 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
93 }
94 }
95 }
96
97 // for better performance we transform map of sets to parallel vectors:
98 // vector< double > and vector< vector< double >* >
99 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
100 {
101 fEtaBins.push_back(itEtaBin->first);
102 phiBins = new vector<double>(itEtaBin->second.size());
103 fPhiBins.push_back(phiBins);
104 phiBins->clear();
105 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
106 {
107 phiBins->push_back(*itPhiBin);
108 }
109 }
110
111 // Eta x Phi smearing to be applied
112 fEtaPhiRes = GetDouble("EtaPhiRes", 0.003);
113
114 fTrackInputArray = ImportArray(GetString("TrackInputArray", "TrackMergerProp/tracks"));
115 fItTrackInputArray = fTrackInputArray->MakeIterator();
116
117 fTrackOutputArray = ExportArray(GetString("TrackOutputArray", "tracks"));
118 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
119 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
120 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
121}
122
123//------------------------------------------------------------------------------
124
125void DenseTrackFilter::Finish()
126{
127 vector<vector<Double_t> *>::iterator itPhiBin;
128 if(fItTrackInputArray) delete fItTrackInputArray;
129 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
130 {
131 delete *itPhiBin;
132 }
133}
134
135//------------------------------------------------------------------------------
136
137void DenseTrackFilter::Process()
138{
139 Candidate *track;
140 TLorentzVector position, momentum;
141 Short_t etaBin, phiBin, flags;
142 Int_t number;
143 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
144 Double_t ptmax;
145
146 vector<Double_t>::iterator itEtaBin;
147 vector<Double_t>::iterator itPhiBin;
148 vector<Double_t> *phiBins;
149
150 vector<Long64_t>::iterator itTowerHits;
151
152 fTowerHits.clear();
153
154 // loop over all tracks
155 fItTrackInputArray->Reset();
156 number = -1;
157 while((track = static_cast<Candidate *>(fItTrackInputArray->Next())))
158 {
159 const TLorentzVector &trackPosition = track->Position;
160 ++number;
161
162 // find eta bin [1, fEtaBins.size - 1]
163 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
164 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
165 etaBin = distance(fEtaBins.begin(), itEtaBin);
166
167 // phi bins for given eta bin
168 phiBins = fPhiBins[etaBin];
169
170 // find phi bin [1, phiBins.size - 1]
171 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
172 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
173 phiBin = distance(phiBins->begin(), itPhiBin);
174
175 flags = 1;
176
177 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
178 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
179
180 fTowerHits.push_back(towerHit);
181 }
182
183 // all hits are sorted first by eta bin number, then by phi bin number,
184 // then by flags and then by particle or track number
185 sort(fTowerHits.begin(), fTowerHits.end());
186
187 // loop over all hits
188 towerEtaPhi = 0;
189 fBestTrack = 0;
190 ptmax = 0.0;
191 fTowerTrackHits = 0;
192
193 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
194 {
195 towerHit = (*itTowerHits);
196 flags = (towerHit >> 24) & 0x00000000000000FFLL;
197 number = (towerHit)&0x0000000000FFFFFFLL;
198 hitEtaPhi = towerHit >> 32;
199
200 if(towerEtaPhi != hitEtaPhi)
201 {
202 // switch to next tower
203 towerEtaPhi = hitEtaPhi;
204
205 // saving track with highest pT that hit the tower
206 FillTrack();
207
208 ptmax = 0.0;
209 fTowerTrackHits = 0;
210 fBestTrack = 0;
211 }
212 // check for track hits
213
214 if(flags & 1)
215 {
216 ++fTowerTrackHits;
217 track = static_cast<Candidate *>(fTrackInputArray->At(number));
218 momentum = track->Momentum;
219
220 if(momentum.Pt() > ptmax)
221 {
222 ptmax = momentum.Pt();
223 fBestTrack = track;
224 }
225 continue;
226 }
227 }
228
229 // here fill last tower
230 FillTrack();
231}
232
233//------------------------------------------------------------------------------
234
235void DenseTrackFilter::FillTrack()
236{
237
238 Candidate *candidate, *track;
239 Double_t pt, eta, phi;
240 Int_t numberOfCandidates;
241
242 // saving track with highest pT that hit the tower
243 if(fTowerTrackHits < 1) return;
244
245 numberOfCandidates = fBestTrack->GetCandidates()->GetEntriesFast();
246 if(numberOfCandidates < 1) return;
247
248 track = static_cast<Candidate *>(fBestTrack->GetCandidates()->At(numberOfCandidates - 1));
249 candidate = static_cast<Candidate *>(track->Clone());
250 pt = candidate->Momentum.Pt();
251 eta = candidate->Momentum.Eta();
252 phi = candidate->Momentum.Phi();
253 eta = gRandom->Gaus(eta, fEtaPhiRes);
254 phi = gRandom->Gaus(phi, fEtaPhiRes);
255 candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt * TMath::CosH(eta));
256 candidate->AddCandidate(track);
257
258 fTrackOutputArray->Add(candidate);
259 switch(TMath::Abs(candidate->PID))
260 {
261 case 11:
262 fElectronOutputArray->Add(candidate);
263 break;
264 case 13:
265 fMuonOutputArray->Add(candidate);
266 break;
267 default:
268 fChargedHadronOutputArray->Add(candidate);
269 }
270}
Note: See TracBrowser for help on using the repository browser.