Fork me on GitHub

source: git/modules/Calorimeter.cc@ e0f8f99

ImprovedOutputFile Timing dual_readout llp
Last change on this file since e0f8f99 was efac6f9, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

added photon PID to ecal eflow towers

  • Property mode set to 100644
File size: 20.5 KB
RevLine 
[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
52using namespace std;
53
54//------------------------------------------------------------------------------
55
56Calorimeter::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
74Calorimeter::~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
90void 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
192void 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
205void 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;
[d7d2da3]214 Int_t pdgCode;
215
216 TFractionMap::iterator itFractionMap;
217
218 vector< Double_t >::iterator itEtaBin;
219 vector< Double_t >::iterator itPhiBin;
220 vector< Double_t > *phiBins;
221
222 vector< Long64_t >::iterator itTowerHits;
223
224 DelphesFactory *factory = GetFactory();
225 fTowerHits.clear();
[00e8dca]226 fECalTowerFractions.clear();
227 fHCalTowerFractions.clear();
228 fECalTrackFractions.clear();
229 fHCalTrackFractions.clear();
[82575a3]230
[d7d2da3]231 // loop over all particles
232 fItParticleInputArray->Reset();
233 number = -1;
234 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
235 {
236 const TLorentzVector &particlePosition = particle->Position;
237 ++number;
238
239 pdgCode = TMath::Abs(particle->PID);
240
241 itFractionMap = fFractionMap.find(pdgCode);
242 if(itFractionMap == fFractionMap.end())
243 {
244 itFractionMap = fFractionMap.find(0);
245 }
246
[82575a3]247 ecalFraction = itFractionMap->second.first;
248 hcalFraction = itFractionMap->second.second;
249
[00e8dca]250 fECalTowerFractions.push_back(ecalFraction);
251 fHCalTowerFractions.push_back(hcalFraction);
[82575a3]252
253 if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
[d7d2da3]254
255 // find eta bin [1, fEtaBins.size - 1]
256 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
257 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
258 etaBin = distance(fEtaBins.begin(), itEtaBin);
259
260 // phi bins for given eta bin
261 phiBins = fPhiBins[etaBin];
262
263 // find phi bin [1, phiBins.size - 1]
264 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
265 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
266 phiBin = distance(phiBins->begin(), itPhiBin);
267
[73e0386]268 flags = 0;
[2dab783]269 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
[d7d2da3]270
271 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
272 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
273
274 fTowerHits.push_back(towerHit);
275 }
276
277 // loop over all tracks
278 fItTrackInputArray->Reset();
279 number = -1;
280 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
281 {
282 const TLorentzVector &trackPosition = track->Position;
283 ++number;
284
[73e0386]285 pdgCode = TMath::Abs(track->PID);
286
287 itFractionMap = fFractionMap.find(pdgCode);
288 if(itFractionMap == fFractionMap.end())
289 {
290 itFractionMap = fFractionMap.find(0);
291 }
292
[82575a3]293 ecalFraction = itFractionMap->second.first;
294 hcalFraction = itFractionMap->second.second;
295
[00e8dca]296 fECalTrackFractions.push_back(ecalFraction);
297 fHCalTrackFractions.push_back(hcalFraction);
[82575a3]298
[d7d2da3]299 // find eta bin [1, fEtaBins.size - 1]
300 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
301 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
302 etaBin = distance(fEtaBins.begin(), itEtaBin);
303
304 // phi bins for given eta bin
305 phiBins = fPhiBins[etaBin];
306
307 // find phi bin [1, phiBins.size - 1]
308 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
309 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
310 phiBin = distance(phiBins->begin(), itPhiBin);
311
[73e0386]312 flags = 1;
313
[d7d2da3]314 // 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]315 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
[d7d2da3]316
317 fTowerHits.push_back(towerHit);
318 }
319
320 // all hits are sorted first by eta bin number, then by phi bin number,
321 // then by flags and then by particle or track number
322 sort(fTowerHits.begin(), fTowerHits.end());
323
324 // loop over all hits
325 towerEtaPhi = 0;
326 fTower = 0;
327 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
328 {
329 towerHit = (*itTowerHits);
330 flags = (towerHit >> 24) & 0x00000000000000FFLL;
331 number = (towerHit) & 0x0000000000FFFFFFLL;
332 hitEtaPhi = towerHit >> 32;
333
334 if(towerEtaPhi != hitEtaPhi)
335 {
336 // switch to next tower
337 towerEtaPhi = hitEtaPhi;
338
339 // finalize previous tower
340 FinalizeTower();
341
342 // create new tower
343 fTower = factory->NewCandidate();
344
345 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
346 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
347
348 // phi bins for given eta bin
349 phiBins = fPhiBins[etaBin];
350
351 // calculate eta and phi of the tower's center
352 fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
353 fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
354
355 fTowerEdges[0] = fEtaBins[etaBin - 1];
356 fTowerEdges[1] = fEtaBins[etaBin];
357 fTowerEdges[2] = (*phiBins)[phiBin - 1];
358 fTowerEdges[3] = (*phiBins)[phiBin];
359
[00e8dca]360 fECalTowerEnergy = 0.0;
361 fHCalTowerEnergy = 0.0;
362
[04290b1]363 fECalTrackEnergy = 0.0;
364 fHCalTrackEnergy = 0.0;
365
366 fECalTrackSigma = 0.0;
367 fHCalTrackSigma = 0.0;
368
[2dab783]369 fTowerTrackHits = 0;
370 fTowerPhotonHits = 0;
[a221d1f]371
[04290b1]372 fECalTowerTrackArray->Clear();
373 fHCalTowerTrackArray->Clear();
374
[d7d2da3]375 }
376
377 // check for track hits
[73e0386]378 if(flags & 1)
[d7d2da3]379 {
[2dab783]380 ++fTowerTrackHits;
[73e0386]381
[2dab783]382 track = static_cast<Candidate*>(fTrackInputArray->At(number));
383 momentum = track->Momentum;
[22dc7fd]384 position = track->Position;
[82575a3]385
[00e8dca]386 ecalEnergy = momentum.E() * fECalTrackFractions[number];
387 hcalEnergy = momentum.E() * fHCalTrackFractions[number];
[a221d1f]388
[839deb7]389 if(ecalEnergy > fTimingEnergyMin && fTower)
390 {
391 if(fElectronsFromTrack)
392 {
393 fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, track->Position.T()));
394 }
[3db5282]395 }
[a221d1f]396
[9da65a5]397 if(fECalTrackFractions[number] > 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
[00e8dca]398 {
[04290b1]399 fECalTrackEnergy += ecalEnergy;
400 fECalTrackSigma += (track->TrackResolution)*momentum.E()*(track->TrackResolution)*momentum.E();
401 fECalTowerTrackArray->Add(track);
[00e8dca]402 }
[04290b1]403
[9da65a5]404 else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] > 1.0E-9)
[00e8dca]405 {
[04290b1]406 fHCalTrackEnergy += hcalEnergy;
407 fHCalTrackSigma += (track->TrackResolution)*momentum.E()*(track->TrackResolution)*momentum.E();
408 fHCalTowerTrackArray->Add(track);
[00e8dca]409 }
[04290b1]410
[9da65a5]411 else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
412 {
413 fEFlowTrackOutputArray->Add(track);
414 }
[73e0386]415
[2dab783]416 continue;
[73e0386]417 }
[839deb7]418
[73e0386]419 // check for photon and electron hits in current tower
[2dab783]420 if(flags & 2) ++fTowerPhotonHits;
[a221d1f]421
[d7d2da3]422 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
423 momentum = particle->Momentum;
[22dc7fd]424 position = particle->Position;
[d7d2da3]425
426 // fill current tower
[00e8dca]427 ecalEnergy = momentum.E() * fECalTowerFractions[number];
428 hcalEnergy = momentum.E() * fHCalTowerFractions[number];
[82575a3]429
[00e8dca]430 fECalTowerEnergy += ecalEnergy;
431 fHCalTowerEnergy += hcalEnergy;
[82575a3]432
[839deb7]433 if(ecalEnergy > fTimingEnergyMin && fTower)
434 {
435 if (abs(particle->PID) != 11 || !fElectronsFromTrack)
436 {
437 fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, particle->Position.T()));
[3db5282]438 }
439 }
[82575a3]440
[3079350]441 fTower->AddCandidate(particle);
[d7d2da3]442 }
443
444 // finalize last tower
445 FinalizeTower();
446}
447
448//------------------------------------------------------------------------------
449
450void Calorimeter::FinalizeTower()
451{
[e2dd4c5]452 Candidate *track, *tower, *mother;
[d7d2da3]453 Double_t energy, pt, eta, phi;
[82575a3]454 Double_t ecalEnergy, hcalEnergy;
[04290b1]455 Double_t ecalNeutralEnergy, hcalNeutralEnergy;
456
[82575a3]457 Double_t ecalSigma, hcalSigma;
[04290b1]458 Double_t ecalNeutralSigma, hcalNeutralSigma;
[e2dd4c5]459
[04290b1]460 Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
461
[a98c7ef]462 TLorentzVector momentum;
463 TFractionMap::iterator itFractionMap;
[e2dd4c5]464
[839deb7]465 Float_t weight, sumWeightedTime, sumWeight;
466
[d7d2da3]467 if(!fTower) return;
468
[00e8dca]469 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fECalTowerEnergy);
470 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fHCalTowerEnergy);
[4600a41]471
[00e8dca]472 ecalEnergy = LogNormal(fECalTowerEnergy, ecalSigma);
473 hcalEnergy = LogNormal(fHCalTowerEnergy, hcalSigma);
[82575a3]474
[4b9a2dc]475 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
476 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
477
[38bf1ae]478 if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
479 if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
[4b9a2dc]480
[82575a3]481 energy = ecalEnergy + hcalEnergy;
[839deb7]482
[4e09c3a]483 if(fSmearTowerCenter)
[a221d1f]484 {
[38bf1ae]485 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
486 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
[a221d1f]487 }
488 else
489 {
[38bf1ae]490 eta = fTowerEta;
491 phi = fTowerPhi;
[a221d1f]492 }
[d7d2da3]493
494 pt = energy / TMath::CosH(eta);
495
[3db5282]496 // Time calculation for tower
[839deb7]497 fTower->NTimeHits = 0;
498 sumWeightedTime = 0.0;
499 sumWeight = 0.0;
500
501 for(size_t i = 0; i < fTower->ECalEnergyTimePairs.size(); ++i)
[3db5282]502 {
[839deb7]503 weight = TMath::Sqrt(fTower->ECalEnergyTimePairs[i].first);
504 sumWeightedTime += weight * fTower->ECalEnergyTimePairs[i].second;
505 sumWeight += weight;
506 fTower->NTimeHits++;
[3db5282]507 }
[839deb7]508
509 if(sumWeight > 0.0)
510 {
511 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, sumWeightedTime/sumWeight);
512 }
513 else
514 {
515 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, 999999.9);
[3db5282]516 }
517
518
[d7d2da3]519 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
[82575a3]520 fTower->Eem = ecalEnergy;
521 fTower->Ehad = hcalEnergy;
522
[d7d2da3]523 fTower->Edges[0] = fTowerEdges[0];
524 fTower->Edges[1] = fTowerEdges[1];
525 fTower->Edges[2] = fTowerEdges[2];
526 fTower->Edges[3] = fTowerEdges[3];
527
[38bf1ae]528 if(energy > 0.0)
[82575a3]529 {
530 if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
531 {
532 fPhotonOutputArray->Add(fTower);
533 }
[a221d1f]534
[82575a3]535 fTowerOutputArray->Add(fTower);
536 }
[04290b1]537
[f6b9fec]538 // fill energy flow candidates
[04290b1]539 fECalTrackSigma = TMath::Sqrt(fECalTrackSigma);
540 fHCalTrackSigma = TMath::Sqrt(fHCalTrackSigma);
541
542 //compute neutral excesses
543 ecalNeutralEnergy = max( (ecalEnergy - fECalTrackEnergy) , 0.0);
544 hcalNeutralEnergy = max( (hcalEnergy - fHCalTrackEnergy) , 0.0);
545
546 ecalNeutralSigma = ecalNeutralEnergy / TMath::Sqrt(fECalTrackSigma*fECalTrackSigma + ecalSigma*ecalSigma);
547 hcalNeutralSigma = hcalNeutralEnergy / TMath::Sqrt(fHCalTrackSigma*fHCalTrackSigma + hcalSigma*hcalSigma);
548
549 // if ecal neutral excess is significant, simply create neutral EflowPhoton tower and clone each track into eflowtrack
550 if(ecalNeutralEnergy > fECalEnergyMin && ecalNeutralSigma > fECalEnergySignificanceMin)
[00e8dca]551 {
[04290b1]552 // create new photon tower
553 tower = static_cast<Candidate*>(fTower->Clone());
554 pt = ecalNeutralEnergy / TMath::CosH(eta);
555
556 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalNeutralEnergy);
557 tower->Eem = ecalNeutralEnergy;
558 tower->Ehad = 0.0;
[efac6f9]559 tower->PID = 22;
[04290b1]560
561 fEFlowPhotonOutputArray->Add(tower);
562
563 //clone tracks
564 fItECalTowerTrackArray->Reset();
565 while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
566 {
567 mother = track;
568 track = static_cast<Candidate*>(track->Clone());
569 track->AddCandidate(mother);
[e2dd4c5]570
[04290b1]571 fEFlowTrackOutputArray->Add(track);
572 }
573
[82575a3]574 }
[04290b1]575
576 // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
577 else if(fECalTrackEnergy > 0.0)
[00e8dca]578 {
[04290b1]579 weightTrack = (fECalTrackSigma > 0.0) ? 1 / (fECalTrackSigma*fECalTrackSigma) : 0.0;
580 weightCalo = (ecalSigma > 0.0) ? 1 / (ecalSigma*ecalSigma) : 0.0;
581
582 bestEnergyEstimate = (weightTrack*fECalTrackEnergy + weightCalo*ecalEnergy) / (weightTrack + weightCalo);
583 rescaleFactor = bestEnergyEstimate/fECalTrackEnergy;
584
585 //rescale tracks
586 fItECalTowerTrackArray->Reset();
587 while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
588 {
589 mother = track;
590 track = static_cast<Candidate*>(track->Clone());
591 track->AddCandidate(mother);
592
593 track->Momentum *= rescaleFactor;
594
595 fEFlowTrackOutputArray->Add(track);
596 }
[00e8dca]597 }
[e2dd4c5]598
[82575a3]599
[04290b1]600 // if hcal neutral excess is significant, simply create neutral EflowNeutralHadron tower and clone each track into eflowtrack
601 if(hcalNeutralEnergy > fHCalEnergyMin && hcalNeutralSigma > fHCalEnergySignificanceMin)
[2dab783]602 {
[d4c4d9d]603 // create new photon tower
[2dab783]604 tower = static_cast<Candidate*>(fTower->Clone());
[04290b1]605 pt = hcalNeutralEnergy / TMath::CosH(eta);
606
607 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalNeutralEnergy);
608 tower->Ehad = hcalNeutralEnergy;
609 tower->Eem = 0.0;
610
611 fEFlowNeutralHadronOutputArray->Add(tower);
612
613 //clone tracks
614 fItHCalTowerTrackArray->Reset();
615 while((track = static_cast<Candidate*>(fItHCalTowerTrackArray->Next())))
616 {
617 mother = track;
618 track = static_cast<Candidate*>(track->Clone());
619 track->AddCandidate(mother);
[0d5f77c]620
[04290b1]621 fEFlowTrackOutputArray->Add(track);
622 }
623
[82575a3]624 }
[04290b1]625
626 // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
627 else if(fHCalTrackEnergy > 0.0)
[82575a3]628 {
[04290b1]629 weightTrack = (fHCalTrackSigma > 0.0) ? 1 / (fHCalTrackSigma*fHCalTrackSigma) : 0.0;
630 weightCalo = (hcalSigma > 0.0) ? 1 / (hcalSigma*hcalSigma) : 0.0;
631
632 bestEnergyEstimate = (weightTrack*fHCalTrackEnergy + weightCalo*hcalEnergy) / (weightTrack + weightCalo);
633 rescaleFactor = bestEnergyEstimate / fHCalTrackEnergy;
634
635 //rescale tracks
636 fItECalTowerTrackArray->Reset();
637 while((track = static_cast<Candidate*>(fItECalTowerTrackArray->Next())))
638 {
639 mother = track;
640 track = static_cast<Candidate*>(track->Clone());
641 track->AddCandidate(mother);
642
643 track->Momentum *= rescaleFactor;
644
645 fEFlowTrackOutputArray->Add(track);
646 }
[d7d2da3]647 }
[04290b1]648
649
[d7d2da3]650}
651
652//------------------------------------------------------------------------------
[39022b6]653
[4600a41]654Double_t Calorimeter::LogNormal(Double_t mean, Double_t sigma)
655{
656 Double_t a, b;
657
658 if(mean > 0.0)
[39022b6]659 {
[4600a41]660 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
661 a = TMath::Log(mean) - 0.5*b*b;
662
[38bf1ae]663 return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
[39022b6]664 }
[4600a41]665 else
666 {
667 return 0.0;
668 }
669}
Note: See TracBrowser for help on using the repository browser.