Fork me on GitHub

source: git/modules/SimpleCalorimeter.cc@ e70228d

Timing
Last change on this file since e70228d was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 15.5 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/** \class SimpleCalorimeter
20 *
21 * Fills SimpleCalorimeter towers, performs SimpleCalorimeter resolution smearing,
22 * and creates energy flow objects (tracks, photons, and neutral hadrons).
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/SimpleCalorimeter.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootClassifier.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootResult.h"
37
38#include "TDatabasePDG.h"
39#include "TFormula.h"
40#include "TLorentzVector.h"
41#include "TMath.h"
42#include "TObjArray.h"
43#include "TRandom3.h"
44#include "TString.h"
45
46#include <algorithm>
47#include <iostream>
48#include <sstream>
49#include <stdexcept>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55SimpleCalorimeter::SimpleCalorimeter() :
56 fResolutionFormula(0),
57 fItParticleInputArray(0), fItTrackInputArray(0)
58{
59
60 fResolutionFormula = new DelphesFormula;
61 fTowerTrackArray = new TObjArray;
62 fItTowerTrackArray = fTowerTrackArray->MakeIterator();
63}
64
65//------------------------------------------------------------------------------
66
67SimpleCalorimeter::~SimpleCalorimeter()
68{
69
70 if(fResolutionFormula) delete fResolutionFormula;
71 if(fTowerTrackArray) delete fTowerTrackArray;
72 if(fItTowerTrackArray) delete fItTowerTrackArray;
73}
74
75//------------------------------------------------------------------------------
76
77void SimpleCalorimeter::Init()
78{
79 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
80 Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
81 Double_t fraction;
82 TBinMap::iterator itEtaBin;
83 set<Double_t>::iterator itPhiBin;
84 vector<Double_t> *phiBins;
85
86 // read eta and phi bins
87 param = GetParam("EtaPhiBins");
88 size = param.GetSize();
89 fBinMap.clear();
90 fEtaBins.clear();
91 fPhiBins.clear();
92 for(i = 0; i < size / 2; ++i)
93 {
94 paramEtaBins = param[i * 2];
95 sizeEtaBins = paramEtaBins.GetSize();
96 paramPhiBins = param[i * 2 + 1];
97 sizePhiBins = paramPhiBins.GetSize();
98
99 for(j = 0; j < sizeEtaBins; ++j)
100 {
101 for(k = 0; k < sizePhiBins; ++k)
102 {
103 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
104 }
105 }
106 }
107
108 // for better performance we transform map of sets to parallel vectors:
109 // vector< double > and vector< vector< double >* >
110 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
111 {
112 fEtaBins.push_back(itEtaBin->first);
113 phiBins = new vector<double>(itEtaBin->second.size());
114 fPhiBins.push_back(phiBins);
115 phiBins->clear();
116 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
117 {
118 phiBins->push_back(*itPhiBin);
119 }
120 }
121
122 // read energy fractions for different particles
123 param = GetParam("EnergyFraction");
124 size = param.GetSize();
125
126 // set default energy fractions values
127 fFractionMap.clear();
128 fFractionMap[0] = 1.0;
129
130 for(i = 0; i < size / 2; ++i)
131 {
132 paramFractions = param[i * 2 + 1];
133 fraction = paramFractions[0].GetDouble();
134 fFractionMap[param[i * 2].GetInt()] = fraction;
135 }
136
137 // read min E value for towers to be saved
138 fEnergyMin = GetDouble("EnergyMin", 0.0);
139
140 fEnergySignificanceMin = GetDouble("EnergySignificanceMin", 0.0);
141
142 // flag that says if current calo is Ecal of Hcal (will then fill correct values of Eem and Ehad)
143 fIsEcal = GetBool("IsEcal", false);
144
145 // switch on or off the dithering of the center of calorimeter towers
146 fSmearTowerCenter = GetBool("SmearTowerCenter", true);
147
148 // read resolution formulas
149 fResolutionFormula->Compile(GetString("ResolutionFormula", "0"));
150
151 // import array with output from other modules
152 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
153 fItParticleInputArray = fParticleInputArray->MakeIterator();
154
155 fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
156 fItTrackInputArray = fTrackInputArray->MakeIterator();
157
158 // create output arrays
159 fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
160
161 fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
162 fEFlowTowerOutputArray = ExportArray(GetString("EFlowTowerOutputArray", "eflowTowers"));
163}
164
165//------------------------------------------------------------------------------
166
167void SimpleCalorimeter::Finish()
168{
169 vector<vector<Double_t> *>::iterator itPhiBin;
170 if(fItParticleInputArray) delete fItParticleInputArray;
171 if(fItTrackInputArray) delete fItTrackInputArray;
172 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
173 {
174 delete *itPhiBin;
175 }
176}
177
178//------------------------------------------------------------------------------
179
180void SimpleCalorimeter::Process()
181{
182 Candidate *particle, *track;
183 TLorentzVector position, momentum;
184 Short_t etaBin, phiBin, flags;
185 Int_t number;
186 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
187 Double_t fraction;
188 Double_t energy;
189 Double_t sigma;
190 Double_t energyGuess;
191
192 Int_t pdgCode;
193
194 TFractionMap::iterator itFractionMap;
195
196 vector<Double_t>::iterator itEtaBin;
197 vector<Double_t>::iterator itPhiBin;
198 vector<Double_t> *phiBins;
199
200 vector<Long64_t>::iterator itTowerHits;
201
202 DelphesFactory *factory = GetFactory();
203 fTowerHits.clear();
204 fTowerFractions.clear();
205 fTrackFractions.clear();
206
207 // loop over all particles
208 fItParticleInputArray->Reset();
209 number = -1;
210 while((particle = static_cast<Candidate *>(fItParticleInputArray->Next())))
211 {
212 const TLorentzVector &particlePosition = particle->Position;
213 ++number;
214
215 pdgCode = TMath::Abs(particle->PID);
216
217 itFractionMap = fFractionMap.find(pdgCode);
218 if(itFractionMap == fFractionMap.end())
219 {
220 itFractionMap = fFractionMap.find(0);
221 }
222
223 fraction = itFractionMap->second;
224 fTowerFractions.push_back(fraction);
225
226 if(fraction < 1.0E-9) continue;
227
228 // find eta bin [1, fEtaBins.size - 1]
229 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
230 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
231 etaBin = distance(fEtaBins.begin(), itEtaBin);
232
233 // phi bins for given eta bin
234 phiBins = fPhiBins[etaBin];
235
236 // find phi bin [1, phiBins.size - 1]
237 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
238 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
239 phiBin = distance(phiBins->begin(), itPhiBin);
240
241 flags = 0;
242 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
243
244 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
245 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
246
247 fTowerHits.push_back(towerHit);
248 }
249
250 // loop over all tracks
251 fItTrackInputArray->Reset();
252 number = -1;
253 while((track = static_cast<Candidate *>(fItTrackInputArray->Next())))
254 {
255 const TLorentzVector &trackPosition = track->Position;
256 ++number;
257
258 pdgCode = TMath::Abs(track->PID);
259
260 itFractionMap = fFractionMap.find(pdgCode);
261 if(itFractionMap == fFractionMap.end())
262 {
263 itFractionMap = fFractionMap.find(0);
264 }
265
266 fraction = itFractionMap->second;
267
268 fTrackFractions.push_back(fraction);
269
270 // find eta bin [1, fEtaBins.size - 1]
271 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.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(), trackPosition.Phi());
280 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
281 phiBin = distance(phiBins->begin(), itPhiBin);
282
283 flags = 1;
284
285 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
286 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
287
288 fTowerHits.push_back(towerHit);
289 }
290
291 // all hits are sorted first by eta bin number, then by phi bin number,
292 // then by flags and then by particle or track number
293 sort(fTowerHits.begin(), fTowerHits.end());
294
295 // loop over all hits
296 towerEtaPhi = 0;
297 fTower = 0;
298 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
299 {
300 towerHit = (*itTowerHits);
301 flags = (towerHit >> 24) & 0x00000000000000FFLL;
302 number = (towerHit)&0x0000000000FFFFFFLL;
303 hitEtaPhi = towerHit >> 32;
304
305 if(towerEtaPhi != hitEtaPhi)
306 {
307 // switch to next tower
308 towerEtaPhi = hitEtaPhi;
309
310 // finalize previous tower
311 FinalizeTower();
312
313 // create new tower
314 fTower = factory->NewCandidate();
315
316 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
317 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
318
319 // phi bins for given eta bin
320 phiBins = fPhiBins[etaBin];
321
322 // calculate eta and phi of the tower's center
323 fTowerEta = 0.5 * (fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
324 fTowerPhi = 0.5 * ((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
325
326 fTowerEdges[0] = fEtaBins[etaBin - 1];
327 fTowerEdges[1] = fEtaBins[etaBin];
328 fTowerEdges[2] = (*phiBins)[phiBin - 1];
329 fTowerEdges[3] = (*phiBins)[phiBin];
330
331 fTowerEnergy = 0.0;
332
333 fTrackEnergy = 0.0;
334 fTrackSigma = 0.0;
335
336 fTowerTime = 0.0;
337 fTrackTime = 0.0;
338
339 fTowerTimeWeight = 0.0;
340
341 fTowerTrackHits = 0;
342 fTowerPhotonHits = 0;
343
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;
355
356 energy = momentum.E() * fTrackFractions[number];
357
358 fTrackTime += TMath::Sqrt(energy) * position.T();
359 fTrackTimeWeight += TMath::Sqrt(energy);
360
361 if(fTrackFractions[number] > 1.0E-9)
362 {
363
364 // compute total charged energy
365 fTrackEnergy += energy;
366 sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, momentum.E());
367 if(sigma / momentum.E() < track->TrackResolution)
368 energyGuess = energy;
369 else
370 energyGuess = momentum.E();
371
372 fTrackSigma += ((track->TrackResolution) * energyGuess) * ((track->TrackResolution) * energyGuess);
373 fTowerTrackArray->Add(track);
374 }
375
376 else
377 {
378 fEFlowTrackOutputArray->Add(track);
379 }
380
381 continue;
382 }
383
384 // check for photon and electron hits in current tower
385 if(flags & 2) ++fTowerPhotonHits;
386
387 particle = static_cast<Candidate *>(fParticleInputArray->At(number));
388 momentum = particle->Momentum;
389 position = particle->Position;
390
391 // fill current tower
392 energy = momentum.E() * fTowerFractions[number];
393
394 fTowerEnergy += energy;
395
396 fTowerTime += energy * position.T();
397 fTowerTimeWeight += energy;
398
399 fTower->AddCandidate(particle);
400 }
401
402 // finalize last tower
403 FinalizeTower();
404}
405
406//------------------------------------------------------------------------------
407
408void SimpleCalorimeter::FinalizeTower()
409{
410 Candidate *tower, *track, *mother;
411 Double_t energy, neutralEnergy, pt, eta, phi;
412 Double_t sigma, neutralSigma;
413 Double_t time;
414
415 Double_t weightTrack, weightCalo, bestEnergyEstimate, rescaleFactor;
416
417 TLorentzVector momentum;
418 TFractionMap::iterator itFractionMap;
419
420 if(!fTower) return;
421
422 sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerEnergy);
423
424 energy = LogNormal(fTowerEnergy, sigma);
425
426 time = (fTowerTimeWeight < 1.0E-09) ? 0.0 : fTowerTime / fTowerTimeWeight;
427
428 sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, energy);
429
430 if(energy < fEnergyMin || energy < fEnergySignificanceMin * sigma) energy = 0.0;
431
432 if(fSmearTowerCenter)
433 {
434 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
435 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
436 }
437 else
438 {
439 eta = fTowerEta;
440 phi = fTowerPhi;
441 }
442
443 pt = energy / TMath::CosH(eta);
444
445 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, time);
446 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
447
448 fTower->Eem = (!fIsEcal) ? 0 : energy;
449 fTower->Ehad = (fIsEcal) ? 0 : energy;
450
451 fTower->Edges[0] = fTowerEdges[0];
452 fTower->Edges[1] = fTowerEdges[1];
453 fTower->Edges[2] = fTowerEdges[2];
454 fTower->Edges[3] = fTowerEdges[3];
455
456 // fill SimpleCalorimeter towers
457 if(energy > 0.0) fTowerOutputArray->Add(fTower);
458
459 // e-flow candidates
460
461 //compute neutral excess
462
463 fTrackSigma = TMath::Sqrt(fTrackSigma);
464 neutralEnergy = max((energy - fTrackEnergy), 0.0);
465
466 //compute sigma_trk total
467 neutralSigma = neutralEnergy / TMath::Sqrt(fTrackSigma * fTrackSigma + sigma * sigma);
468
469 // if neutral excess is significant, simply create neutral Eflow tower and clone each track into eflowtrack
470 if(neutralEnergy > fEnergyMin && neutralSigma > fEnergySignificanceMin)
471 {
472 // create new photon tower
473 tower = static_cast<Candidate *>(fTower->Clone());
474 pt = neutralEnergy / TMath::CosH(eta);
475
476 tower->Eem = (!fIsEcal) ? 0 : neutralEnergy;
477 tower->Ehad = (fIsEcal) ? 0 : neutralEnergy;
478 tower->PID = (fIsEcal) ? 22 : 0;
479
480 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, neutralEnergy);
481 fEFlowTowerOutputArray->Add(tower);
482
483 fItTowerTrackArray->Reset();
484 while((track = static_cast<Candidate *>(fItTowerTrackArray->Next())))
485 {
486 mother = track;
487 track = static_cast<Candidate *>(track->Clone());
488 track->AddCandidate(mother);
489
490 fEFlowTrackOutputArray->Add(track);
491 }
492 }
493
494 // if neutral excess is not significant, rescale eflow tracks, such that the total charged equals the best measurement given by the calorimeter and tracking
495 else if(fTrackEnergy > 0.0)
496 {
497 weightTrack = (fTrackSigma > 0.0) ? 1 / (fTrackSigma * fTrackSigma) : 0.0;
498 weightCalo = (sigma > 0.0) ? 1 / (sigma * sigma) : 0.0;
499
500 bestEnergyEstimate = (weightTrack * fTrackEnergy + weightCalo * energy) / (weightTrack + weightCalo);
501 rescaleFactor = bestEnergyEstimate / fTrackEnergy;
502
503 fItTowerTrackArray->Reset();
504 while((track = static_cast<Candidate *>(fItTowerTrackArray->Next())))
505 {
506 mother = track;
507 track = static_cast<Candidate *>(track->Clone());
508 track->AddCandidate(mother);
509
510 track->Momentum *= rescaleFactor;
511
512 fEFlowTrackOutputArray->Add(track);
513 }
514 }
515}
516
517//------------------------------------------------------------------------------
518
519Double_t SimpleCalorimeter::LogNormal(Double_t mean, Double_t sigma)
520{
521 Double_t a, b;
522
523 if(mean > 0.0)
524 {
525 b = TMath::Sqrt(TMath::Log((1.0 + (sigma * sigma) / (mean * mean))));
526 a = TMath::Log(mean) - 0.5 * b * b;
527
528 return TMath::Exp(a + b * gRandom->Gaus(0.0, 1.0));
529 }
530 else
531 {
532 return 0.0;
533 }
534}
Note: See TracBrowser for help on using the repository browser.