[01f457a] | 1 | /*
|
---|
| 2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
| 3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
[1fa50c2] | 4 | *
|
---|
[01f457a] | 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.
|
---|
[1fa50c2] | 9 | *
|
---|
[01f457a] | 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.
|
---|
[1fa50c2] | 14 | *
|
---|
[01f457a] | 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 |
|
---|
[d7d2da3] | 19 |
|
---|
| 20 | /** \class Calorimeter
|
---|
| 21 | *
|
---|
| 22 | * Fills calorimeter towers, performs calorimeter resolution smearing,
|
---|
[d4c4d9d] | 23 | * and creates energy flow objects (tracks, photons, and neutral hadrons).
|
---|
[d7d2da3] | 24 | *
|
---|
| 25 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 26 | *
|
---|
| 27 | */
|
---|
| 28 |
|
---|
| 29 | #include "modules/Calorimeter.h"
|
---|
| 30 |
|
---|
| 31 | #include "classes/DelphesClasses.h"
|
---|
| 32 | #include "classes/DelphesFactory.h"
|
---|
| 33 | #include "classes/DelphesFormula.h"
|
---|
| 34 |
|
---|
| 35 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 36 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 37 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 38 |
|
---|
| 39 | #include "TMath.h"
|
---|
| 40 | #include "TString.h"
|
---|
| 41 | #include "TFormula.h"
|
---|
| 42 | #include "TRandom3.h"
|
---|
| 43 | #include "TObjArray.h"
|
---|
| 44 | #include "TDatabasePDG.h"
|
---|
| 45 | #include "TLorentzVector.h"
|
---|
| 46 |
|
---|
| 47 | #include <algorithm>
|
---|
| 48 | #include <stdexcept>
|
---|
| 49 | #include <iostream>
|
---|
| 50 | #include <sstream>
|
---|
| 51 |
|
---|
| 52 | using namespace std;
|
---|
| 53 |
|
---|
| 54 | //------------------------------------------------------------------------------
|
---|
| 55 |
|
---|
| 56 | Calorimeter::Calorimeter() :
|
---|
[82575a3] | 57 | fECalResolutionFormula(0), fHCalResolutionFormula(0),
|
---|
[d7d2da3] | 58 | fItParticleInputArray(0), fItTrackInputArray(0),
|
---|
[2dab783] | 59 | fTowerTrackArray(0), fItTowerTrackArray(0)
|
---|
[d7d2da3] | 60 | {
|
---|
[82575a3] | 61 | fECalResolutionFormula = new DelphesFormula;
|
---|
| 62 | fHCalResolutionFormula = new DelphesFormula;
|
---|
| 63 |
|
---|
[9493a0f] | 64 | fTowerTrackArray = new TObjArray;
|
---|
| 65 | fItTowerTrackArray = fTowerTrackArray->MakeIterator();
|
---|
[d7d2da3] | 66 | }
|
---|
| 67 |
|
---|
| 68 | //------------------------------------------------------------------------------
|
---|
| 69 |
|
---|
| 70 | Calorimeter::~Calorimeter()
|
---|
| 71 | {
|
---|
[82575a3] | 72 | if(fECalResolutionFormula) delete fECalResolutionFormula;
|
---|
| 73 | if(fHCalResolutionFormula) delete fHCalResolutionFormula;
|
---|
| 74 |
|
---|
[9493a0f] | 75 | if(fTowerTrackArray) delete fTowerTrackArray;
|
---|
| 76 | if(fItTowerTrackArray) delete fItTowerTrackArray;
|
---|
[2dab783] | 77 | }
|
---|
[d7d2da3] | 78 |
|
---|
| 79 | //------------------------------------------------------------------------------
|
---|
| 80 |
|
---|
| 81 | void Calorimeter::Init()
|
---|
| 82 | {
|
---|
| 83 | ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
|
---|
[a221d1f] | 84 | Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
|
---|
[82575a3] | 85 | Double_t ecalFraction, hcalFraction;
|
---|
[d7d2da3] | 86 | TBinMap::iterator itEtaBin;
|
---|
| 87 | set< Double_t >::iterator itPhiBin;
|
---|
| 88 | vector< Double_t > *phiBins;
|
---|
| 89 |
|
---|
| 90 | // read eta and phi bins
|
---|
| 91 | param = GetParam("EtaPhiBins");
|
---|
| 92 | size = param.GetSize();
|
---|
| 93 | fBinMap.clear();
|
---|
| 94 | fEtaBins.clear();
|
---|
| 95 | fPhiBins.clear();
|
---|
| 96 | for(i = 0; i < size/2; ++i)
|
---|
| 97 | {
|
---|
| 98 | paramEtaBins = param[i*2];
|
---|
| 99 | sizeEtaBins = paramEtaBins.GetSize();
|
---|
| 100 | paramPhiBins = param[i*2 + 1];
|
---|
| 101 | sizePhiBins = paramPhiBins.GetSize();
|
---|
| 102 |
|
---|
| 103 | for(j = 0; j < sizeEtaBins; ++j)
|
---|
| 104 | {
|
---|
| 105 | for(k = 0; k < sizePhiBins; ++k)
|
---|
| 106 | {
|
---|
| 107 | fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
|
---|
| 108 | }
|
---|
| 109 | }
|
---|
| 110 | }
|
---|
| 111 |
|
---|
| 112 | // for better performance we transform map of sets to parallel vectors:
|
---|
| 113 | // vector< double > and vector< vector< double >* >
|
---|
| 114 | for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
|
---|
| 115 | {
|
---|
| 116 | fEtaBins.push_back(itEtaBin->first);
|
---|
| 117 | phiBins = new vector< double >(itEtaBin->second.size());
|
---|
| 118 | fPhiBins.push_back(phiBins);
|
---|
| 119 | phiBins->clear();
|
---|
| 120 | for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
|
---|
| 121 | {
|
---|
| 122 | phiBins->push_back(*itPhiBin);
|
---|
| 123 | }
|
---|
| 124 | }
|
---|
| 125 |
|
---|
| 126 | // read energy fractions for different particles
|
---|
| 127 | param = GetParam("EnergyFraction");
|
---|
| 128 | size = param.GetSize();
|
---|
| 129 |
|
---|
| 130 | // set default energy fractions values
|
---|
| 131 | fFractionMap.clear();
|
---|
[82575a3] | 132 | fFractionMap[0] = make_pair(0.0, 1.0);
|
---|
[d7d2da3] | 133 |
|
---|
| 134 | for(i = 0; i < size/2; ++i)
|
---|
| 135 | {
|
---|
| 136 | paramFractions = param[i*2 + 1];
|
---|
[82575a3] | 137 |
|
---|
| 138 | ecalFraction = paramFractions[0].GetDouble();
|
---|
| 139 | hcalFraction = paramFractions[1].GetDouble();
|
---|
| 140 |
|
---|
| 141 | fFractionMap[param[i*2].GetInt()] = make_pair(ecalFraction, hcalFraction);
|
---|
[d7d2da3] | 142 | }
|
---|
[8624f58] | 143 |
|
---|
[3db5282] | 144 | // read min E value for timing measurement in ECAL
|
---|
[839deb7] | 145 | fTimingEnergyMin = GetDouble("TimingEnergyMin",4.);
|
---|
[3db5282] | 146 | // For timing
|
---|
| 147 | // So far this flag needs to be false
|
---|
| 148 | // Curved extrapolation not supported
|
---|
| 149 | fElectronsFromTrack = false;
|
---|
| 150 |
|
---|
[4b9a2dc] | 151 | // read min E value for towers to be saved
|
---|
[38bf1ae] | 152 | fECalEnergyMin = GetDouble("ECalEnergyMin", 0.0);
|
---|
| 153 | fHCalEnergyMin = GetDouble("HCalEnergyMin", 0.0);
|
---|
[a221d1f] | 154 |
|
---|
[38bf1ae] | 155 | fECalEnergySignificanceMin = GetDouble("ECalEnergySignificanceMin", 0.0);
|
---|
| 156 | fHCalEnergySignificanceMin = GetDouble("HCalEnergySignificanceMin", 0.0);
|
---|
[a221d1f] | 157 |
|
---|
[8624f58] | 158 | // switch on or off the dithering of the center of calorimeter towers
|
---|
[4e09c3a] | 159 | fSmearTowerCenter = GetBool("SmearTowerCenter", true);
|
---|
[4b9a2dc] | 160 |
|
---|
[d7d2da3] | 161 | // read resolution formulas
|
---|
[82575a3] | 162 | fECalResolutionFormula->Compile(GetString("ECalResolutionFormula", "0"));
|
---|
| 163 | fHCalResolutionFormula->Compile(GetString("HCalResolutionFormula", "0"));
|
---|
| 164 |
|
---|
[d7d2da3] | 165 | // import array with output from other modules
|
---|
| 166 | fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
|
---|
| 167 | fItParticleInputArray = fParticleInputArray->MakeIterator();
|
---|
| 168 |
|
---|
| 169 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
|
---|
| 170 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
| 171 |
|
---|
| 172 | // create output arrays
|
---|
| 173 | fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
|
---|
[82575a3] | 174 | fPhotonOutputArray = ExportArray(GetString("PhotonOutputArray", "photons"));
|
---|
[a221d1f] | 175 |
|
---|
[82575a3] | 176 | fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
|
---|
| 177 | fEFlowPhotonOutputArray = ExportArray(GetString("EFlowPhotonOutputArray", "eflowPhotons"));
|
---|
| 178 | fEFlowNeutralHadronOutputArray = ExportArray(GetString("EFlowNeutralHadronOutputArray", "eflowNeutralHadrons"));
|
---|
[d7d2da3] | 179 | }
|
---|
| 180 |
|
---|
| 181 | //------------------------------------------------------------------------------
|
---|
| 182 |
|
---|
| 183 | void Calorimeter::Finish()
|
---|
| 184 | {
|
---|
[2dab783] | 185 | vector< vector< Double_t >* >::iterator itPhiBin;
|
---|
[d7d2da3] | 186 | if(fItParticleInputArray) delete fItParticleInputArray;
|
---|
| 187 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
| 188 | for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
|
---|
| 189 | {
|
---|
| 190 | delete *itPhiBin;
|
---|
| 191 | }
|
---|
| 192 | }
|
---|
| 193 |
|
---|
| 194 | //------------------------------------------------------------------------------
|
---|
| 195 |
|
---|
| 196 | void Calorimeter::Process()
|
---|
| 197 | {
|
---|
| 198 | Candidate *particle, *track;
|
---|
| 199 | TLorentzVector position, momentum;
|
---|
| 200 | Short_t etaBin, phiBin, flags;
|
---|
| 201 | Int_t number;
|
---|
| 202 | Long64_t towerHit, towerEtaPhi, hitEtaPhi;
|
---|
[82575a3] | 203 | Double_t ecalFraction, hcalFraction;
|
---|
| 204 | Double_t ecalEnergy, hcalEnergy;
|
---|
[d7d2da3] | 205 | Int_t pdgCode;
|
---|
| 206 |
|
---|
| 207 | TFractionMap::iterator itFractionMap;
|
---|
| 208 |
|
---|
| 209 | vector< Double_t >::iterator itEtaBin;
|
---|
| 210 | vector< Double_t >::iterator itPhiBin;
|
---|
| 211 | vector< Double_t > *phiBins;
|
---|
| 212 |
|
---|
| 213 | vector< Long64_t >::iterator itTowerHits;
|
---|
| 214 |
|
---|
| 215 | DelphesFactory *factory = GetFactory();
|
---|
| 216 | fTowerHits.clear();
|
---|
[82575a3] | 217 | fTowerECalFractions.clear();
|
---|
| 218 | fTowerHCalFractions.clear();
|
---|
| 219 | fTrackECalFractions.clear();
|
---|
| 220 | fTrackHCalFractions.clear();
|
---|
| 221 |
|
---|
[d7d2da3] | 222 | // loop over all particles
|
---|
| 223 | fItParticleInputArray->Reset();
|
---|
| 224 | number = -1;
|
---|
| 225 | while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
|
---|
| 226 | {
|
---|
| 227 | const TLorentzVector &particlePosition = particle->Position;
|
---|
| 228 | ++number;
|
---|
| 229 |
|
---|
| 230 | pdgCode = TMath::Abs(particle->PID);
|
---|
| 231 |
|
---|
| 232 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 233 | if(itFractionMap == fFractionMap.end())
|
---|
| 234 | {
|
---|
| 235 | itFractionMap = fFractionMap.find(0);
|
---|
| 236 | }
|
---|
| 237 |
|
---|
[82575a3] | 238 | ecalFraction = itFractionMap->second.first;
|
---|
| 239 | hcalFraction = itFractionMap->second.second;
|
---|
| 240 |
|
---|
| 241 | fTowerECalFractions.push_back(ecalFraction);
|
---|
| 242 | fTowerHCalFractions.push_back(hcalFraction);
|
---|
| 243 |
|
---|
| 244 | if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
|
---|
[d7d2da3] | 245 |
|
---|
| 246 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 247 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
|
---|
| 248 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 249 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 250 |
|
---|
| 251 | // phi bins for given eta bin
|
---|
| 252 | phiBins = fPhiBins[etaBin];
|
---|
| 253 |
|
---|
| 254 | // find phi bin [1, phiBins.size - 1]
|
---|
| 255 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
|
---|
| 256 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 257 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 258 |
|
---|
[73e0386] | 259 | flags = 0;
|
---|
[2dab783] | 260 | flags |= (pdgCode == 11 || pdgCode == 22) << 1;
|
---|
[d7d2da3] | 261 |
|
---|
| 262 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
|
---|
| 263 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
| 264 |
|
---|
| 265 | fTowerHits.push_back(towerHit);
|
---|
| 266 | }
|
---|
| 267 |
|
---|
| 268 | // loop over all tracks
|
---|
| 269 | fItTrackInputArray->Reset();
|
---|
| 270 | number = -1;
|
---|
| 271 | while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
| 272 | {
|
---|
| 273 | const TLorentzVector &trackPosition = track->Position;
|
---|
| 274 | ++number;
|
---|
| 275 |
|
---|
[73e0386] | 276 | pdgCode = TMath::Abs(track->PID);
|
---|
| 277 |
|
---|
| 278 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 279 | if(itFractionMap == fFractionMap.end())
|
---|
| 280 | {
|
---|
| 281 | itFractionMap = fFractionMap.find(0);
|
---|
| 282 | }
|
---|
| 283 |
|
---|
[82575a3] | 284 | ecalFraction = itFractionMap->second.first;
|
---|
| 285 | hcalFraction = itFractionMap->second.second;
|
---|
| 286 |
|
---|
| 287 | fTrackECalFractions.push_back(ecalFraction);
|
---|
| 288 | fTrackHCalFractions.push_back(hcalFraction);
|
---|
| 289 |
|
---|
[d7d2da3] | 290 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 291 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
|
---|
| 292 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 293 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 294 |
|
---|
| 295 | // phi bins for given eta bin
|
---|
| 296 | phiBins = fPhiBins[etaBin];
|
---|
| 297 |
|
---|
| 298 | // find phi bin [1, phiBins.size - 1]
|
---|
| 299 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
|
---|
| 300 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 301 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 302 |
|
---|
[73e0386] | 303 | flags = 1;
|
---|
| 304 |
|
---|
[d7d2da3] | 305 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
|
---|
[73e0386] | 306 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
[d7d2da3] | 307 |
|
---|
| 308 | fTowerHits.push_back(towerHit);
|
---|
| 309 | }
|
---|
| 310 |
|
---|
| 311 | // all hits are sorted first by eta bin number, then by phi bin number,
|
---|
| 312 | // then by flags and then by particle or track number
|
---|
| 313 | sort(fTowerHits.begin(), fTowerHits.end());
|
---|
| 314 |
|
---|
| 315 | // loop over all hits
|
---|
| 316 | towerEtaPhi = 0;
|
---|
| 317 | fTower = 0;
|
---|
| 318 | for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
|
---|
| 319 | {
|
---|
| 320 | towerHit = (*itTowerHits);
|
---|
| 321 | flags = (towerHit >> 24) & 0x00000000000000FFLL;
|
---|
| 322 | number = (towerHit) & 0x0000000000FFFFFFLL;
|
---|
| 323 | hitEtaPhi = towerHit >> 32;
|
---|
| 324 |
|
---|
| 325 | if(towerEtaPhi != hitEtaPhi)
|
---|
| 326 | {
|
---|
| 327 | // switch to next tower
|
---|
| 328 | towerEtaPhi = hitEtaPhi;
|
---|
| 329 |
|
---|
| 330 | // finalize previous tower
|
---|
| 331 | FinalizeTower();
|
---|
| 332 |
|
---|
| 333 | // create new tower
|
---|
| 334 | fTower = factory->NewCandidate();
|
---|
| 335 |
|
---|
| 336 | phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
|
---|
| 337 | etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
|
---|
| 338 |
|
---|
| 339 | // phi bins for given eta bin
|
---|
| 340 | phiBins = fPhiBins[etaBin];
|
---|
| 341 |
|
---|
| 342 | // calculate eta and phi of the tower's center
|
---|
| 343 | fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
|
---|
| 344 | fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
|
---|
| 345 |
|
---|
| 346 | fTowerEdges[0] = fEtaBins[etaBin - 1];
|
---|
| 347 | fTowerEdges[1] = fEtaBins[etaBin];
|
---|
| 348 | fTowerEdges[2] = (*phiBins)[phiBin - 1];
|
---|
| 349 | fTowerEdges[3] = (*phiBins)[phiBin];
|
---|
| 350 |
|
---|
[82575a3] | 351 | fTowerECalEnergy = 0.0;
|
---|
| 352 | fTowerHCalEnergy = 0.0;
|
---|
| 353 |
|
---|
| 354 | fTrackECalEnergy = 0.0;
|
---|
| 355 | fTrackHCalEnergy = 0.0;
|
---|
| 356 |
|
---|
[2dab783] | 357 | fTowerTrackHits = 0;
|
---|
| 358 | fTowerPhotonHits = 0;
|
---|
[a221d1f] | 359 |
|
---|
[9493a0f] | 360 | fTowerTrackArray->Clear();
|
---|
[d7d2da3] | 361 | }
|
---|
| 362 |
|
---|
| 363 | // check for track hits
|
---|
[73e0386] | 364 | if(flags & 1)
|
---|
[d7d2da3] | 365 | {
|
---|
[2dab783] | 366 | ++fTowerTrackHits;
|
---|
[73e0386] | 367 |
|
---|
[2dab783] | 368 | track = static_cast<Candidate*>(fTrackInputArray->At(number));
|
---|
| 369 | momentum = track->Momentum;
|
---|
[22dc7fd] | 370 | position = track->Position;
|
---|
[82575a3] | 371 |
|
---|
| 372 | ecalEnergy = momentum.E() * fTrackECalFractions[number];
|
---|
| 373 | hcalEnergy = momentum.E() * fTrackHCalFractions[number];
|
---|
| 374 |
|
---|
| 375 | fTrackECalEnergy += ecalEnergy;
|
---|
| 376 | fTrackHCalEnergy += hcalEnergy;
|
---|
[a221d1f] | 377 |
|
---|
[839deb7] | 378 | if(ecalEnergy > fTimingEnergyMin && fTower)
|
---|
| 379 | {
|
---|
| 380 | if(fElectronsFromTrack)
|
---|
| 381 | {
|
---|
| 382 | fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, track->Position.T()));
|
---|
| 383 | }
|
---|
[3db5282] | 384 | }
|
---|
[a221d1f] | 385 |
|
---|
[2dab783] | 386 | fTowerTrackArray->Add(track);
|
---|
[73e0386] | 387 |
|
---|
[2dab783] | 388 | continue;
|
---|
[73e0386] | 389 | }
|
---|
[839deb7] | 390 |
|
---|
[73e0386] | 391 | // check for photon and electron hits in current tower
|
---|
[2dab783] | 392 | if(flags & 2) ++fTowerPhotonHits;
|
---|
[a221d1f] | 393 |
|
---|
[d7d2da3] | 394 | particle = static_cast<Candidate*>(fParticleInputArray->At(number));
|
---|
| 395 | momentum = particle->Momentum;
|
---|
[22dc7fd] | 396 | position = particle->Position;
|
---|
[d7d2da3] | 397 |
|
---|
| 398 | // fill current tower
|
---|
[82575a3] | 399 | ecalEnergy = momentum.E() * fTowerECalFractions[number];
|
---|
| 400 | hcalEnergy = momentum.E() * fTowerHCalFractions[number];
|
---|
| 401 |
|
---|
| 402 | fTowerECalEnergy += ecalEnergy;
|
---|
| 403 | fTowerHCalEnergy += hcalEnergy;
|
---|
| 404 |
|
---|
[839deb7] | 405 | if(ecalEnergy > fTimingEnergyMin && fTower)
|
---|
| 406 | {
|
---|
| 407 | if (abs(particle->PID) != 11 || !fElectronsFromTrack)
|
---|
| 408 | {
|
---|
| 409 | fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, particle->Position.T()));
|
---|
[3db5282] | 410 | }
|
---|
| 411 | }
|
---|
[82575a3] | 412 |
|
---|
[3079350] | 413 | fTower->AddCandidate(particle);
|
---|
[d7d2da3] | 414 | }
|
---|
| 415 |
|
---|
| 416 | // finalize last tower
|
---|
| 417 | FinalizeTower();
|
---|
| 418 | }
|
---|
| 419 |
|
---|
| 420 | //------------------------------------------------------------------------------
|
---|
| 421 |
|
---|
| 422 | void Calorimeter::FinalizeTower()
|
---|
| 423 | {
|
---|
[e2dd4c5] | 424 | Candidate *track, *tower, *mother;
|
---|
[d7d2da3] | 425 | Double_t energy, pt, eta, phi;
|
---|
[82575a3] | 426 | Double_t ecalEnergy, hcalEnergy;
|
---|
| 427 | Double_t ecalSigma, hcalSigma;
|
---|
[a98c7ef] | 428 | Double_t ecalTrkSigma, hcalTrkSigma;
|
---|
| 429 | Double_t ecalFraction, hcalFraction;
|
---|
[e2dd4c5] | 430 |
|
---|
[a98c7ef] | 431 | Int_t pdgCode;
|
---|
| 432 | TLorentzVector momentum;
|
---|
| 433 | TFractionMap::iterator itFractionMap;
|
---|
[e2dd4c5] | 434 |
|
---|
[839deb7] | 435 | Float_t weight, sumWeightedTime, sumWeight;
|
---|
| 436 |
|
---|
[d7d2da3] | 437 | if(!fTower) return;
|
---|
| 438 |
|
---|
[e2dd4c5] | 439 |
|
---|
[82575a3] | 440 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerECalEnergy);
|
---|
| 441 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerHCalEnergy);
|
---|
[4600a41] | 442 |
|
---|
[38bf1ae] | 443 | ecalEnergy = LogNormal(fTowerECalEnergy, ecalSigma);
|
---|
[82575a3] | 444 | hcalEnergy = LogNormal(fTowerHCalEnergy, hcalSigma);
|
---|
| 445 |
|
---|
[4b9a2dc] | 446 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
|
---|
| 447 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
|
---|
| 448 |
|
---|
[38bf1ae] | 449 | if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
|
---|
| 450 | if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
|
---|
[4b9a2dc] | 451 |
|
---|
[82575a3] | 452 | energy = ecalEnergy + hcalEnergy;
|
---|
[839deb7] | 453 |
|
---|
[4e09c3a] | 454 | if(fSmearTowerCenter)
|
---|
[a221d1f] | 455 | {
|
---|
[38bf1ae] | 456 | eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
|
---|
| 457 | phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
|
---|
[a221d1f] | 458 | }
|
---|
| 459 | else
|
---|
| 460 | {
|
---|
[38bf1ae] | 461 | eta = fTowerEta;
|
---|
| 462 | phi = fTowerPhi;
|
---|
[a221d1f] | 463 | }
|
---|
[d7d2da3] | 464 |
|
---|
| 465 | pt = energy / TMath::CosH(eta);
|
---|
| 466 |
|
---|
[3db5282] | 467 | // Time calculation for tower
|
---|
[839deb7] | 468 | fTower->NTimeHits = 0;
|
---|
| 469 | sumWeightedTime = 0.0;
|
---|
| 470 | sumWeight = 0.0;
|
---|
| 471 |
|
---|
| 472 | for(size_t i = 0; i < fTower->ECalEnergyTimePairs.size(); ++i)
|
---|
[3db5282] | 473 | {
|
---|
[839deb7] | 474 | weight = TMath::Sqrt(fTower->ECalEnergyTimePairs[i].first);
|
---|
| 475 | sumWeightedTime += weight * fTower->ECalEnergyTimePairs[i].second;
|
---|
| 476 | sumWeight += weight;
|
---|
| 477 | fTower->NTimeHits++;
|
---|
[3db5282] | 478 | }
|
---|
[839deb7] | 479 |
|
---|
| 480 | if(sumWeight > 0.0)
|
---|
| 481 | {
|
---|
| 482 | fTower->Position.SetPtEtaPhiE(1.0, eta, phi, sumWeightedTime/sumWeight);
|
---|
| 483 | }
|
---|
| 484 | else
|
---|
| 485 | {
|
---|
| 486 | fTower->Position.SetPtEtaPhiE(1.0, eta, phi, 999999.9);
|
---|
[3db5282] | 487 | }
|
---|
| 488 |
|
---|
| 489 |
|
---|
[d7d2da3] | 490 | fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
|
---|
[82575a3] | 491 | fTower->Eem = ecalEnergy;
|
---|
| 492 | fTower->Ehad = hcalEnergy;
|
---|
| 493 |
|
---|
[d7d2da3] | 494 | fTower->Edges[0] = fTowerEdges[0];
|
---|
| 495 | fTower->Edges[1] = fTowerEdges[1];
|
---|
| 496 | fTower->Edges[2] = fTowerEdges[2];
|
---|
| 497 | fTower->Edges[3] = fTowerEdges[3];
|
---|
| 498 |
|
---|
[38bf1ae] | 499 | if(energy > 0.0)
|
---|
[82575a3] | 500 | {
|
---|
| 501 | if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
|
---|
| 502 | {
|
---|
| 503 | fPhotonOutputArray->Add(fTower);
|
---|
| 504 | }
|
---|
[a221d1f] | 505 |
|
---|
[82575a3] | 506 | fTowerOutputArray->Add(fTower);
|
---|
| 507 | }
|
---|
[0d5f77c] | 508 |
|
---|
[f6b9fec] | 509 | // fill energy flow candidates
|
---|
[82575a3] | 510 |
|
---|
[e2dd4c5] | 511 | // save as eflowtracks only tracks that have better resolution than calo
|
---|
| 512 |
|
---|
| 513 | fItTowerTrackArray->Reset();
|
---|
| 514 | while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
|
---|
| 515 | {
|
---|
| 516 | pdgCode = TMath::Abs(track->PID);
|
---|
| 517 |
|
---|
| 518 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 519 | if(itFractionMap == fFractionMap.end())
|
---|
| 520 | {
|
---|
| 521 | itFractionMap = fFractionMap.find(0);
|
---|
| 522 | }
|
---|
| 523 |
|
---|
| 524 | ecalFraction = itFractionMap->second.first;
|
---|
| 525 | hcalFraction = itFractionMap->second.second;
|
---|
| 526 |
|
---|
| 527 | mother = track;
|
---|
| 528 | track = static_cast<Candidate*>(track->Clone());
|
---|
| 529 | track->AddCandidate(mother);
|
---|
| 530 | momentum = track->Momentum;
|
---|
| 531 |
|
---|
| 532 | if(ecalFraction > 1.0E-9 && hcalFraction < 1.0E-9 )
|
---|
| 533 | {
|
---|
| 534 | ecalTrkSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
|
---|
| 535 | if(ecalTrkSigma/momentum.E() < track->TrackResolution)
|
---|
| 536 | {
|
---|
| 537 | momentum *= ecalEnergy/fTrackECalEnergy;
|
---|
| 538 | }
|
---|
| 539 | }
|
---|
| 540 | else if(ecalFraction < 1.0E-9 && hcalFraction > 1.0E-9 )
|
---|
| 541 | {
|
---|
| 542 | hcalTrkSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
|
---|
| 543 | if(hcalTrkSigma/momentum.E() < track->TrackResolution)
|
---|
| 544 | {
|
---|
| 545 | momentum *= hcalEnergy/fTrackHCalEnergy;
|
---|
| 546 | }
|
---|
| 547 | }
|
---|
| 548 |
|
---|
| 549 | fEFlowTrackOutputArray->Add(track);
|
---|
[82575a3] | 550 | }
|
---|
| 551 |
|
---|
[e2dd4c5] | 552 | ecalEnergy -= fTrackECalEnergy;
|
---|
| 553 | hcalEnergy -= fTrackHCalEnergy;
|
---|
| 554 |
|
---|
[38bf1ae] | 555 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
|
---|
| 556 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
|
---|
| 557 |
|
---|
| 558 | if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
|
---|
| 559 | if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
|
---|
[82575a3] | 560 |
|
---|
| 561 | energy = ecalEnergy + hcalEnergy;
|
---|
| 562 |
|
---|
| 563 | if(ecalEnergy > 0.0)
|
---|
[2dab783] | 564 | {
|
---|
[d4c4d9d] | 565 | // create new photon tower
|
---|
[2dab783] | 566 | tower = static_cast<Candidate*>(fTower->Clone());
|
---|
[0d5f77c] | 567 |
|
---|
[82575a3] | 568 | pt = ecalEnergy / TMath::CosH(eta);
|
---|
| 569 |
|
---|
| 570 | tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalEnergy);
|
---|
| 571 | tower->Eem = ecalEnergy;
|
---|
[38bf1ae] | 572 | tower->Ehad = 0.0;
|
---|
[82575a3] | 573 |
|
---|
| 574 | fEFlowPhotonOutputArray->Add(tower);
|
---|
| 575 | }
|
---|
| 576 | if(hcalEnergy > 0.0)
|
---|
| 577 | {
|
---|
| 578 | // create new neutral hadron tower
|
---|
| 579 | tower = static_cast<Candidate*>(fTower->Clone());
|
---|
| 580 |
|
---|
| 581 | pt = hcalEnergy / TMath::CosH(eta);
|
---|
| 582 |
|
---|
| 583 | tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalEnergy);
|
---|
[38bf1ae] | 584 | tower->Eem = 0.0;
|
---|
[82575a3] | 585 | tower->Ehad = hcalEnergy;
|
---|
| 586 |
|
---|
| 587 | fEFlowNeutralHadronOutputArray->Add(tower);
|
---|
[d7d2da3] | 588 | }
|
---|
| 589 | }
|
---|
| 590 |
|
---|
| 591 | //------------------------------------------------------------------------------
|
---|
[39022b6] | 592 |
|
---|
[4600a41] | 593 | Double_t Calorimeter::LogNormal(Double_t mean, Double_t sigma)
|
---|
| 594 | {
|
---|
| 595 | Double_t a, b;
|
---|
| 596 |
|
---|
| 597 | if(mean > 0.0)
|
---|
[39022b6] | 598 | {
|
---|
[4600a41] | 599 | b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
|
---|
| 600 | a = TMath::Log(mean) - 0.5*b*b;
|
---|
| 601 |
|
---|
[38bf1ae] | 602 | return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
|
---|
[39022b6] | 603 | }
|
---|
[4600a41] | 604 | else
|
---|
| 605 | {
|
---|
| 606 | return 0.0;
|
---|
| 607 | }
|
---|
| 608 | }
|
---|