[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),
|
---|
[00e8dca] | 58 | fItParticleInputArray(0), fItTrackInputArray(0)
|
---|
[d7d2da3] | 59 | {
|
---|
[04290b1] | 60 |
|
---|
[82575a3] | 61 | fECalResolutionFormula = new DelphesFormula;
|
---|
| 62 | fHCalResolutionFormula = new DelphesFormula;
|
---|
| 63 |
|
---|
[04290b1] | 64 | fECalTowerTrackArray = new TObjArray;
|
---|
| 65 | fItECalTowerTrackArray = fECalTowerTrackArray->MakeIterator();
|
---|
[00e8dca] | 66 |
|
---|
[04290b1] | 67 | fHCalTowerTrackArray = new TObjArray;
|
---|
| 68 | fItHCalTowerTrackArray = fHCalTowerTrackArray->MakeIterator();
|
---|
| 69 |
|
---|
[d7d2da3] | 70 | }
|
---|
| 71 |
|
---|
| 72 | //------------------------------------------------------------------------------
|
---|
| 73 |
|
---|
| 74 | Calorimeter::~Calorimeter()
|
---|
| 75 | {
|
---|
[04290b1] | 76 |
|
---|
[82575a3] | 77 | if(fECalResolutionFormula) delete fECalResolutionFormula;
|
---|
| 78 | if(fHCalResolutionFormula) delete fHCalResolutionFormula;
|
---|
| 79 |
|
---|
[04290b1] | 80 | if(fECalTowerTrackArray) delete fECalTowerTrackArray;
|
---|
| 81 | if(fItECalTowerTrackArray) delete fItECalTowerTrackArray;
|
---|
[00e8dca] | 82 |
|
---|
[04290b1] | 83 | if(fHCalTowerTrackArray) delete fHCalTowerTrackArray;
|
---|
| 84 | if(fItHCalTowerTrackArray) delete fItHCalTowerTrackArray;
|
---|
| 85 |
|
---|
[2dab783] | 86 | }
|
---|
[d7d2da3] | 87 |
|
---|
| 88 | //------------------------------------------------------------------------------
|
---|
| 89 |
|
---|
| 90 | void Calorimeter::Init()
|
---|
| 91 | {
|
---|
| 92 | ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
|
---|
[a221d1f] | 93 | Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
|
---|
[82575a3] | 94 | Double_t ecalFraction, hcalFraction;
|
---|
[d7d2da3] | 95 | TBinMap::iterator itEtaBin;
|
---|
| 96 | set< Double_t >::iterator itPhiBin;
|
---|
| 97 | vector< Double_t > *phiBins;
|
---|
| 98 |
|
---|
| 99 | // read eta and phi bins
|
---|
| 100 | param = GetParam("EtaPhiBins");
|
---|
| 101 | size = param.GetSize();
|
---|
| 102 | fBinMap.clear();
|
---|
| 103 | fEtaBins.clear();
|
---|
| 104 | fPhiBins.clear();
|
---|
| 105 | for(i = 0; i < size/2; ++i)
|
---|
| 106 | {
|
---|
| 107 | paramEtaBins = param[i*2];
|
---|
| 108 | sizeEtaBins = paramEtaBins.GetSize();
|
---|
| 109 | paramPhiBins = param[i*2 + 1];
|
---|
| 110 | sizePhiBins = paramPhiBins.GetSize();
|
---|
| 111 |
|
---|
| 112 | for(j = 0; j < sizeEtaBins; ++j)
|
---|
| 113 | {
|
---|
| 114 | for(k = 0; k < sizePhiBins; ++k)
|
---|
| 115 | {
|
---|
| 116 | fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
|
---|
| 117 | }
|
---|
| 118 | }
|
---|
| 119 | }
|
---|
| 120 |
|
---|
| 121 | // for better performance we transform map of sets to parallel vectors:
|
---|
| 122 | // vector< double > and vector< vector< double >* >
|
---|
| 123 | for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
|
---|
| 124 | {
|
---|
| 125 | fEtaBins.push_back(itEtaBin->first);
|
---|
| 126 | phiBins = new vector< double >(itEtaBin->second.size());
|
---|
| 127 | fPhiBins.push_back(phiBins);
|
---|
| 128 | phiBins->clear();
|
---|
| 129 | for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
|
---|
| 130 | {
|
---|
| 131 | phiBins->push_back(*itPhiBin);
|
---|
| 132 | }
|
---|
| 133 | }
|
---|
| 134 |
|
---|
| 135 | // read energy fractions for different particles
|
---|
| 136 | param = GetParam("EnergyFraction");
|
---|
| 137 | size = param.GetSize();
|
---|
| 138 |
|
---|
| 139 | // set default energy fractions values
|
---|
| 140 | fFractionMap.clear();
|
---|
[82575a3] | 141 | fFractionMap[0] = make_pair(0.0, 1.0);
|
---|
[d7d2da3] | 142 |
|
---|
| 143 | for(i = 0; i < size/2; ++i)
|
---|
| 144 | {
|
---|
| 145 | paramFractions = param[i*2 + 1];
|
---|
[82575a3] | 146 |
|
---|
| 147 | ecalFraction = paramFractions[0].GetDouble();
|
---|
| 148 | hcalFraction = paramFractions[1].GetDouble();
|
---|
| 149 |
|
---|
| 150 | fFractionMap[param[i*2].GetInt()] = make_pair(ecalFraction, hcalFraction);
|
---|
[d7d2da3] | 151 | }
|
---|
[8624f58] | 152 |
|
---|
[3db5282] | 153 | // read min E value for timing measurement in ECAL
|
---|
[839deb7] | 154 | fTimingEnergyMin = GetDouble("TimingEnergyMin",4.);
|
---|
[3db5282] | 155 | // For timing
|
---|
| 156 | // So far this flag needs to be false
|
---|
| 157 | // Curved extrapolation not supported
|
---|
| 158 | fElectronsFromTrack = false;
|
---|
| 159 |
|
---|
[4b9a2dc] | 160 | // read min E value for towers to be saved
|
---|
[38bf1ae] | 161 | fECalEnergyMin = GetDouble("ECalEnergyMin", 0.0);
|
---|
| 162 | fHCalEnergyMin = GetDouble("HCalEnergyMin", 0.0);
|
---|
[a221d1f] | 163 |
|
---|
[38bf1ae] | 164 | fECalEnergySignificanceMin = GetDouble("ECalEnergySignificanceMin", 0.0);
|
---|
| 165 | fHCalEnergySignificanceMin = GetDouble("HCalEnergySignificanceMin", 0.0);
|
---|
[a221d1f] | 166 |
|
---|
[8624f58] | 167 | // switch on or off the dithering of the center of calorimeter towers
|
---|
[4e09c3a] | 168 | fSmearTowerCenter = GetBool("SmearTowerCenter", true);
|
---|
[4b9a2dc] | 169 |
|
---|
[d7d2da3] | 170 | // read resolution formulas
|
---|
[82575a3] | 171 | fECalResolutionFormula->Compile(GetString("ECalResolutionFormula", "0"));
|
---|
| 172 | fHCalResolutionFormula->Compile(GetString("HCalResolutionFormula", "0"));
|
---|
| 173 |
|
---|
[d7d2da3] | 174 | // import array with output from other modules
|
---|
| 175 | fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
|
---|
| 176 | fItParticleInputArray = fParticleInputArray->MakeIterator();
|
---|
| 177 |
|
---|
| 178 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
|
---|
| 179 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
| 180 |
|
---|
| 181 | // create output arrays
|
---|
| 182 | fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
|
---|
[82575a3] | 183 | fPhotonOutputArray = ExportArray(GetString("PhotonOutputArray", "photons"));
|
---|
[a221d1f] | 184 |
|
---|
[82575a3] | 185 | fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
|
---|
| 186 | fEFlowPhotonOutputArray = ExportArray(GetString("EFlowPhotonOutputArray", "eflowPhotons"));
|
---|
| 187 | fEFlowNeutralHadronOutputArray = ExportArray(GetString("EFlowNeutralHadronOutputArray", "eflowNeutralHadrons"));
|
---|
[d7d2da3] | 188 | }
|
---|
| 189 |
|
---|
| 190 | //------------------------------------------------------------------------------
|
---|
| 191 |
|
---|
| 192 | void Calorimeter::Finish()
|
---|
| 193 | {
|
---|
[2dab783] | 194 | vector< vector< Double_t >* >::iterator itPhiBin;
|
---|
[d7d2da3] | 195 | if(fItParticleInputArray) delete fItParticleInputArray;
|
---|
| 196 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
| 197 | for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
|
---|
| 198 | {
|
---|
| 199 | delete *itPhiBin;
|
---|
| 200 | }
|
---|
| 201 | }
|
---|
| 202 |
|
---|
| 203 | //------------------------------------------------------------------------------
|
---|
| 204 |
|
---|
| 205 | void Calorimeter::Process()
|
---|
| 206 | {
|
---|
| 207 | Candidate *particle, *track;
|
---|
| 208 | TLorentzVector position, momentum;
|
---|
| 209 | Short_t etaBin, phiBin, flags;
|
---|
| 210 | Int_t number;
|
---|
| 211 | Long64_t towerHit, towerEtaPhi, hitEtaPhi;
|
---|
[82575a3] | 212 | Double_t ecalFraction, hcalFraction;
|
---|
| 213 | Double_t ecalEnergy, hcalEnergy;
|
---|
[298734e] | 214 | Double_t ecalSigma, hcalSigma;
|
---|
| 215 | Double_t energyGuess;
|
---|
[d7d2da3] | 216 | Int_t pdgCode;
|
---|
| 217 |
|
---|
| 218 | TFractionMap::iterator itFractionMap;
|
---|
| 219 |
|
---|
| 220 | vector< Double_t >::iterator itEtaBin;
|
---|
| 221 | vector< Double_t >::iterator itPhiBin;
|
---|
| 222 | vector< Double_t > *phiBins;
|
---|
| 223 |
|
---|
| 224 | vector< Long64_t >::iterator itTowerHits;
|
---|
| 225 |
|
---|
| 226 | DelphesFactory *factory = GetFactory();
|
---|
| 227 | fTowerHits.clear();
|
---|
[00e8dca] | 228 | fECalTowerFractions.clear();
|
---|
| 229 | fHCalTowerFractions.clear();
|
---|
| 230 | fECalTrackFractions.clear();
|
---|
| 231 | fHCalTrackFractions.clear();
|
---|
[82575a3] | 232 |
|
---|
[d7d2da3] | 233 | // loop over all particles
|
---|
| 234 | fItParticleInputArray->Reset();
|
---|
| 235 | number = -1;
|
---|
| 236 | while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
|
---|
| 237 | {
|
---|
| 238 | const TLorentzVector &particlePosition = particle->Position;
|
---|
| 239 | ++number;
|
---|
| 240 |
|
---|
| 241 | pdgCode = TMath::Abs(particle->PID);
|
---|
| 242 |
|
---|
| 243 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 244 | if(itFractionMap == fFractionMap.end())
|
---|
| 245 | {
|
---|
| 246 | itFractionMap = fFractionMap.find(0);
|
---|
| 247 | }
|
---|
| 248 |
|
---|
[82575a3] | 249 | ecalFraction = itFractionMap->second.first;
|
---|
| 250 | hcalFraction = itFractionMap->second.second;
|
---|
| 251 |
|
---|
[00e8dca] | 252 | fECalTowerFractions.push_back(ecalFraction);
|
---|
| 253 | fHCalTowerFractions.push_back(hcalFraction);
|
---|
[82575a3] | 254 |
|
---|
| 255 | if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
|
---|
[d7d2da3] | 256 |
|
---|
| 257 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 258 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
|
---|
| 259 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 260 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 261 |
|
---|
| 262 | // phi bins for given eta bin
|
---|
| 263 | phiBins = fPhiBins[etaBin];
|
---|
| 264 |
|
---|
| 265 | // find phi bin [1, phiBins.size - 1]
|
---|
| 266 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
|
---|
| 267 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 268 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 269 |
|
---|
[73e0386] | 270 | flags = 0;
|
---|
[2dab783] | 271 | flags |= (pdgCode == 11 || pdgCode == 22) << 1;
|
---|
[d7d2da3] | 272 |
|
---|
| 273 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
|
---|
| 274 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
| 275 |
|
---|
| 276 | fTowerHits.push_back(towerHit);
|
---|
| 277 | }
|
---|
| 278 |
|
---|
| 279 | // loop over all tracks
|
---|
| 280 | fItTrackInputArray->Reset();
|
---|
| 281 | number = -1;
|
---|
| 282 | while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
| 283 | {
|
---|
| 284 | const TLorentzVector &trackPosition = track->Position;
|
---|
| 285 | ++number;
|
---|
| 286 |
|
---|
[73e0386] | 287 | pdgCode = TMath::Abs(track->PID);
|
---|
| 288 |
|
---|
| 289 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 290 | if(itFractionMap == fFractionMap.end())
|
---|
| 291 | {
|
---|
| 292 | itFractionMap = fFractionMap.find(0);
|
---|
| 293 | }
|
---|
| 294 |
|
---|
[82575a3] | 295 | ecalFraction = itFractionMap->second.first;
|
---|
| 296 | hcalFraction = itFractionMap->second.second;
|
---|
| 297 |
|
---|
[00e8dca] | 298 | fECalTrackFractions.push_back(ecalFraction);
|
---|
| 299 | fHCalTrackFractions.push_back(hcalFraction);
|
---|
[82575a3] | 300 |
|
---|
[d7d2da3] | 301 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 302 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
|
---|
| 303 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 304 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 305 |
|
---|
| 306 | // phi bins for given eta bin
|
---|
| 307 | phiBins = fPhiBins[etaBin];
|
---|
| 308 |
|
---|
| 309 | // find phi bin [1, phiBins.size - 1]
|
---|
| 310 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
|
---|
| 311 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 312 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 313 |
|
---|
[73e0386] | 314 | flags = 1;
|
---|
| 315 |
|
---|
[d7d2da3] | 316 | // 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] | 317 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
[d7d2da3] | 318 |
|
---|
| 319 | fTowerHits.push_back(towerHit);
|
---|
| 320 | }
|
---|
| 321 |
|
---|
| 322 | // all hits are sorted first by eta bin number, then by phi bin number,
|
---|
| 323 | // then by flags and then by particle or track number
|
---|
| 324 | sort(fTowerHits.begin(), fTowerHits.end());
|
---|
| 325 |
|
---|
| 326 | // loop over all hits
|
---|
| 327 | towerEtaPhi = 0;
|
---|
| 328 | fTower = 0;
|
---|
| 329 | for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
|
---|
| 330 | {
|
---|
| 331 | towerHit = (*itTowerHits);
|
---|
| 332 | flags = (towerHit >> 24) & 0x00000000000000FFLL;
|
---|
| 333 | number = (towerHit) & 0x0000000000FFFFFFLL;
|
---|
| 334 | hitEtaPhi = towerHit >> 32;
|
---|
| 335 |
|
---|
| 336 | if(towerEtaPhi != hitEtaPhi)
|
---|
| 337 | {
|
---|
| 338 | // switch to next tower
|
---|
| 339 | towerEtaPhi = hitEtaPhi;
|
---|
| 340 |
|
---|
| 341 | // finalize previous tower
|
---|
| 342 | FinalizeTower();
|
---|
| 343 |
|
---|
| 344 | // create new tower
|
---|
| 345 | fTower = factory->NewCandidate();
|
---|
| 346 |
|
---|
| 347 | phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
|
---|
| 348 | etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
|
---|
| 349 |
|
---|
| 350 | // phi bins for given eta bin
|
---|
| 351 | phiBins = fPhiBins[etaBin];
|
---|
| 352 |
|
---|
| 353 | // calculate eta and phi of the tower's center
|
---|
| 354 | fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
|
---|
| 355 | fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
|
---|
| 356 |
|
---|
| 357 | fTowerEdges[0] = fEtaBins[etaBin - 1];
|
---|
| 358 | fTowerEdges[1] = fEtaBins[etaBin];
|
---|
| 359 | fTowerEdges[2] = (*phiBins)[phiBin - 1];
|
---|
| 360 | fTowerEdges[3] = (*phiBins)[phiBin];
|
---|
| 361 |
|
---|
[00e8dca] | 362 | fECalTowerEnergy = 0.0;
|
---|
| 363 | fHCalTowerEnergy = 0.0;
|
---|
| 364 |
|
---|
[04290b1] | 365 | fECalTrackEnergy = 0.0;
|
---|
| 366 | fHCalTrackEnergy = 0.0;
|
---|
| 367 |
|
---|
| 368 | fECalTrackSigma = 0.0;
|
---|
| 369 | fHCalTrackSigma = 0.0;
|
---|
| 370 |
|
---|
[2dab783] | 371 | fTowerTrackHits = 0;
|
---|
| 372 | fTowerPhotonHits = 0;
|
---|
[a221d1f] | 373 |
|
---|
[04290b1] | 374 | fECalTowerTrackArray->Clear();
|
---|
| 375 | fHCalTowerTrackArray->Clear();
|
---|
| 376 |
|
---|
[d7d2da3] | 377 | }
|
---|
| 378 |
|
---|
| 379 | // check for track hits
|
---|
[73e0386] | 380 | if(flags & 1)
|
---|
[d7d2da3] | 381 | {
|
---|
[2dab783] | 382 | ++fTowerTrackHits;
|
---|
[73e0386] | 383 |
|
---|
[2dab783] | 384 | track = static_cast<Candidate*>(fTrackInputArray->At(number));
|
---|
| 385 | momentum = track->Momentum;
|
---|
[22dc7fd] | 386 | position = track->Position;
|
---|
[82575a3] | 387 |
|
---|
[00e8dca] | 388 | ecalEnergy = momentum.E() * fECalTrackFractions[number];
|
---|
| 389 | hcalEnergy = momentum.E() * fHCalTrackFractions[number];
|
---|
[a221d1f] | 390 |
|
---|
[839deb7] | 391 | if(ecalEnergy > fTimingEnergyMin && fTower)
|
---|
| 392 | {
|
---|
| 393 | if(fElectronsFromTrack)
|
---|
| 394 | {
|
---|
| 395 | fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, track->Position.T()));
|
---|
| 396 | }
|
---|
[3db5282] | 397 | }
|
---|
[a221d1f] | 398 |
|
---|
[9da65a5] | 399 | if(fECalTrackFractions[number] > 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
|
---|
[00e8dca] | 400 | {
|
---|
[04290b1] | 401 | fECalTrackEnergy += ecalEnergy;
|
---|
[298734e] | 402 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
|
---|
| 403 | if(ecalSigma/momentum.E() < track->TrackResolution) energyGuess = ecalEnergy;
|
---|
| 404 | else energyGuess = momentum.E();
|
---|
| 405 |
|
---|
| 406 | fECalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
|
---|
[04290b1] | 407 | fECalTowerTrackArray->Add(track);
|
---|
[00e8dca] | 408 | }
|
---|
[04290b1] | 409 |
|
---|
[9da65a5] | 410 | else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] > 1.0E-9)
|
---|
[00e8dca] | 411 | {
|
---|
[04290b1] | 412 | fHCalTrackEnergy += hcalEnergy;
|
---|
[298734e] | 413 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
|
---|
| 414 | if(hcalSigma/momentum.E() < track->TrackResolution) energyGuess = hcalEnergy;
|
---|
| 415 | else energyGuess = momentum.E();
|
---|
| 416 |
|
---|
| 417 | fHCalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
|
---|
[04290b1] | 418 | fHCalTowerTrackArray->Add(track);
|
---|
[00e8dca] | 419 | }
|
---|
[04290b1] | 420 |
|
---|
[9da65a5] | 421 | else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
|
---|
| 422 | {
|
---|
| 423 | fEFlowTrackOutputArray->Add(track);
|
---|
| 424 | }
|
---|
[73e0386] | 425 |
|
---|
[2dab783] | 426 | continue;
|
---|
[73e0386] | 427 | }
|
---|
[839deb7] | 428 |
|
---|
[73e0386] | 429 | // check for photon and electron hits in current tower
|
---|
[2dab783] | 430 | if(flags & 2) ++fTowerPhotonHits;
|
---|
[a221d1f] | 431 |
|
---|
[d7d2da3] | 432 | particle = static_cast<Candidate*>(fParticleInputArray->At(number));
|
---|
| 433 | momentum = particle->Momentum;
|
---|
[22dc7fd] | 434 | position = particle->Position;
|
---|
[d7d2da3] | 435 |
|
---|
| 436 | // fill current tower
|
---|
[00e8dca] | 437 | ecalEnergy = momentum.E() * fECalTowerFractions[number];
|
---|
| 438 | hcalEnergy = momentum.E() * fHCalTowerFractions[number];
|
---|
[82575a3] | 439 |
|
---|
[00e8dca] | 440 | fECalTowerEnergy += ecalEnergy;
|
---|
| 441 | fHCalTowerEnergy += hcalEnergy;
|
---|
[82575a3] | 442 |
|
---|
[839deb7] | 443 | if(ecalEnergy > fTimingEnergyMin && fTower)
|
---|
| 444 | {
|
---|
| 445 | if (abs(particle->PID) != 11 || !fElectronsFromTrack)
|
---|
| 446 | {
|
---|
| 447 | fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, particle->Position.T()));
|
---|
[3db5282] | 448 | }
|
---|
| 449 | }
|
---|
[82575a3] | 450 |
|
---|
[3079350] | 451 | fTower->AddCandidate(particle);
|
---|
[d7d2da3] | 452 | }
|
---|
| 453 |
|
---|
| 454 | // finalize last tower
|
---|
| 455 | FinalizeTower();
|
---|
| 456 | }
|
---|
| 457 |
|
---|
| 458 | //------------------------------------------------------------------------------
|
---|
| 459 |
|
---|
| 460 | void Calorimeter::FinalizeTower()
|
---|
| 461 | {
|
---|
[e2dd4c5] | 462 | Candidate *track, *tower, *mother;
|
---|
[d7d2da3] | 463 | Double_t energy, pt, eta, phi;
|
---|
[82575a3] | 464 | Double_t ecalEnergy, hcalEnergy;
|
---|
[04290b1] | 465 | Double_t ecalNeutralEnergy, hcalNeutralEnergy;
|
---|
| 466 |
|
---|
[82575a3] | 467 | Double_t ecalSigma, hcalSigma;
|
---|
[04290b1] | 468 | Double_t ecalNeutralSigma, hcalNeutralSigma;
|
---|
[e2dd4c5] | 469 |
|
---|
[04290b1] | 470 | Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
|
---|
| 471 |
|
---|
[a98c7ef] | 472 | TLorentzVector momentum;
|
---|
| 473 | TFractionMap::iterator itFractionMap;
|
---|
[e2dd4c5] | 474 |
|
---|
[839deb7] | 475 | Float_t weight, sumWeightedTime, sumWeight;
|
---|
| 476 |
|
---|
[d7d2da3] | 477 | if(!fTower) return;
|
---|
| 478 |
|
---|
[00e8dca] | 479 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fECalTowerEnergy);
|
---|
| 480 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fHCalTowerEnergy);
|
---|
[4600a41] | 481 |
|
---|
[00e8dca] | 482 | ecalEnergy = LogNormal(fECalTowerEnergy, ecalSigma);
|
---|
| 483 | hcalEnergy = LogNormal(fHCalTowerEnergy, hcalSigma);
|
---|
[82575a3] | 484 |
|
---|
[4b9a2dc] | 485 | ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
|
---|
| 486 | hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
|
---|
| 487 |
|
---|
[38bf1ae] | 488 | if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
|
---|
| 489 | if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
|
---|
[4b9a2dc] | 490 |
|
---|
[82575a3] | 491 | energy = ecalEnergy + hcalEnergy;
|
---|
[839deb7] | 492 |
|
---|
[4e09c3a] | 493 | if(fSmearTowerCenter)
|
---|
[a221d1f] | 494 | {
|
---|
[38bf1ae] | 495 | eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
|
---|
| 496 | phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
|
---|
[a221d1f] | 497 | }
|
---|
| 498 | else
|
---|
| 499 | {
|
---|
[38bf1ae] | 500 | eta = fTowerEta;
|
---|
| 501 | phi = fTowerPhi;
|
---|
[a221d1f] | 502 | }
|
---|
[d7d2da3] | 503 |
|
---|
| 504 | pt = energy / TMath::CosH(eta);
|
---|
| 505 |
|
---|
[3db5282] | 506 | // Time calculation for tower
|
---|
[839deb7] | 507 | fTower->NTimeHits = 0;
|
---|
| 508 | sumWeightedTime = 0.0;
|
---|
| 509 | sumWeight = 0.0;
|
---|
| 510 |
|
---|
| 511 | for(size_t i = 0; i < fTower->ECalEnergyTimePairs.size(); ++i)
|
---|
[3db5282] | 512 | {
|
---|
[839deb7] | 513 | weight = TMath::Sqrt(fTower->ECalEnergyTimePairs[i].first);
|
---|
| 514 | sumWeightedTime += weight * fTower->ECalEnergyTimePairs[i].second;
|
---|
| 515 | sumWeight += weight;
|
---|
| 516 | fTower->NTimeHits++;
|
---|
[3db5282] | 517 | }
|
---|
[839deb7] | 518 |
|
---|
| 519 | if(sumWeight > 0.0)
|
---|
| 520 | {
|
---|
| 521 | fTower->Position.SetPtEtaPhiE(1.0, eta, phi, sumWeightedTime/sumWeight);
|
---|
| 522 | }
|
---|
| 523 | else
|
---|
| 524 | {
|
---|
| 525 | fTower->Position.SetPtEtaPhiE(1.0, eta, phi, 999999.9);
|
---|
[3db5282] | 526 | }
|
---|
| 527 |
|
---|
| 528 |
|
---|
[d7d2da3] | 529 | fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
|
---|
[82575a3] | 530 | fTower->Eem = ecalEnergy;
|
---|
| 531 | fTower->Ehad = hcalEnergy;
|
---|
| 532 |
|
---|
[d7d2da3] | 533 | fTower->Edges[0] = fTowerEdges[0];
|
---|
| 534 | fTower->Edges[1] = fTowerEdges[1];
|
---|
| 535 | fTower->Edges[2] = fTowerEdges[2];
|
---|
| 536 | fTower->Edges[3] = fTowerEdges[3];
|
---|
| 537 |
|
---|
[38bf1ae] | 538 | if(energy > 0.0)
|
---|
[82575a3] | 539 | {
|
---|
| 540 | if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
|
---|
| 541 | {
|
---|
| 542 | fPhotonOutputArray->Add(fTower);
|
---|
| 543 | }
|
---|
[a221d1f] | 544 |
|
---|
[82575a3] | 545 | fTowerOutputArray->Add(fTower);
|
---|
| 546 | }
|
---|
[04290b1] | 547 |
|
---|
[f6b9fec] | 548 | // fill energy flow candidates
|
---|
[04290b1] | 549 | fECalTrackSigma = TMath::Sqrt(fECalTrackSigma);
|
---|
| 550 | fHCalTrackSigma = TMath::Sqrt(fHCalTrackSigma);
|
---|
| 551 |
|
---|
| 552 | //compute neutral excesses
|
---|
| 553 | ecalNeutralEnergy = max( (ecalEnergy - fECalTrackEnergy) , 0.0);
|
---|
| 554 | hcalNeutralEnergy = max( (hcalEnergy - fHCalTrackEnergy) , 0.0);
|
---|
| 555 |
|
---|
| 556 | ecalNeutralSigma = ecalNeutralEnergy / TMath::Sqrt(fECalTrackSigma*fECalTrackSigma + ecalSigma*ecalSigma);
|
---|
| 557 | hcalNeutralSigma = hcalNeutralEnergy / TMath::Sqrt(fHCalTrackSigma*fHCalTrackSigma + hcalSigma*hcalSigma);
|
---|
[b6e6d36] | 558 |
|
---|
[04290b1] | 559 | // if ecal neutral excess is significant, simply create neutral EflowPhoton tower and clone each track into eflowtrack
|
---|
| 560 | if(ecalNeutralEnergy > fECalEnergyMin && ecalNeutralSigma > fECalEnergySignificanceMin)
|
---|
[00e8dca] | 561 | {
|
---|
[04290b1] | 562 | // create new photon tower
|
---|
| 563 | tower = static_cast<Candidate*>(fTower->Clone());
|
---|
| 564 | pt = ecalNeutralEnergy / TMath::CosH(eta);
|
---|
| 565 |
|
---|
| 566 | tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalNeutralEnergy);
|
---|
| 567 | tower->Eem = ecalNeutralEnergy;
|
---|
| 568 | tower->Ehad = 0.0;
|
---|
[efac6f9] | 569 | tower->PID = 22;
|
---|
[04290b1] | 570 |
|
---|
| 571 | fEFlowPhotonOutputArray->Add(tower);
|
---|
| 572 |
|
---|
| 573 | //clone tracks
|
---|
| 574 | fItECalTowerTrackArray->Reset();
|
---|
| 575 | while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
|
---|
| 576 | {
|
---|
| 577 | mother = track;
|
---|
| 578 | track = static_cast<Candidate*>(track->Clone());
|
---|
| 579 | track->AddCandidate(mother);
|
---|
[e2dd4c5] | 580 |
|
---|
[04290b1] | 581 | fEFlowTrackOutputArray->Add(track);
|
---|
| 582 | }
|
---|
| 583 |
|
---|
[82575a3] | 584 | }
|
---|
[04290b1] | 585 |
|
---|
| 586 | // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
|
---|
| 587 | else if(fECalTrackEnergy > 0.0)
|
---|
[00e8dca] | 588 | {
|
---|
[04290b1] | 589 | weightTrack = (fECalTrackSigma > 0.0) ? 1 / (fECalTrackSigma*fECalTrackSigma) : 0.0;
|
---|
| 590 | weightCalo = (ecalSigma > 0.0) ? 1 / (ecalSigma*ecalSigma) : 0.0;
|
---|
| 591 |
|
---|
| 592 | bestEnergyEstimate = (weightTrack*fECalTrackEnergy + weightCalo*ecalEnergy) / (weightTrack + weightCalo);
|
---|
| 593 | rescaleFactor = bestEnergyEstimate/fECalTrackEnergy;
|
---|
| 594 |
|
---|
| 595 | //rescale tracks
|
---|
| 596 | fItECalTowerTrackArray->Reset();
|
---|
| 597 | while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
|
---|
| 598 | {
|
---|
| 599 | mother = track;
|
---|
| 600 | track = static_cast<Candidate*>(track->Clone());
|
---|
| 601 | track->AddCandidate(mother);
|
---|
| 602 |
|
---|
| 603 | track->Momentum *= rescaleFactor;
|
---|
| 604 |
|
---|
| 605 | fEFlowTrackOutputArray->Add(track);
|
---|
| 606 | }
|
---|
[00e8dca] | 607 | }
|
---|
[e2dd4c5] | 608 |
|
---|
[82575a3] | 609 |
|
---|
[04290b1] | 610 | // if hcal neutral excess is significant, simply create neutral EflowNeutralHadron tower and clone each track into eflowtrack
|
---|
| 611 | if(hcalNeutralEnergy > fHCalEnergyMin && hcalNeutralSigma > fHCalEnergySignificanceMin)
|
---|
[2dab783] | 612 | {
|
---|
[d4c4d9d] | 613 | // create new photon tower
|
---|
[2dab783] | 614 | tower = static_cast<Candidate*>(fTower->Clone());
|
---|
[04290b1] | 615 | pt = hcalNeutralEnergy / TMath::CosH(eta);
|
---|
| 616 |
|
---|
| 617 | tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalNeutralEnergy);
|
---|
| 618 | tower->Ehad = hcalNeutralEnergy;
|
---|
| 619 | tower->Eem = 0.0;
|
---|
| 620 |
|
---|
| 621 | fEFlowNeutralHadronOutputArray->Add(tower);
|
---|
| 622 |
|
---|
| 623 | //clone tracks
|
---|
| 624 | fItHCalTowerTrackArray->Reset();
|
---|
| 625 | while((track = static_cast<Candidate*>(fItHCalTowerTrackArray->Next())))
|
---|
| 626 | {
|
---|
| 627 | mother = track;
|
---|
| 628 | track = static_cast<Candidate*>(track->Clone());
|
---|
| 629 | track->AddCandidate(mother);
|
---|
[0d5f77c] | 630 |
|
---|
[04290b1] | 631 | fEFlowTrackOutputArray->Add(track);
|
---|
| 632 | }
|
---|
| 633 |
|
---|
[82575a3] | 634 | }
|
---|
[04290b1] | 635 |
|
---|
| 636 | // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
|
---|
| 637 | else if(fHCalTrackEnergy > 0.0)
|
---|
[82575a3] | 638 | {
|
---|
[04290b1] | 639 | weightTrack = (fHCalTrackSigma > 0.0) ? 1 / (fHCalTrackSigma*fHCalTrackSigma) : 0.0;
|
---|
| 640 | weightCalo = (hcalSigma > 0.0) ? 1 / (hcalSigma*hcalSigma) : 0.0;
|
---|
| 641 |
|
---|
| 642 | bestEnergyEstimate = (weightTrack*fHCalTrackEnergy + weightCalo*hcalEnergy) / (weightTrack + weightCalo);
|
---|
| 643 | rescaleFactor = bestEnergyEstimate / fHCalTrackEnergy;
|
---|
| 644 |
|
---|
| 645 | //rescale tracks
|
---|
[b6e6d36] | 646 | fItHCalTowerTrackArray->Reset();
|
---|
| 647 | while((track = static_cast<Candidate*>(fItHCalTowerTrackArray->Next())))
|
---|
[04290b1] | 648 | {
|
---|
| 649 | mother = track;
|
---|
| 650 | track = static_cast<Candidate*>(track->Clone());
|
---|
| 651 | track->AddCandidate(mother);
|
---|
| 652 |
|
---|
| 653 | track->Momentum *= rescaleFactor;
|
---|
| 654 |
|
---|
| 655 | fEFlowTrackOutputArray->Add(track);
|
---|
| 656 | }
|
---|
[d7d2da3] | 657 | }
|
---|
[04290b1] | 658 |
|
---|
| 659 |
|
---|
[d7d2da3] | 660 | }
|
---|
| 661 |
|
---|
| 662 | //------------------------------------------------------------------------------
|
---|
[39022b6] | 663 |
|
---|
[4600a41] | 664 | Double_t Calorimeter::LogNormal(Double_t mean, Double_t sigma)
|
---|
| 665 | {
|
---|
| 666 | Double_t a, b;
|
---|
| 667 |
|
---|
| 668 | if(mean > 0.0)
|
---|
[39022b6] | 669 | {
|
---|
[4600a41] | 670 | b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
|
---|
| 671 | a = TMath::Log(mean) - 0.5*b*b;
|
---|
| 672 |
|
---|
[38bf1ae] | 673 | return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
|
---|
[39022b6] | 674 | }
|
---|
[4600a41] | 675 | else
|
---|
| 676 | {
|
---|
| 677 | return 0.0;
|
---|
| 678 | }
|
---|
| 679 | }
|
---|