Fork me on GitHub

source: git/modules/DualReadoutCalorimeter.cc@ 4e8e72b

Last change on this file since 4e8e72b was 9a7ea36, checked in by michele <michele.selvaggi@…>, 3 years ago

removed trailing comments

  • Property mode set to 100644
File size: 21.0 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 fTowerRmax=0.;
250
251 //cout<<"--------- new event ---------- "<<endl;
252
253 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
254 {
255 const TLorentzVector &particlePosition = particle->Position;
256 ++number;
257
258 // compute maximum radius (needed in FinalizeTower to assess whether barrel or endcap tower)
259 if (particlePosition.Perp() > fTowerRmax)
260 fTowerRmax=particlePosition.Perp();
261
262 pdgCode = TMath::Abs(particle->PID);
263
264 itFractionMap = fFractionMap.find(pdgCode);
265 if(itFractionMap == fFractionMap.end())
266 {
267 itFractionMap = fFractionMap.find(0);
268 }
269
270 ecalFraction = itFractionMap->second.first;
271 hcalFraction = itFractionMap->second.second;
272
273 fECalTowerFractions.push_back(ecalFraction);
274 fHCalTowerFractions.push_back(hcalFraction);
275
276 if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
277
278 // find eta bin [1, fEtaBins.size - 1]
279 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
280 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
281 etaBin = distance(fEtaBins.begin(), itEtaBin);
282
283 // phi bins for given eta bin
284 phiBins = fPhiBins[etaBin];
285
286 // find phi bin [1, phiBins.size - 1]
287 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
288 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
289 phiBin = distance(phiBins->begin(), itPhiBin);
290
291 flags = 0;
292 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
293
294 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
295 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
296
297 fTowerHits.push_back(towerHit);
298 }
299
300 // loop over all tracks
301 fItTrackInputArray->Reset();
302 number = -1;
303 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
304 {
305 const TLorentzVector &trackPosition = track->Position;
306 ++number;
307
308 pdgCode = TMath::Abs(track->PID);
309
310 itFractionMap = fFractionMap.find(pdgCode);
311 if(itFractionMap == fFractionMap.end())
312 {
313 itFractionMap = fFractionMap.find(0);
314 }
315
316 ecalFraction = itFractionMap->second.first;
317 hcalFraction = itFractionMap->second.second;
318
319 fECalTrackFractions.push_back(ecalFraction);
320 fHCalTrackFractions.push_back(hcalFraction);
321
322 // find eta bin [1, fEtaBins.size - 1]
323 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
324 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
325 etaBin = distance(fEtaBins.begin(), itEtaBin);
326
327 // phi bins for given eta bin
328 phiBins = fPhiBins[etaBin];
329
330 // find phi bin [1, phiBins.size - 1]
331 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
332 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
333 phiBin = distance(phiBins->begin(), itPhiBin);
334
335 flags = 1;
336
337 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
338 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
339
340 fTowerHits.push_back(towerHit);
341 }
342
343 // all hits are sorted first by eta bin number, then by phi bin number,
344 // then by flags and then by particle or track number
345 sort(fTowerHits.begin(), fTowerHits.end());
346
347 // loop over all hits
348 towerEtaPhi = 0;
349 fTower = 0;
350 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
351 {
352 towerHit = (*itTowerHits);
353 flags = (towerHit >> 24) & 0x00000000000000FFLL;
354 number = (towerHit) & 0x0000000000FFFFFFLL;
355 hitEtaPhi = towerHit >> 32;
356
357 if(towerEtaPhi != hitEtaPhi)
358 {
359 // switch to next tower
360 towerEtaPhi = hitEtaPhi;
361
362 // finalize previous tower
363 FinalizeTower();
364
365 // create new tower
366 fTower = factory->NewCandidate();
367
368 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
369 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
370
371 // phi bins for given eta bin
372 phiBins = fPhiBins[etaBin];
373
374 // calculate eta and phi of the tower's center
375 fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
376 fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
377
378 fTowerEdges[0] = fEtaBins[etaBin - 1];
379 fTowerEdges[1] = fEtaBins[etaBin];
380 fTowerEdges[2] = (*phiBins)[phiBin - 1];
381 fTowerEdges[3] = (*phiBins)[phiBin];
382
383 fECalTowerEnergy = 0.0;
384 fHCalTowerEnergy = 0.0;
385
386 fECalTrackEnergy = 0.0;
387 fHCalTrackEnergy = 0.0;
388 fTrackEnergy = 0.0;
389
390 fECalTrackSigma = 0.0;
391 fHCalTrackSigma = 0.0;
392 fTrackSigma = 0.0;
393
394 fTowerTrackHits = 0;
395 fTowerPhotonHits = 0;
396
397 fTowerTime = 0.0;
398 fTowerTimeWeight = 0.0;
399
400 fECalTowerTrackArray->Clear();
401 fHCalTowerTrackArray->Clear();
402 fTowerTrackArray->Clear();
403
404 }
405
406 // check for track hits
407 if(flags & 1)
408 {
409 ++fTowerTrackHits;
410
411 track = static_cast<Candidate*>(fTrackInputArray->At(number));
412 momentum = track->Momentum;
413 position = track->Position;
414
415 ecalEnergy = momentum.E() * fECalTrackFractions[number];
416 hcalEnergy = momentum.E() * fHCalTrackFractions[number];
417 energy = ecalEnergy + hcalEnergy;
418
419 if(ecalEnergy > fTimingEnergyMin && fTower)
420 {
421 if(fElectronsFromTrack)
422 {
423 fTower->ECalEnergyTimePairs.push_back(make_pair<Float_t, Float_t>(ecalEnergy, track->Position.T()));
424 }
425 }
426
427 // in Dual Readout we do not care if tracks are ECAL of HCAL
428 if(fECalTrackFractions[number] > 1.0E-9 || fHCalTrackFractions[number] > 1.0E-9)
429 {
430 fTrackEnergy += energy;
431 // 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)
432 sigma = 0.0;
433 if(fHCalTrackFractions[number] > 0)
434 sigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
435 else
436 sigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
437
438 if(sigma/momentum.E() < track->TrackResolution)
439 energyGuess = ecalEnergy + hcalEnergy;
440 else
441 energyGuess = momentum.E();
442
443 fTrackSigma += (track->TrackResolution)*energyGuess*(track->TrackResolution)*energyGuess;
444 fTowerTrackArray->Add(track);
445 }
446 else
447 {
448 fEFlowTrackOutputArray->Add(track);
449 }
450
451 continue;
452 }
453
454 // check for photon and electron hits in current tower
455 if(flags & 2) ++fTowerPhotonHits;
456
457 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
458 momentum = particle->Momentum;
459 position = particle->Position;
460
461
462 // fill current tower
463 ecalEnergy = momentum.E() * fECalTowerFractions[number];
464 hcalEnergy = momentum.E() * fHCalTowerFractions[number];
465
466 fECalTowerEnergy += ecalEnergy;
467 fHCalTowerEnergy += hcalEnergy;
468
469 // assume combined timing measurements in ECAL/HCAL sections
470 fTowerTime += (ecalEnergy + hcalEnergy) * position.T(); //sigma_t ~ 1/sqrt(E)
471 fTowerTimeWeight += ecalEnergy + hcalEnergy;
472
473 fTower->AddCandidate(particle);
474 fTower->Position = position;
475 }
476
477 // finalize last tower
478 FinalizeTower();
479}
480
481//------------------------------------------------------------------------------
482
483void DualReadoutCalorimeter::FinalizeTower()
484{
485
486 Candidate *track, *tower, *mother;
487 Double_t energy, pt, eta, phi, r, time;
488 Double_t ecalEnergy, hcalEnergy;
489 Double_t ecalNeutralEnergy, hcalNeutralEnergy, neutralEnergy;
490
491 Double_t ecalSigma, hcalSigma, sigma;
492 Double_t ecalNeutralSigma, hcalNeutralSigma, neutralSigma;
493
494 Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
495
496 TLorentzVector momentum;
497 TFractionMap::iterator itFractionMap;
498
499 Float_t weight, sumWeightedTime, sumWeight;
500
501 if(!fTower) return;
502
503 // if no hadronic energy, use ECAL resolution
504 if (fHCalTowerEnergy <= fHCalEnergyMin)
505 {
506 energy = fECalTowerEnergy;
507 sigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
508 }
509
510 // if hadronic fraction > 0, use HCAL resolution
511 else
512 {
513 energy = fECalTowerEnergy + fHCalTowerEnergy;
514 sigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
515 }
516
517 energy = LogNormal(energy, sigma);
518
519 if(energy < fEnergyMin || energy < fEnergySignificanceMin*sigma) energy = 0.0;
520
521 // for now keep this the same
522 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fECalTowerEnergy);
523 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fHCalTowerEnergy);
524
525 ecalEnergy = LogNormal(fECalTowerEnergy, ecalSigma);
526 hcalEnergy = LogNormal(fHCalTowerEnergy, hcalSigma);
527
528 time = (fTowerTimeWeight < 1.0E-09) ? 0.0 : fTowerTime / fTowerTimeWeight;
529
530 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
531 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
532
533 if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
534 if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
535
536 if(fSmearTowerCenter)
537 {
538 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
539 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
540 }
541 else
542 {
543 eta = fTowerEta;
544 phi = fTowerPhi;
545 }
546
547 pt = energy / TMath::CosH(eta);
548
549 // check whether barrel or endcap tower
550
551 // endcap
552 if (TMath::Abs(fTower->Position.Pt() - fTowerRmax) > 1.e-06 && TMath::Abs(eta) > 0.){
553 r = fTower->Position.Z()/TMath::SinH(eta);
554 }
555 // barrel
556 else {
557 r = fTower->Position.Pt();
558 }
559
560 fTower->Position.SetPtEtaPhiE(r, eta, phi, time);
561 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
562 fTower->L = fTower->Position.Vect().Mag();
563
564 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
565 fTower->Eem = ecalEnergy;
566 fTower->Ehad = hcalEnergy;
567 fTower->Etrk = fTrackEnergy;
568 fTower->Edges[0] = fTowerEdges[0];
569 fTower->Edges[1] = fTowerEdges[1];
570 fTower->Edges[2] = fTowerEdges[2];
571 fTower->Edges[3] = fTowerEdges[3];
572
573 if(energy > 0.0)
574 {
575 if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
576 {
577 fPhotonOutputArray->Add(fTower);
578 }
579 fTowerOutputArray->Add(fTower);
580 }
581
582 // fill energy flow candidates
583
584 fTrackSigma = TMath::Sqrt(fTrackSigma);
585 neutralEnergy = max( (energy - fTrackEnergy) , 0.0);
586 neutralSigma = neutralEnergy / TMath::Sqrt(fTrackSigma*fTrackSigma + sigma*sigma);
587
588 if(neutralEnergy > fEnergyMin && neutralSigma > fEnergySignificanceMin)
589 {
590 // create new photon tower
591 tower = static_cast<Candidate*>(fTower->Clone());
592 pt = neutralEnergy / TMath::CosH(eta);
593 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, neutralEnergy);
594
595 // if no hadronic energy, use ECAL resolution
596 if (fHCalTowerEnergy <= fHCalEnergyMin)
597 {
598 tower->Eem = neutralEnergy;
599 tower->Ehad = 0.0;
600 tower->PID = 22;
601 fEFlowPhotonOutputArray->Add(tower);
602 }
603 // if hadronic fraction > 0, use HCAL resolution
604 else
605 {
606 tower->Eem = 0;
607 tower->Ehad = neutralEnergy;
608 tower->PID = 130;
609 fEFlowNeutralHadronOutputArray->Add(tower);
610 }
611
612 //clone tracks
613 fItTowerTrackArray->Reset();
614 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
615 {
616 mother = track;
617 track = static_cast<Candidate*>(track->Clone());
618 track->AddCandidate(mother);
619 fEFlowTrackOutputArray->Add(track);
620 }
621 }
622
623
624 // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the DualReadoutCalorimeter and tracking
625 else if(fTrackEnergy > 0.0)
626 {
627 //cout<<"no significant neutral excess found:"<<endl;
628 weightTrack = (fTrackSigma > 0.0) ? 1 / (fTrackSigma*fTrackSigma) : 0.0;
629 weightCalo = (sigma > 0.0) ? 1 / (sigma*sigma) : 0.0;
630
631 bestEnergyEstimate = (weightTrack*fTrackEnergy + weightCalo*energy) / (weightTrack + weightCalo);
632 rescaleFactor = bestEnergyEstimate/fTrackEnergy;
633
634 //rescale tracks
635 fItTowerTrackArray->Reset();
636 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
637 {
638 mother = track;
639 track = static_cast<Candidate *>(track->Clone());
640 track->AddCandidate(mother);
641 track->Momentum.SetPtEtaPhiM(track->Momentum.Pt()*rescaleFactor, track->Momentum.Eta(), track->Momentum.Phi(), track->Momentum.M());
642 fEFlowTrackOutputArray->Add(track);
643 }
644 }
645
646
647}
648
649//------------------------------------------------------------------------------
650
651Double_t DualReadoutCalorimeter::LogNormal(Double_t mean, Double_t sigma)
652{
653 Double_t a, b;
654
655 if(mean > 0.0)
656 {
657 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
658 a = TMath::Log(mean) - 0.5*b*b;
659
660 return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
661 }
662 else
663 {
664 return 0.0;
665 }
666}
Note: See TracBrowser for help on using the repository browser.