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