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