[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 |
|
---|
| 51 | using namespace std;
|
---|
| 52 |
|
---|
| 53 | //------------------------------------------------------------------------------
|
---|
| 54 |
|
---|
| 55 | DenseTrackFilter::DenseTrackFilter() :
|
---|
| 56 | fItTrackInputArray(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 = GetBool("EtaPhiRes", 0.005);
|
---|
| 114 |
|
---|
| 115 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
|
---|
| 116 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
| 117 |
|
---|
| 118 | fTrackOutputArray = ExportArray(GetString("TrackOutputArray", "tracks"));
|
---|
| 119 | }
|
---|
| 120 |
|
---|
| 121 | //------------------------------------------------------------------------------
|
---|
| 122 |
|
---|
| 123 | void DenseTrackFilter::Finish()
|
---|
| 124 | {
|
---|
| 125 | vector< vector< Double_t >* >::iterator itPhiBin;
|
---|
| 126 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
| 127 | for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
|
---|
| 128 | {
|
---|
| 129 | delete *itPhiBin;
|
---|
| 130 | }
|
---|
| 131 | }
|
---|
| 132 |
|
---|
| 133 | //------------------------------------------------------------------------------
|
---|
| 134 |
|
---|
| 135 | void DenseTrackFilter::Process()
|
---|
| 136 | {
|
---|
| 137 | Candidate *candidate, *track, *bestTrack, *mother;
|
---|
| 138 | TLorentzVector position, momentum;
|
---|
| 139 | Short_t etaBin, phiBin, flags;
|
---|
[f8a2939] | 140 | Int_t number, towerTrackHits;
|
---|
[196994a] | 141 | Long64_t towerHit, towerEtaPhi, hitEtaPhi;
|
---|
| 142 |
|
---|
| 143 | Double_t pt, ptmax, eta, phi;
|
---|
| 144 |
|
---|
| 145 | vector< Double_t >::iterator itEtaBin;
|
---|
| 146 | vector< Double_t >::iterator itPhiBin;
|
---|
| 147 | vector< Double_t > *phiBins;
|
---|
| 148 |
|
---|
| 149 | vector< Long64_t >::iterator itTowerHits;
|
---|
| 150 |
|
---|
| 151 | fTowerHits.clear();
|
---|
| 152 |
|
---|
| 153 | // loop over all tracks
|
---|
| 154 | fItTrackInputArray->Reset();
|
---|
| 155 | number = -1;
|
---|
| 156 | while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
| 157 | {
|
---|
| 158 | const TLorentzVector &trackPosition = track->Position;
|
---|
| 159 | ++number;
|
---|
| 160 |
|
---|
| 161 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 162 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
|
---|
| 163 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 164 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 165 |
|
---|
| 166 | // phi bins for given eta bin
|
---|
| 167 | phiBins = fPhiBins[etaBin];
|
---|
| 168 |
|
---|
| 169 | // find phi bin [1, phiBins.size - 1]
|
---|
| 170 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
|
---|
| 171 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 172 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 173 |
|
---|
| 174 | flags = 1;
|
---|
| 175 |
|
---|
| 176 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
|
---|
| 177 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
| 178 |
|
---|
| 179 | fTowerHits.push_back(towerHit);
|
---|
| 180 | }
|
---|
| 181 |
|
---|
| 182 | // all hits are sorted first by eta bin number, then by phi bin number,
|
---|
| 183 | // then by flags and then by particle or track number
|
---|
| 184 | sort(fTowerHits.begin(), fTowerHits.end());
|
---|
| 185 |
|
---|
| 186 | // loop over all hits
|
---|
| 187 | towerEtaPhi = 0;
|
---|
| 188 | bestTrack = 0;
|
---|
| 189 | fTower = 0;
|
---|
| 190 | ptmax = 0.0;
|
---|
[f8a2939] | 191 | towerTrackHits = 0;
|
---|
[196994a] | 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 | if(towerEtaPhi != hitEtaPhi)
|
---|
| 200 | {
|
---|
| 201 | // switch to next tower
|
---|
| 202 | towerEtaPhi = hitEtaPhi;
|
---|
| 203 |
|
---|
| 204 | // saving track with highest pT that hit the tower
|
---|
[f8a2939] | 205 | if(towerTrackHits > 0)
|
---|
[196994a] | 206 | {
|
---|
| 207 | mother = bestTrack;
|
---|
| 208 | candidate = static_cast<Candidate*>(bestTrack->Clone());
|
---|
| 209 | pt = candidate->Momentum.Pt();
|
---|
| 210 | eta = candidate->Momentum.Eta();
|
---|
| 211 | phi = candidate->Momentum.Phi();
|
---|
| 212 | eta = gRandom->Gaus(eta, fEtaPhiRes);
|
---|
| 213 | phi = gRandom->Gaus(phi, fEtaPhiRes);
|
---|
| 214 | candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt*TMath::CosH(eta));
|
---|
| 215 | candidate->AddCandidate(mother);
|
---|
| 216 | fTrackOutputArray->Add(candidate);
|
---|
| 217 | }
|
---|
| 218 |
|
---|
| 219 | ptmax = 0.0;
|
---|
[f8a2939] | 220 | towerTrackHits = 0;
|
---|
[196994a] | 221 | bestTrack = 0;
|
---|
| 222 | }
|
---|
| 223 | // check for track hits
|
---|
| 224 |
|
---|
| 225 | if(flags & 1)
|
---|
| 226 | {
|
---|
[f8a2939] | 227 | ++towerTrackHits;
|
---|
[196994a] | 228 | track = static_cast<Candidate*>(fTrackInputArray->At(number));
|
---|
| 229 | momentum = track->Momentum;
|
---|
| 230 |
|
---|
| 231 | if (momentum.Pt() > ptmax)
|
---|
| 232 | {
|
---|
| 233 | ptmax = momentum.Pt();
|
---|
| 234 | bestTrack = track;
|
---|
| 235 | }
|
---|
| 236 | continue;
|
---|
| 237 | }
|
---|
| 238 | }
|
---|
| 239 | }
|
---|