Fork me on GitHub

source: git/modules/DenseTrackFilter.cc@ 8a58fff

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 8a58fff was b24f05a, checked in by Michele Selvaggi <michele.selvaggi@…>, 6 years ago

push PhaseII final card

  • Property mode set to 100644
File size: 7.8 KB
RevLine 
[196994a]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
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55DenseTrackFilter::DenseTrackFilter() :
[ec37bc3]56 fItTrackInputArray(0)
[196994a]57{
58}
59
60//------------------------------------------------------------------------------
61
62DenseTrackFilter::~DenseTrackFilter()
63{
64}
65
66//------------------------------------------------------------------------------
67
68void 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
[4170544]113 fEtaPhiRes = GetDouble("EtaPhiRes", 0.003);
[196994a]114
[4170544]115 fTrackInputArray = ImportArray(GetString("TrackInputArray", "TrackMergerProp/tracks"));
[196994a]116 fItTrackInputArray = fTrackInputArray->MakeIterator();
117
118 fTrackOutputArray = ExportArray(GetString("TrackOutputArray", "tracks"));
[ec37bc3]119 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
120 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
121 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
[196994a]122}
123
124//------------------------------------------------------------------------------
125
126void DenseTrackFilter::Finish()
127{
128 vector< vector< Double_t >* >::iterator itPhiBin;
129 if(fItTrackInputArray) delete fItTrackInputArray;
130 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
131 {
132 delete *itPhiBin;
133 }
134}
135
136//------------------------------------------------------------------------------
137
138void DenseTrackFilter::Process()
139{
[4170544]140 Candidate *track;
[196994a]141 TLorentzVector position, momentum;
142 Short_t etaBin, phiBin, flags;
[4170544]143 Int_t number;
[196994a]144 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
[4170544]145 Double_t ptmax;
[196994a]146
147 vector< Double_t >::iterator itEtaBin;
148 vector< Double_t >::iterator itPhiBin;
149 vector< Double_t > *phiBins;
150
151 vector< Long64_t >::iterator itTowerHits;
152
153 fTowerHits.clear();
154
155 // loop over all tracks
[ec37bc3]156 fItTrackInputArray->Reset();
[196994a]157 number = -1;
[ec37bc3]158 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
[196994a]159 {
160 const TLorentzVector &trackPosition = track->Position;
161 ++number;
162
163 // find eta bin [1, fEtaBins.size - 1]
164 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
165 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
166 etaBin = distance(fEtaBins.begin(), itEtaBin);
167
168 // phi bins for given eta bin
169 phiBins = fPhiBins[etaBin];
170
171 // find phi bin [1, phiBins.size - 1]
172 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
173 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
174 phiBin = distance(phiBins->begin(), itPhiBin);
175
176 flags = 1;
177
178 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
179 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
180
181 fTowerHits.push_back(towerHit);
182 }
183
184 // all hits are sorted first by eta bin number, then by phi bin number,
185 // then by flags and then by particle or track number
186 sort(fTowerHits.begin(), fTowerHits.end());
187
188 // loop over all hits
189 towerEtaPhi = 0;
[4170544]190 fBestTrack = 0;
[196994a]191 ptmax = 0.0;
[4170544]192 fTowerTrackHits = 0;
[196994a]193
194 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
195 {
196 towerHit = (*itTowerHits);
197 flags = (towerHit >> 24) & 0x00000000000000FFLL;
198 number = (towerHit) & 0x0000000000FFFFFFLL;
199 hitEtaPhi = towerHit >> 32;
[4170544]200
[196994a]201 if(towerEtaPhi != hitEtaPhi)
202 {
203 // switch to next tower
204 towerEtaPhi = hitEtaPhi;
[60e1de6]205
[196994a]206 // saving track with highest pT that hit the tower
[4170544]207 FillTrack();
[196994a]208
209 ptmax = 0.0;
[4170544]210 fTowerTrackHits = 0;
211 fBestTrack = 0;
[196994a]212 }
213 // check for track hits
[60e1de6]214
[196994a]215 if(flags & 1)
216 {
[4170544]217 ++fTowerTrackHits;
[ec37bc3]218 track = static_cast<Candidate*>(fTrackInputArray->At(number));
[196994a]219 momentum = track->Momentum;
220
[60e1de6]221 if(momentum.Pt() > ptmax)
[196994a]222 {
[60e1de6]223 ptmax = momentum.Pt();
[4170544]224 fBestTrack = track;
[196994a]225 }
226 continue;
227 }
[60e1de6]228
[196994a]229 }
[60e1de6]230
[4170544]231 // here fill last tower
232 FillTrack();
233
[60e1de6]234}
[4170544]235
236//------------------------------------------------------------------------------
237
238void DenseTrackFilter::FillTrack()
239{
[60e1de6]240
241 Candidate *candidate, *track;
[4170544]242 Double_t pt, eta, phi;
[60e1de6]243 Int_t numberOfCandidates;
[4170544]244
[60e1de6]245 // saving track with highest pT that hit the tower
[ec37bc3]246 if(fTowerTrackHits > 0)
[60e1de6]247 {
[ec37bc3]248
249 numberOfCandidates = fBestTrack->GetCandidates()->GetEntriesFast();
250 if (numberOfCandidates > 1)
251 {
252
253 track = static_cast<Candidate*>(fBestTrack->GetCandidates()->At(numberOfCandidates - 1));
254 candidate = static_cast<Candidate*>(track->Clone());
255 pt = candidate->Momentum.Pt();
256 eta = candidate->Momentum.Eta();
257 phi = candidate->Momentum.Phi();
258 eta = gRandom->Gaus(eta, fEtaPhiRes);
259 phi = gRandom->Gaus(phi, fEtaPhiRes);
260 candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt*TMath::CosH(eta));
261 candidate->AddCandidate(track);
262 fTrackOutputArray->Add(candidate);
263 switch(TMath::Abs(candidate->PID))
264 {
265 case 11:
266 fElectronOutputArray->Add(candidate);
267 break;
268 case 13:
269 fMuonOutputArray->Add(candidate);
270 break;
271 default:
272 fChargedHadronOutputArray->Add(candidate);
273 }
[b24f05a]274
[ec37bc3]275 }
[60e1de6]276 }
[ec37bc3]277
[196994a]278}
Note: See TracBrowser for help on using the repository browser.