[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 |
|
---|
[82575a3] | 19 |
|
---|
| 20 | /** \class SimpleCalorimeter
|
---|
| 21 | *
|
---|
| 22 | * Fills SimpleCalorimeter towers, performs SimpleCalorimeter resolution smearing,
|
---|
| 23 | * and creates energy flow objects (tracks, photons, and neutral hadrons).
|
---|
| 24 | *
|
---|
| 25 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 26 | *
|
---|
| 27 | */
|
---|
| 28 |
|
---|
| 29 | #include "modules/SimpleCalorimeter.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 | SimpleCalorimeter::SimpleCalorimeter() :
|
---|
| 57 | fResolutionFormula(0),
|
---|
| 58 | fItParticleInputArray(0), fItTrackInputArray(0),
|
---|
| 59 | fTowerTrackArray(0), fItTowerTrackArray(0)
|
---|
| 60 | {
|
---|
| 61 | fResolutionFormula = new DelphesFormula;
|
---|
[e5767b57] | 62 |
|
---|
[82575a3] | 63 | fTowerTrackArray = new TObjArray;
|
---|
| 64 | fItTowerTrackArray = fTowerTrackArray->MakeIterator();
|
---|
| 65 | }
|
---|
| 66 |
|
---|
| 67 | //------------------------------------------------------------------------------
|
---|
| 68 |
|
---|
| 69 | SimpleCalorimeter::~SimpleCalorimeter()
|
---|
| 70 | {
|
---|
| 71 | if(fResolutionFormula) delete fResolutionFormula;
|
---|
[e5767b57] | 72 |
|
---|
[82575a3] | 73 | if(fTowerTrackArray) delete fTowerTrackArray;
|
---|
| 74 | if(fItTowerTrackArray) delete fItTowerTrackArray;
|
---|
| 75 | }
|
---|
| 76 |
|
---|
| 77 | //------------------------------------------------------------------------------
|
---|
| 78 |
|
---|
| 79 | void SimpleCalorimeter::Init()
|
---|
| 80 | {
|
---|
| 81 | ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
|
---|
| 82 | Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
|
---|
| 83 | Double_t fraction;
|
---|
| 84 | TBinMap::iterator itEtaBin;
|
---|
| 85 | set< Double_t >::iterator itPhiBin;
|
---|
| 86 | vector< Double_t > *phiBins;
|
---|
| 87 |
|
---|
| 88 | // read eta and phi bins
|
---|
| 89 | param = GetParam("EtaPhiBins");
|
---|
| 90 | size = param.GetSize();
|
---|
| 91 | fBinMap.clear();
|
---|
| 92 | fEtaBins.clear();
|
---|
| 93 | fPhiBins.clear();
|
---|
| 94 | for(i = 0; i < size/2; ++i)
|
---|
| 95 | {
|
---|
| 96 | paramEtaBins = param[i*2];
|
---|
| 97 | sizeEtaBins = paramEtaBins.GetSize();
|
---|
| 98 | paramPhiBins = param[i*2 + 1];
|
---|
| 99 | sizePhiBins = paramPhiBins.GetSize();
|
---|
| 100 |
|
---|
| 101 | for(j = 0; j < sizeEtaBins; ++j)
|
---|
| 102 | {
|
---|
| 103 | for(k = 0; k < sizePhiBins; ++k)
|
---|
| 104 | {
|
---|
| 105 | fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
|
---|
| 106 | }
|
---|
| 107 | }
|
---|
| 108 | }
|
---|
| 109 |
|
---|
| 110 | // for better performance we transform map of sets to parallel vectors:
|
---|
| 111 | // vector< double > and vector< vector< double >* >
|
---|
| 112 | for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
|
---|
| 113 | {
|
---|
| 114 | fEtaBins.push_back(itEtaBin->first);
|
---|
| 115 | phiBins = new vector< double >(itEtaBin->second.size());
|
---|
| 116 | fPhiBins.push_back(phiBins);
|
---|
| 117 | phiBins->clear();
|
---|
| 118 | for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
|
---|
| 119 | {
|
---|
| 120 | phiBins->push_back(*itPhiBin);
|
---|
| 121 | }
|
---|
| 122 | }
|
---|
| 123 |
|
---|
| 124 | // read energy fractions for different particles
|
---|
| 125 | param = GetParam("EnergyFraction");
|
---|
| 126 | size = param.GetSize();
|
---|
| 127 |
|
---|
| 128 | // set default energy fractions values
|
---|
| 129 | fFractionMap.clear();
|
---|
| 130 | fFractionMap[0] = 1.0;
|
---|
| 131 |
|
---|
| 132 | for(i = 0; i < size/2; ++i)
|
---|
| 133 | {
|
---|
| 134 | paramFractions = param[i*2 + 1];
|
---|
| 135 | fraction = paramFractions[0].GetDouble();
|
---|
| 136 | fFractionMap[param[i*2].GetInt()] = fraction;
|
---|
| 137 | }
|
---|
[e5767b57] | 138 |
|
---|
[82575a3] | 139 | /*
|
---|
| 140 | TFractionMap::iterator itFractionMap;
|
---|
| 141 | for(itFractionMap = fFractionMap.begin(); itFractionMap != fFractionMap.end(); ++itFractionMap)
|
---|
| 142 | {
|
---|
| 143 | cout << itFractionMap->first << " " << itFractionMap->second.first << " " << itFractionMap->second.second << endl;
|
---|
| 144 | }
|
---|
| 145 | */
|
---|
[4b9a2dc] | 146 |
|
---|
| 147 | // read min E value for towers to be saved
|
---|
[e5767b57] | 148 | fEnergyMin = GetDouble("EnergyMin", 0.0);
|
---|
| 149 |
|
---|
| 150 | fEnergySignificanceMin = GetDouble("EnergySignificanceMin", 0.0);
|
---|
| 151 |
|
---|
| 152 | // switch on or off the dithering of the center of calorimeter towers
|
---|
| 153 | fDitherTowerCenter = GetBool("DitherTowerCenter", true);
|
---|
| 154 |
|
---|
[82575a3] | 155 | // read resolution formulas
|
---|
| 156 | fResolutionFormula->Compile(GetString("ResolutionFormula", "0"));
|
---|
[e5767b57] | 157 |
|
---|
[82575a3] | 158 | // import array with output from other modules
|
---|
| 159 | fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
|
---|
| 160 | fItParticleInputArray = fParticleInputArray->MakeIterator();
|
---|
| 161 |
|
---|
| 162 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
|
---|
| 163 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
| 164 |
|
---|
| 165 | // create output arrays
|
---|
| 166 | fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
|
---|
| 167 | fEFlowTowerOutputArray = ExportArray(GetString("EFlowTowerOutputArray", "eflowTowers"));
|
---|
| 168 | }
|
---|
| 169 |
|
---|
| 170 | //------------------------------------------------------------------------------
|
---|
| 171 |
|
---|
| 172 | void SimpleCalorimeter::Finish()
|
---|
| 173 | {
|
---|
| 174 | vector< vector< Double_t >* >::iterator itPhiBin;
|
---|
| 175 | if(fItParticleInputArray) delete fItParticleInputArray;
|
---|
| 176 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
| 177 | for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
|
---|
| 178 | {
|
---|
| 179 | delete *itPhiBin;
|
---|
| 180 | }
|
---|
| 181 | }
|
---|
| 182 |
|
---|
| 183 | //------------------------------------------------------------------------------
|
---|
| 184 |
|
---|
| 185 | void SimpleCalorimeter::Process()
|
---|
| 186 | {
|
---|
| 187 | Candidate *particle, *track;
|
---|
| 188 | TLorentzVector position, momentum;
|
---|
| 189 | Short_t etaBin, phiBin, flags;
|
---|
| 190 | Int_t number;
|
---|
| 191 | Long64_t towerHit, towerEtaPhi, hitEtaPhi;
|
---|
| 192 | Double_t fraction;
|
---|
| 193 | Double_t energy;
|
---|
| 194 | Int_t pdgCode;
|
---|
| 195 |
|
---|
| 196 | TFractionMap::iterator itFractionMap;
|
---|
| 197 |
|
---|
| 198 | vector< Double_t >::iterator itEtaBin;
|
---|
| 199 | vector< Double_t >::iterator itPhiBin;
|
---|
| 200 | vector< Double_t > *phiBins;
|
---|
| 201 |
|
---|
| 202 | vector< Long64_t >::iterator itTowerHits;
|
---|
| 203 |
|
---|
| 204 | DelphesFactory *factory = GetFactory();
|
---|
| 205 | fTowerHits.clear();
|
---|
| 206 | fTowerFractions.clear();
|
---|
| 207 | fTrackFractions.clear();
|
---|
[e5767b57] | 208 |
|
---|
[82575a3] | 209 | // loop over all particles
|
---|
| 210 | fItParticleInputArray->Reset();
|
---|
| 211 | number = -1;
|
---|
| 212 | while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
|
---|
| 213 | {
|
---|
| 214 | const TLorentzVector &particlePosition = particle->Position;
|
---|
| 215 | ++number;
|
---|
| 216 |
|
---|
| 217 | pdgCode = TMath::Abs(particle->PID);
|
---|
| 218 |
|
---|
| 219 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 220 | if(itFractionMap == fFractionMap.end())
|
---|
| 221 | {
|
---|
| 222 | itFractionMap = fFractionMap.find(0);
|
---|
| 223 | }
|
---|
| 224 |
|
---|
| 225 | fraction = itFractionMap->second;
|
---|
| 226 | fTowerFractions.push_back(fraction);
|
---|
[e5767b57] | 227 |
|
---|
[82575a3] | 228 | if(fraction < 1.0E-9) continue;
|
---|
| 229 |
|
---|
| 230 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 231 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
|
---|
| 232 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 233 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 234 |
|
---|
| 235 | // phi bins for given eta bin
|
---|
| 236 | phiBins = fPhiBins[etaBin];
|
---|
| 237 |
|
---|
| 238 | // find phi bin [1, phiBins.size - 1]
|
---|
| 239 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
|
---|
| 240 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 241 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 242 |
|
---|
| 243 | flags = 0;
|
---|
| 244 | flags |= (pdgCode == 11 || pdgCode == 22) << 1;
|
---|
| 245 |
|
---|
| 246 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
|
---|
| 247 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
| 248 |
|
---|
| 249 | fTowerHits.push_back(towerHit);
|
---|
| 250 | }
|
---|
| 251 |
|
---|
| 252 | // loop over all tracks
|
---|
| 253 | fItTrackInputArray->Reset();
|
---|
| 254 | number = -1;
|
---|
| 255 | while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
| 256 | {
|
---|
| 257 | const TLorentzVector &trackPosition = track->Position;
|
---|
| 258 | ++number;
|
---|
| 259 |
|
---|
| 260 | pdgCode = TMath::Abs(track->PID);
|
---|
| 261 |
|
---|
| 262 | itFractionMap = fFractionMap.find(pdgCode);
|
---|
| 263 | if(itFractionMap == fFractionMap.end())
|
---|
| 264 | {
|
---|
| 265 | itFractionMap = fFractionMap.find(0);
|
---|
| 266 | }
|
---|
| 267 |
|
---|
| 268 | fraction = itFractionMap->second;
|
---|
[e5767b57] | 269 |
|
---|
[82575a3] | 270 | fTrackFractions.push_back(fraction);
|
---|
[e5767b57] | 271 |
|
---|
[82575a3] | 272 | // find eta bin [1, fEtaBins.size - 1]
|
---|
| 273 | itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
|
---|
| 274 | if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
|
---|
| 275 | etaBin = distance(fEtaBins.begin(), itEtaBin);
|
---|
| 276 |
|
---|
| 277 | // phi bins for given eta bin
|
---|
| 278 | phiBins = fPhiBins[etaBin];
|
---|
| 279 |
|
---|
| 280 | // find phi bin [1, phiBins.size - 1]
|
---|
| 281 | itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
|
---|
| 282 | if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
|
---|
| 283 | phiBin = distance(phiBins->begin(), itPhiBin);
|
---|
| 284 |
|
---|
| 285 | flags = 1;
|
---|
| 286 |
|
---|
| 287 | // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
|
---|
| 288 | towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
|
---|
| 289 |
|
---|
| 290 | fTowerHits.push_back(towerHit);
|
---|
| 291 | }
|
---|
| 292 |
|
---|
| 293 | // all hits are sorted first by eta bin number, then by phi bin number,
|
---|
| 294 | // then by flags and then by particle or track number
|
---|
| 295 | sort(fTowerHits.begin(), fTowerHits.end());
|
---|
| 296 |
|
---|
| 297 | // loop over all hits
|
---|
| 298 | towerEtaPhi = 0;
|
---|
| 299 | fTower = 0;
|
---|
| 300 | for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
|
---|
| 301 | {
|
---|
| 302 | towerHit = (*itTowerHits);
|
---|
| 303 | flags = (towerHit >> 24) & 0x00000000000000FFLL;
|
---|
| 304 | number = (towerHit) & 0x0000000000FFFFFFLL;
|
---|
| 305 | hitEtaPhi = towerHit >> 32;
|
---|
| 306 |
|
---|
| 307 | if(towerEtaPhi != hitEtaPhi)
|
---|
| 308 | {
|
---|
| 309 | // switch to next tower
|
---|
| 310 | towerEtaPhi = hitEtaPhi;
|
---|
| 311 |
|
---|
| 312 | // finalize previous tower
|
---|
| 313 | FinalizeTower();
|
---|
| 314 |
|
---|
| 315 | // create new tower
|
---|
| 316 | fTower = factory->NewCandidate();
|
---|
| 317 |
|
---|
| 318 | phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
|
---|
| 319 | etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
|
---|
| 320 |
|
---|
| 321 | // phi bins for given eta bin
|
---|
| 322 | phiBins = fPhiBins[etaBin];
|
---|
| 323 |
|
---|
| 324 | // calculate eta and phi of the tower's center
|
---|
| 325 | fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
|
---|
| 326 | fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
|
---|
| 327 |
|
---|
| 328 | fTowerEdges[0] = fEtaBins[etaBin - 1];
|
---|
| 329 | fTowerEdges[1] = fEtaBins[etaBin];
|
---|
| 330 | fTowerEdges[2] = (*phiBins)[phiBin - 1];
|
---|
| 331 | fTowerEdges[3] = (*phiBins)[phiBin];
|
---|
| 332 |
|
---|
| 333 | fTowerEnergy = 0.0;
|
---|
| 334 | fTrackEnergy = 0.0;
|
---|
[e5767b57] | 335 |
|
---|
[82575a3] | 336 | fTowerTime = 0.0;
|
---|
| 337 | fTrackTime = 0.0;
|
---|
[e5767b57] | 338 |
|
---|
| 339 | fTowerTimeWeight = 0.0;
|
---|
| 340 |
|
---|
[82575a3] | 341 | fTowerTrackHits = 0;
|
---|
| 342 | fTowerPhotonHits = 0;
|
---|
[e5767b57] | 343 |
|
---|
[82575a3] | 344 | fTowerTrackArray->Clear();
|
---|
| 345 | }
|
---|
| 346 |
|
---|
| 347 | // check for track hits
|
---|
| 348 | if(flags & 1)
|
---|
| 349 | {
|
---|
| 350 | ++fTowerTrackHits;
|
---|
| 351 |
|
---|
| 352 | track = static_cast<Candidate*>(fTrackInputArray->At(number));
|
---|
| 353 | momentum = track->Momentum;
|
---|
| 354 | position = track->Position;
|
---|
[e5767b57] | 355 |
|
---|
[82575a3] | 356 | energy = momentum.E() * fTrackFractions[number];
|
---|
[e5767b57] | 357 |
|
---|
[82575a3] | 358 | fTrackEnergy += energy;
|
---|
[e5767b57] | 359 |
|
---|
[82575a3] | 360 | fTrackTime += TMath::Sqrt(energy)*position.T();
|
---|
[e5767b57] | 361 | fTrackTimeWeight += TMath::Sqrt(energy);
|
---|
| 362 |
|
---|
[82575a3] | 363 | fTowerTrackArray->Add(track);
|
---|
| 364 |
|
---|
| 365 | continue;
|
---|
| 366 | }
|
---|
[e5767b57] | 367 |
|
---|
[82575a3] | 368 | // check for photon and electron hits in current tower
|
---|
| 369 | if(flags & 2) ++fTowerPhotonHits;
|
---|
[e5767b57] | 370 |
|
---|
[82575a3] | 371 | particle = static_cast<Candidate*>(fParticleInputArray->At(number));
|
---|
| 372 | momentum = particle->Momentum;
|
---|
| 373 | position = particle->Position;
|
---|
| 374 |
|
---|
| 375 | // fill current tower
|
---|
| 376 | energy = momentum.E() * fTowerFractions[number];
|
---|
[e5767b57] | 377 |
|
---|
[82575a3] | 378 | fTowerEnergy += energy;
|
---|
[e5767b57] | 379 |
|
---|
[82575a3] | 380 | fTowerTime += TMath::Sqrt(energy)*position.T();
|
---|
[e5767b57] | 381 | fTowerTimeWeight += TMath::Sqrt(energy);
|
---|
| 382 |
|
---|
[82575a3] | 383 | fTower->AddCandidate(particle);
|
---|
| 384 | }
|
---|
| 385 |
|
---|
| 386 | // finalize last tower
|
---|
| 387 | FinalizeTower();
|
---|
| 388 | }
|
---|
| 389 |
|
---|
| 390 | //------------------------------------------------------------------------------
|
---|
| 391 |
|
---|
| 392 | void SimpleCalorimeter::FinalizeTower()
|
---|
| 393 | {
|
---|
| 394 | Candidate *tower;
|
---|
| 395 | Double_t energy, pt, eta, phi;
|
---|
| 396 | Double_t sigma;
|
---|
| 397 | Double_t time;
|
---|
| 398 |
|
---|
| 399 | if(!fTower) return;
|
---|
| 400 |
|
---|
| 401 | sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerEnergy);
|
---|
| 402 |
|
---|
| 403 | energy = LogNormal(fTowerEnergy, sigma);
|
---|
[e5767b57] | 404 |
|
---|
| 405 | time = (fTowerTimeWeight < 1.0E-09 ) ? 0.0 : fTowerTime/fTowerTimeWeight;
|
---|
[82575a3] | 406 |
|
---|
[4b9a2dc] | 407 | sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
|
---|
[e5767b57] | 408 |
|
---|
| 409 | if(energy < fEnergyMin || energy < fEnergySignificanceMin*sigma) energy = 0.0;
|
---|
| 410 |
|
---|
| 411 | if(fDitherTowerCenter)
|
---|
| 412 | {
|
---|
| 413 | eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
|
---|
| 414 | phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
|
---|
| 415 | }
|
---|
| 416 | else
|
---|
| 417 | {
|
---|
| 418 | eta = fTowerEta;
|
---|
| 419 | phi = fTowerPhi;
|
---|
| 420 | }
|
---|
[82575a3] | 421 |
|
---|
| 422 | pt = energy / TMath::CosH(eta);
|
---|
| 423 |
|
---|
| 424 | fTower->Position.SetPtEtaPhiE(1.0, eta, phi, time);
|
---|
| 425 | fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
|
---|
[e5767b57] | 426 |
|
---|
[82575a3] | 427 | fTower->Edges[0] = fTowerEdges[0];
|
---|
| 428 | fTower->Edges[1] = fTowerEdges[1];
|
---|
| 429 | fTower->Edges[2] = fTowerEdges[2];
|
---|
| 430 | fTower->Edges[3] = fTowerEdges[3];
|
---|
| 431 |
|
---|
| 432 | // fill SimpleCalorimeter towers
|
---|
| 433 | if(energy > 0.0) fTowerOutputArray->Add(fTower);
|
---|
| 434 |
|
---|
| 435 | // fill energy flow candidates
|
---|
| 436 | energy -= fTrackEnergy;
|
---|
[e5767b57] | 437 |
|
---|
| 438 | sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
|
---|
| 439 |
|
---|
| 440 | if(energy < fEnergyMin || energy < fEnergySignificanceMin*sigma) energy = 0.0;
|
---|
| 441 |
|
---|
[82575a3] | 442 | // save energy excess as an energy flow tower
|
---|
| 443 | if(energy > 0.0)
|
---|
| 444 | {
|
---|
| 445 | // create new photon tower
|
---|
| 446 | tower = static_cast<Candidate*>(fTower->Clone());
|
---|
| 447 | pt = energy / TMath::CosH(eta);
|
---|
| 448 |
|
---|
| 449 | tower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
|
---|
| 450 | fEFlowTowerOutputArray->Add(tower);
|
---|
| 451 | }
|
---|
| 452 | }
|
---|
| 453 |
|
---|
| 454 | //------------------------------------------------------------------------------
|
---|
| 455 |
|
---|
| 456 | Double_t SimpleCalorimeter::LogNormal(Double_t mean, Double_t sigma)
|
---|
| 457 | {
|
---|
| 458 | Double_t a, b;
|
---|
| 459 |
|
---|
| 460 | if(mean > 0.0)
|
---|
| 461 | {
|
---|
| 462 | b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
|
---|
| 463 | a = TMath::Log(mean) - 0.5*b*b;
|
---|
| 464 |
|
---|
[e5767b57] | 465 | return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
|
---|
[82575a3] | 466 | }
|
---|
| 467 | else
|
---|
| 468 | {
|
---|
| 469 | return 0.0;
|
---|
| 470 | }
|
---|
| 471 | }
|
---|