Fork me on GitHub

source: git/modules/DualReadoutCalorimeter.cc@ e7cccc2

3.4.3pre05
Last change on this file since e7cccc2 was f382b7e, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

now filling neutral hadrons in DR towers

  • Property mode set to 100644
File size: 23.4 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
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.
9 *
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.
14 *
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
19
20/** \class DualReadoutCalorimeter
21 *
22
23 // ============ TODO =========================================================
24 // This implementation of dual calorimetry relies on several approximations:
25 // - If hadronic energy is found in the tower the energy resolution then the full tower enrgy is smeared according to hadronic resolution (pessimistic for (e,n) or (pi+,gamma))
26 // - While e+ vs pi+ (or gamma vs n) separation is in principle possible for single particles (using C/S, PMT timing, lateral shower profile) it is not obvious it can be done overlapping particles.
27 // Now we assume that regarless of the number of particle hits per tower we can always distinguish e+ vs pi+, which is probably not true in the case (e+,n) vs (pi+,gamma) without longitudinal segmentation.
28
29 *
30 * \author M. Selvaggi - CERN
31 *
32 */
33
34#include "modules/DualReadoutCalorimeter.h"
35
36#include "classes/DelphesClasses.h"
37#include "classes/DelphesFactory.h"
38#include "classes/DelphesFormula.h"
39
40#include "ExRootAnalysis/ExRootResult.h"
41#include "ExRootAnalysis/ExRootFilter.h"
42#include "ExRootAnalysis/ExRootClassifier.h"
43
44#include "TMath.h"
45#include "TString.h"
46#include "TFormula.h"
47#include "TRandom3.h"
48#include "TObjArray.h"
49#include "TDatabasePDG.h"
50#include "TLorentzVector.h"
51
52#include <algorithm>
53#include <stdexcept>
54#include <iostream>
55#include <sstream>
56
57using namespace std;
58
59//------------------------------------------------------------------------------
60
61DualReadoutCalorimeter::DualReadoutCalorimeter() :
62 fECalResolutionFormula(0), fHCalResolutionFormula(0),
63 fItParticleInputArray(0), fItTrackInputArray(0)
64{
65
66 fECalResolutionFormula = new DelphesFormula;
67 fHCalResolutionFormula = new DelphesFormula;
68
69 fECalTowerTrackArray = new TObjArray;
70 fItECalTowerTrackArray = fECalTowerTrackArray->MakeIterator();
71
72 fHCalTowerTrackArray = new TObjArray;
73 fItHCalTowerTrackArray = fHCalTowerTrackArray->MakeIterator();
74
75 fTowerTrackArray = new TObjArray;
76 fItTowerTrackArray = fTowerTrackArray->MakeIterator();
77
78}
79
80//------------------------------------------------------------------------------
81
82DualReadoutCalorimeter::~DualReadoutCalorimeter()
83{
84
85 if(fECalResolutionFormula) delete fECalResolutionFormula;
86 if(fHCalResolutionFormula) delete fHCalResolutionFormula;
87
88 if(fECalTowerTrackArray) delete fECalTowerTrackArray;
89 if(fItECalTowerTrackArray) delete fItECalTowerTrackArray;
90
91 if(fHCalTowerTrackArray) delete fHCalTowerTrackArray;
92 if(fItHCalTowerTrackArray) delete fItHCalTowerTrackArray;
93
94 if(fTowerTrackArray) delete fTowerTrackArray;
95 if(fItTowerTrackArray) delete fItTowerTrackArray;
96
97}
98
99//------------------------------------------------------------------------------
100
101void DualReadoutCalorimeter::Init()
102{
103 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
104 Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
105 Double_t ecalFraction, hcalFraction;
106 TBinMap::iterator itEtaBin;
107 set< Double_t >::iterator itPhiBin;
108 vector< Double_t > *phiBins;
109
110 // read eta and phi bins
111 param = GetParam("EtaPhiBins");
112 size = param.GetSize();
113 fBinMap.clear();
114 fEtaBins.clear();
115 fPhiBins.clear();
116 for(i = 0; i < size/2; ++i)
117 {
118 paramEtaBins = param[i*2];
119 sizeEtaBins = paramEtaBins.GetSize();
120 paramPhiBins = param[i*2 + 1];
121 sizePhiBins = paramPhiBins.GetSize();
122
123 for(j = 0; j < sizeEtaBins; ++j)
124 {
125 for(k = 0; k < sizePhiBins; ++k)
126 {
127 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
128 }
129 }
130 }
131
132 // for better performance we transform map of sets to parallel vectors:
133 // vector< double > and vector< vector< double >* >
134 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
135 {
136 fEtaBins.push_back(itEtaBin->first);
137 phiBins = new vector< double >(itEtaBin->second.size());
138 fPhiBins.push_back(phiBins);
139 phiBins->clear();
140 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
141 {
142 phiBins->push_back(*itPhiBin);
143 }
144 }
145
146 // read energy fractions for different particles
147 param = GetParam("EnergyFraction");
148 size = param.GetSize();
149
150 // set default energy fractions values
151 fFractionMap.clear();
152 fFractionMap[0] = make_pair(0.0, 1.0);
153
154 for(i = 0; i < size/2; ++i)
155 {
156 paramFractions = param[i*2 + 1];
157
158 ecalFraction = paramFractions[0].GetDouble();
159 hcalFraction = paramFractions[1].GetDouble();
160
161 fFractionMap[param[i*2].GetInt()] = make_pair(ecalFraction, hcalFraction);
162 }
163
164 // read min E value for timing measurement in ECAL
165 fTimingEnergyMin = GetDouble("TimingEnergyMin",4.);
166 // For timing
167 // So far this flag needs to be false
168 // Curved extrapolation not supported
169 fElectronsFromTrack = false;
170
171 // read min E value for towers to be saved
172 fECalEnergyMin = GetDouble("ECalEnergyMin", 0.0);
173 fHCalEnergyMin = GetDouble("HCalEnergyMin", 0.0);
174 fEnergyMin = GetDouble("EnergyMin", 0.0);
175
176 fECalEnergySignificanceMin = GetDouble("ECalEnergySignificanceMin", 0.0);
177 fHCalEnergySignificanceMin = GetDouble("HCalEnergySignificanceMin", 0.0);
178 fEnergySignificanceMin = GetDouble("EnergySignificanceMin", 0.0);
179
180 // switch on or off the dithering of the center of DualReadoutCalorimeter towers
181 fSmearTowerCenter = GetBool("SmearTowerCenter", true);
182
183 // read resolution formulas
184 fECalResolutionFormula->Compile(GetString("ECalResolutionFormula", "0"));
185 fHCalResolutionFormula->Compile(GetString("HCalResolutionFormula", "0"));
186
187 // import array with output from other modules
188 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
189 fItParticleInputArray = fParticleInputArray->MakeIterator();
190
191 fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
192 fItTrackInputArray = fTrackInputArray->MakeIterator();
193
194 // create output arrays
195 fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
196 fPhotonOutputArray = ExportArray(GetString("PhotonOutputArray", "photons"));
197
198 fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
199 fEFlowPhotonOutputArray = ExportArray(GetString("EFlowPhotonOutputArray", "eflowPhotons"));
200 fEFlowNeutralHadronOutputArray = ExportArray(GetString("EFlowNeutralHadronOutputArray", "eflowNeutralHadrons"));
201}
202
203//------------------------------------------------------------------------------
204
205void DualReadoutCalorimeter::Finish()
206{
207 vector< vector< Double_t >* >::iterator itPhiBin;
208 if(fItParticleInputArray) delete fItParticleInputArray;
209 if(fItTrackInputArray) delete fItTrackInputArray;
210 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
211 {
212 delete *itPhiBin;
213 }
214}
215
216//------------------------------------------------------------------------------
217
218void DualReadoutCalorimeter::Process()
219{
220 Candidate *particle, *track;
221 TLorentzVector position, momentum;
222 Short_t etaBin, phiBin, flags;
223 Int_t number;
224 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
225 Double_t ecalFraction, hcalFraction;
226 Double_t ecalEnergy, hcalEnergy;
227 Double_t ecalSigma, hcalSigma, sigma;
228 Double_t energyGuess, energy;
229 Int_t pdgCode;
230
231 TFractionMap::iterator itFractionMap;
232
233 vector< Double_t >::iterator itEtaBin;
234 vector< Double_t >::iterator itPhiBin;
235 vector< Double_t > *phiBins;
236
237 vector< Long64_t >::iterator itTowerHits;
238
239 DelphesFactory *factory = GetFactory();
240 fTowerHits.clear();
241 fECalTowerFractions.clear();
242 fHCalTowerFractions.clear();
243 fECalTrackFractions.clear();
244 fHCalTrackFractions.clear();
245
246 // loop over all particles
247 fItParticleInputArray->Reset();
248 number = -1;
249 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
250 {
251 const TLorentzVector &particlePosition = particle->Position;
252 ++number;
253
254 pdgCode = TMath::Abs(particle->PID);
255
256 itFractionMap = fFractionMap.find(pdgCode);
257 if(itFractionMap == fFractionMap.end())
258 {
259 itFractionMap = fFractionMap.find(0);
260 }
261
262 ecalFraction = itFractionMap->second.first;
263 hcalFraction = itFractionMap->second.second;
264
265 fECalTowerFractions.push_back(ecalFraction);
266 fHCalTowerFractions.push_back(hcalFraction);
267
268 if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
269
270 // find eta bin [1, fEtaBins.size - 1]
271 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
272 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
273 etaBin = distance(fEtaBins.begin(), itEtaBin);
274
275 // phi bins for given eta bin
276 phiBins = fPhiBins[etaBin];
277
278 // find phi bin [1, phiBins.size - 1]
279 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
280 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
281 phiBin = distance(phiBins->begin(), itPhiBin);
282
283 flags = 0;
284 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
285
286 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
287 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
288
289 fTowerHits.push_back(towerHit);
290 }
291
292 // loop over all tracks
293 fItTrackInputArray->Reset();
294 number = -1;
295 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
296 {
297 const TLorentzVector &trackPosition = track->Position;
298 ++number;
299
300 pdgCode = TMath::Abs(track->PID);
301
302 itFractionMap = fFractionMap.find(pdgCode);
303 if(itFractionMap == fFractionMap.end())
304 {
305 itFractionMap = fFractionMap.find(0);
306 }
307
308 ecalFraction = itFractionMap->second.first;
309 hcalFraction = itFractionMap->second.second;
310
311 fECalTrackFractions.push_back(ecalFraction);
312 fHCalTrackFractions.push_back(hcalFraction);
313
314 // find eta bin [1, fEtaBins.size - 1]
315 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
316 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
317 etaBin = distance(fEtaBins.begin(), itEtaBin);
318
319 // phi bins for given eta bin
320 phiBins = fPhiBins[etaBin];
321
322 // find phi bin [1, phiBins.size - 1]
323 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
324 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
325 phiBin = distance(phiBins->begin(), itPhiBin);
326
327 flags = 1;
328
329 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
330 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
331
332 fTowerHits.push_back(towerHit);
333 }
334
335 // all hits are sorted first by eta bin number, then by phi bin number,
336 // then by flags and then by particle or track number
337 sort(fTowerHits.begin(), fTowerHits.end());
338
339 // loop over all hits
340 towerEtaPhi = 0;
341 fTower = 0;
342 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
343 {
344 towerHit = (*itTowerHits);
345 flags = (towerHit >> 24) & 0x00000000000000FFLL;
346 number = (towerHit) & 0x0000000000FFFFFFLL;
347 hitEtaPhi = towerHit >> 32;
348
349 if(towerEtaPhi != hitEtaPhi)
350 {
351 // switch to next tower
352 towerEtaPhi = hitEtaPhi;
353
354 // finalize previous tower
355 FinalizeTower();
356
357 // create new tower
358 fTower = factory->NewCandidate();
359
360 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
361 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
362
363 // phi bins for given eta bin
364 phiBins = fPhiBins[etaBin];
365
366 // calculate eta and phi of the tower's center
367 fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
368 fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
369
370 fTowerEdges[0] = fEtaBins[etaBin - 1];
371 fTowerEdges[1] = fEtaBins[etaBin];
372 fTowerEdges[2] = (*phiBins)[phiBin - 1];
373 fTowerEdges[3] = (*phiBins)[phiBin];
374
375 fECalTowerEnergy = 0.0;
376 fHCalTowerEnergy = 0.0;
377
378 fECalTrackEnergy = 0.0;
379 fHCalTrackEnergy = 0.0;
380 fTrackEnergy = 0.0;
381
382 fECalTrackSigma = 0.0;
383 fHCalTrackSigma = 0.0;
384 fTrackSigma = 0.0;
385
386 fTowerTrackHits = 0;
387 fTowerPhotonHits = 0;
388
389 fECalTowerTrackArray->Clear();
390 fHCalTowerTrackArray->Clear();
391 fTowerTrackArray->Clear();
392
393 }
394
395 // check for track hits
396 if(flags & 1)
397 {
398 ++fTowerTrackHits;
399
400 track = static_cast<Candidate*>(fTrackInputArray->At(number));
401 momentum = track->Momentum;
402 position = track->Position;
403
404 ecalEnergy = momentum.E() * fECalTrackFractions[number];
405 hcalEnergy = momentum.E() * fHCalTrackFractions[number];
406 energy = ecalEnergy + hcalEnergy;
407
408 if(ecalEnergy > fTimingEnergyMin && fTower)
409 {
410 if(fElectronsFromTrack)
411 {
412 fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, track->Position.T()));
413 }
414 }
415
416
417 /*
418 if(fECalTrackFractions[number] > 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
419 {
420 fECalTrackEnergy += ecalEnergy;
421 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
422 if(ecalSigma/momentum.E() < track->TrackResolution) energyGuess = ecalEnergy;
423 else energyGuess = momentum.E();
424
425 fECalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
426 fECalTowerTrackArray->Add(track);
427 }
428
429 else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] > 1.0E-9)
430 {
431 fHCalTrackEnergy += hcalEnergy;
432 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
433 if(hcalSigma/momentum.E() < track->TrackResolution) energyGuess = hcalEnergy;
434 else energyGuess = momentum.E();
435
436 fHCalTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
437 fHCalTowerTrackArray->Add(track);
438 }
439
440 // muons
441 else if(fECalTrackFractions[number] < 1.0E-9 && fHCalTrackFractions[number] < 1.0E-9)
442 {
443 fEFlowTrackOutputArray->Add(track);
444 }
445 */
446
447 // in Dual Readout we do not care if tracks are ECAL of HCAL
448 if(fECalTrackFractions[number] > 1.0E-9 || fHCalTrackFractions[number] > 1.0E-9)
449 {
450 fTrackEnergy += energy;
451 // this sigma will be used to determine whether neutral excess is significant. We choose the resolution according to bthe higest deposited fraction (in practice had for charged hadrons and em for electrons)
452 sigma = 0.0;
453 if(fHCalTrackFractions[number] > 0)
454 sigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
455 else
456 sigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
457
458 if(sigma/momentum.E() < track->TrackResolution)
459 energyGuess = ecalEnergy + hcalEnergy;
460 else
461 energyGuess = momentum.E();
462
463 fTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
464 fTowerTrackArray->Add(track);
465
466 }
467 else
468 {
469 fEFlowTrackOutputArray->Add(track);
470 }
471
472 continue;
473 }
474
475 // check for photon and electron hits in current tower
476 if(flags & 2) ++fTowerPhotonHits;
477
478 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
479 momentum = particle->Momentum;
480 position = particle->Position;
481
482 // fill current tower
483 ecalEnergy = momentum.E() * fECalTowerFractions[number];
484 hcalEnergy = momentum.E() * fHCalTowerFractions[number];
485
486 fECalTowerEnergy += ecalEnergy;
487 fHCalTowerEnergy += hcalEnergy;
488
489 if(ecalEnergy > fTimingEnergyMin && fTower)
490 {
491 if (abs(particle->PID) != 11 || !fElectronsFromTrack)
492 {
493 fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, particle->Position.T()));
494 }
495 }
496
497 fTower->AddCandidate(particle);
498 }
499
500 // finalize last tower
501 FinalizeTower();
502}
503
504//------------------------------------------------------------------------------
505
506void DualReadoutCalorimeter::FinalizeTower()
507{
508
509 Candidate *track, *tower, *mother;
510 Double_t energy, pt, eta, phi;
511 Double_t ecalEnergy, hcalEnergy;
512 Double_t ecalNeutralEnergy, hcalNeutralEnergy, neutralEnergy;
513
514 Double_t ecalSigma, hcalSigma, sigma;
515 Double_t ecalNeutralSigma, hcalNeutralSigma, neutralSigma;
516
517 Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
518
519 TLorentzVector momentum;
520 TFractionMap::iterator itFractionMap;
521
522 Float_t weight, sumWeightedTime, sumWeight;
523
524 if(!fTower) return;
525
526
527 //if (fHCalTowerEnergy < 30 && fECalTowerEnergy < 30) return;
528 //cout<<"----------- New tower ---------"<<endl;
529
530
531 // here we change behaviour w.r.t to standard calorimeter. Start with worse case scenario. If fHCalTowerEnergy > 0, assume total energy smeared by HCAL resolution.
532 // For example, if overlapping charged pions and photons take hadronic resolution as combined measurement
533
534 // if no hadronic fraction at all, then use ECAL resolution
535
536 //cout<<"fECalTowerEnergy: "<<fECalTowerEnergy<<", fHCalTowerEnergy: "<<fHCalTowerEnergy<<", Eta: "<<fTowerEta<<endl;
537
538 // if no hadronic energy, use ECAL resolution
539 if (fHCalTowerEnergy <= fHCalEnergyMin)
540 {
541 energy = fECalTowerEnergy;
542 sigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
543 //cout<<"em energy"<<energy<<", sigma: "<<sigma<<endl;
544 }
545
546 // if hadronic fraction > 0, use HCAL resolution
547 else
548 {
549 energy = fECalTowerEnergy + fHCalTowerEnergy;
550 sigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
551 //cout<<"had energy: "<<energy<<", sigma: "<<sigma<<endl;
552 }
553
554 energy = LogNormal(energy, sigma);
555 //cout<<energy<<","<<ecalEnergy<<","<<hcalEnergy<<endl;
556
557 if(energy < fEnergyMin || energy < fEnergySignificanceMin*sigma) energy = 0.0;
558
559 // for now keep this the same
560 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fECalTowerEnergy);
561 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fHCalTowerEnergy);
562
563 ecalEnergy = LogNormal(fECalTowerEnergy, ecalSigma);
564 hcalEnergy = LogNormal(fHCalTowerEnergy, hcalSigma);
565
566 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
567 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
568
569 if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
570 if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
571
572 //cout<<"Measured energy: "<<energy<<endl;
573
574 if(fSmearTowerCenter)
575 {
576 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
577 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
578 }
579 else
580 {
581 eta = fTowerEta;
582 phi = fTowerPhi;
583 }
584
585 pt = energy / TMath::CosH(eta);
586
587 // Time calculation for tower
588 fTower->NTimeHits = 0;
589 sumWeightedTime = 0.0;
590 sumWeight = 0.0;
591
592 for(size_t i = 0; i < fTower->ECalEnergyTimePairs.size(); ++i)
593 {
594 weight = TMath::Sqrt(fTower->ECalEnergyTimePairs[i].first);
595 sumWeightedTime += weight * fTower->ECalEnergyTimePairs[i].second;
596 sumWeight += weight;
597 fTower->NTimeHits++;
598 }
599
600 if(sumWeight > 0.0)
601 {
602 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, sumWeightedTime/sumWeight);
603 }
604 else
605 {
606 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, 999999.9);
607 }
608
609 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
610 fTower->Eem = ecalEnergy;
611 fTower->Ehad = hcalEnergy;
612
613 fTower->Edges[0] = fTowerEdges[0];
614 fTower->Edges[1] = fTowerEdges[1];
615 fTower->Edges[2] = fTowerEdges[2];
616 fTower->Edges[3] = fTowerEdges[3];
617
618 if(energy > 0.0)
619 {
620 if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
621 {
622 fPhotonOutputArray->Add(fTower);
623 }
624 fTowerOutputArray->Add(fTower);
625 }
626
627 // fill energy flow candidates
628
629 fTrackSigma = TMath::Sqrt(fTrackSigma);
630 neutralEnergy = max( (energy - fTrackEnergy) , 0.0);
631 neutralSigma = neutralEnergy / TMath::Sqrt(fTrackSigma*fTrackSigma + sigma*sigma);
632
633 //cout<<"trackEnergy: "<<fTrackEnergy<<", trackSigma: "<<fTrackSigma<<", Ntracks: "<<fTowerTrackArray->GetEntries()<<endl;
634
635 //cout<<"neutralEnergy: "<<neutralEnergy<<", neutralSigma: "<<neutralSigma<<", :fEnergyMin "<<fEnergyMin<<", fEnergySignificanceMin: "<<fEnergySignificanceMin<<endl;
636
637 // For now, if neutral excess is significant, simply create neutral EflowPhoton tower and clone each track into eflowtrack !!! -> Creating only photons !! EFlowNeutralHadron collection will be empy!!! TO BE FIXED
638 if(neutralEnergy > fEnergyMin && neutralSigma > fEnergySignificanceMin)
639 {
640
641 //cout<<"significant neutral excess found:"<<endl;
642 // create new photon tower
643 tower = static_cast<Candidate*>(fTower->Clone());
644 pt = neutralEnergy / TMath::CosH(eta);
645 //cout<<"Creating tower with Pt, Eta, Phi, Energy: "<<pt<<","<<eta<<","<<phi<<","<<neutralEnergy<<endl;
646 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, neutralEnergy);
647
648 // if no hadronic energy, use ECAL resolution
649 if (fHCalTowerEnergy <= fHCalEnergyMin)
650 {
651 tower->Eem = neutralEnergy;
652 tower->Ehad = 0.0;
653 tower->PID = 22;
654 }
655
656 // if hadronic fraction > 0, use HCAL resolution
657 else
658 {
659 tower->Eem = 0;
660 tower->Ehad = neutralEnergy;
661 tower->PID = 130;
662 }
663
664 fEFlowPhotonOutputArray->Add(tower);
665
666
667 //clone tracks
668 fItTowerTrackArray->Reset();
669 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
670 {
671 //cout<<"looping over tracks"<<endl;
672 mother = track;
673 track = static_cast<Candidate*>(track->Clone());
674 track->AddCandidate(mother);
675 fEFlowTrackOutputArray->Add(track);
676 }
677 }
678
679
680 // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the DualReadoutCalorimeter and tracking
681 else if(fTrackEnergy > 0.0)
682 {
683 //cout<<"no significant neutral excess found:"<<endl;
684 weightTrack = (fTrackSigma > 0.0) ? 1 / (fTrackSigma*fTrackSigma) : 0.0;
685 weightCalo = (sigma > 0.0) ? 1 / (sigma*sigma) : 0.0;
686
687 bestEnergyEstimate = (weightTrack*fTrackEnergy + weightCalo*energy) / (weightTrack + weightCalo);
688 rescaleFactor = bestEnergyEstimate/fTrackEnergy;
689
690 //rescale tracks
691 fItTowerTrackArray->Reset();
692 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
693 {
694 mother = track;
695 track = static_cast<Candidate*>(track->Clone());
696 track->AddCandidate(mother);
697
698 track->Momentum *= rescaleFactor;
699
700 fEFlowTrackOutputArray->Add(track);
701 }
702 }
703
704
705}
706
707//------------------------------------------------------------------------------
708
709Double_t DualReadoutCalorimeter::LogNormal(Double_t mean, Double_t sigma)
710{
711 Double_t a, b;
712
713 if(mean > 0.0)
714 {
715 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
716 a = TMath::Log(mean) - 0.5*b*b;
717
718 return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
719 }
720 else
721 {
722 return 0.0;
723 }
724}
Note: See TracBrowser for help on using the repository browser.