Fork me on GitHub

source: git/modules/Calorimeter.cc@ a24426f

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a24426f was 4e09c3a, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

replace DitherTowerCenter with SmearTowerCenter

  • Property mode set to 100644
File size: 16.8 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 Calorimeter
21 *
22 * Fills calorimeter towers, performs calorimeter 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/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() :
57 fECalResolutionFormula(0), fHCalResolutionFormula(0),
58 fItParticleInputArray(0), fItTrackInputArray(0),
59 fTowerTrackArray(0), fItTowerTrackArray(0)
60{
61 fECalResolutionFormula = new DelphesFormula;
62 fHCalResolutionFormula = new DelphesFormula;
63
64 fTowerTrackArray = new TObjArray;
65 fItTowerTrackArray = fTowerTrackArray->MakeIterator();
66}
67
68//------------------------------------------------------------------------------
69
70Calorimeter::~Calorimeter()
71{
72 if(fECalResolutionFormula) delete fECalResolutionFormula;
73 if(fHCalResolutionFormula) delete fHCalResolutionFormula;
74
75 if(fTowerTrackArray) delete fTowerTrackArray;
76 if(fItTowerTrackArray) delete fItTowerTrackArray;
77}
78
79//------------------------------------------------------------------------------
80
81void Calorimeter::Init()
82{
83 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
84 Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
85 Double_t ecalFraction, hcalFraction;
86 TBinMap::iterator itEtaBin;
87 set< Double_t >::iterator itPhiBin;
88 vector< Double_t > *phiBins;
89
90 // read eta and phi bins
91 param = GetParam("EtaPhiBins");
92 size = param.GetSize();
93 fBinMap.clear();
94 fEtaBins.clear();
95 fPhiBins.clear();
96 for(i = 0; i < size/2; ++i)
97 {
98 paramEtaBins = param[i*2];
99 sizeEtaBins = paramEtaBins.GetSize();
100 paramPhiBins = param[i*2 + 1];
101 sizePhiBins = paramPhiBins.GetSize();
102
103 for(j = 0; j < sizeEtaBins; ++j)
104 {
105 for(k = 0; k < sizePhiBins; ++k)
106 {
107 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
108 }
109 }
110 }
111
112 // for better performance we transform map of sets to parallel vectors:
113 // vector< double > and vector< vector< double >* >
114 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
115 {
116 fEtaBins.push_back(itEtaBin->first);
117 phiBins = new vector< double >(itEtaBin->second.size());
118 fPhiBins.push_back(phiBins);
119 phiBins->clear();
120 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
121 {
122 phiBins->push_back(*itPhiBin);
123 }
124 }
125
126 // read energy fractions for different particles
127 param = GetParam("EnergyFraction");
128 size = param.GetSize();
129
130 // set default energy fractions values
131 fFractionMap.clear();
132 fFractionMap[0] = make_pair(0.0, 1.0);
133
134 for(i = 0; i < size/2; ++i)
135 {
136 paramFractions = param[i*2 + 1];
137
138 ecalFraction = paramFractions[0].GetDouble();
139 hcalFraction = paramFractions[1].GetDouble();
140
141 fFractionMap[param[i*2].GetInt()] = make_pair(ecalFraction, hcalFraction);
142 }
143
144/*
145 TFractionMap::iterator itFractionMap;
146 for(itFractionMap = fFractionMap.begin(); itFractionMap != fFractionMap.end(); ++itFractionMap)
147 {
148 cout << itFractionMap->first << " " << itFractionMap->second.first << " " << itFractionMap->second.second << endl;
149 }
150*/
151
152 // read min E value for towers to be saved
153 fECalEnergyMin = GetDouble("ECalEnergyMin", 0.0);
154 fHCalEnergyMin = GetDouble("HCalEnergyMin", 0.0);
155
156 fECalEnergySignificanceMin = GetDouble("ECalEnergySignificanceMin", 0.0);
157 fHCalEnergySignificanceMin = GetDouble("HCalEnergySignificanceMin", 0.0);
158
159 // switch on or off the dithering of the center of calorimeter towers
160 fSmearTowerCenter = GetBool("SmearTowerCenter", true);
161
162 // read resolution formulas
163 fECalResolutionFormula->Compile(GetString("ECalResolutionFormula", "0"));
164 fHCalResolutionFormula->Compile(GetString("HCalResolutionFormula", "0"));
165
166 // import array with output from other modules
167 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
168 fItParticleInputArray = fParticleInputArray->MakeIterator();
169
170 fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
171 fItTrackInputArray = fTrackInputArray->MakeIterator();
172
173 // create output arrays
174 fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
175 fPhotonOutputArray = ExportArray(GetString("PhotonOutputArray", "photons"));
176
177 fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
178 fEFlowPhotonOutputArray = ExportArray(GetString("EFlowPhotonOutputArray", "eflowPhotons"));
179 fEFlowNeutralHadronOutputArray = ExportArray(GetString("EFlowNeutralHadronOutputArray", "eflowNeutralHadrons"));
180}
181
182//------------------------------------------------------------------------------
183
184void Calorimeter::Finish()
185{
186 vector< vector< Double_t >* >::iterator itPhiBin;
187 if(fItParticleInputArray) delete fItParticleInputArray;
188 if(fItTrackInputArray) delete fItTrackInputArray;
189 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
190 {
191 delete *itPhiBin;
192 }
193}
194
195//------------------------------------------------------------------------------
196
197void Calorimeter::Process()
198{
199 Candidate *particle, *track;
200 TLorentzVector position, momentum;
201 Short_t etaBin, phiBin, flags;
202 Int_t number;
203 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
204 Double_t ecalFraction, hcalFraction;
205 Double_t ecalEnergy, hcalEnergy;
206 Int_t pdgCode;
207
208 TFractionMap::iterator itFractionMap;
209
210 vector< Double_t >::iterator itEtaBin;
211 vector< Double_t >::iterator itPhiBin;
212 vector< Double_t > *phiBins;
213
214 vector< Long64_t >::iterator itTowerHits;
215
216 DelphesFactory *factory = GetFactory();
217 fTowerHits.clear();
218 fTowerECalFractions.clear();
219 fTowerHCalFractions.clear();
220 fTrackECalFractions.clear();
221 fTrackHCalFractions.clear();
222
223 // loop over all particles
224 fItParticleInputArray->Reset();
225 number = -1;
226 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
227 {
228 const TLorentzVector &particlePosition = particle->Position;
229 ++number;
230
231 pdgCode = TMath::Abs(particle->PID);
232
233 itFractionMap = fFractionMap.find(pdgCode);
234 if(itFractionMap == fFractionMap.end())
235 {
236 itFractionMap = fFractionMap.find(0);
237 }
238
239 ecalFraction = itFractionMap->second.first;
240 hcalFraction = itFractionMap->second.second;
241
242 fTowerECalFractions.push_back(ecalFraction);
243 fTowerHCalFractions.push_back(hcalFraction);
244
245 if(ecalFraction < 1.0E-9 && hcalFraction < 1.0E-9) continue;
246
247 // find eta bin [1, fEtaBins.size - 1]
248 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
249 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
250 etaBin = distance(fEtaBins.begin(), itEtaBin);
251
252 // phi bins for given eta bin
253 phiBins = fPhiBins[etaBin];
254
255 // find phi bin [1, phiBins.size - 1]
256 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
257 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
258 phiBin = distance(phiBins->begin(), itPhiBin);
259
260 flags = 0;
261 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
262
263 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
264 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
265
266 fTowerHits.push_back(towerHit);
267 }
268
269 // loop over all tracks
270 fItTrackInputArray->Reset();
271 number = -1;
272 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
273 {
274 const TLorentzVector &trackPosition = track->Position;
275 ++number;
276
277 pdgCode = TMath::Abs(track->PID);
278
279 itFractionMap = fFractionMap.find(pdgCode);
280 if(itFractionMap == fFractionMap.end())
281 {
282 itFractionMap = fFractionMap.find(0);
283 }
284
285 ecalFraction = itFractionMap->second.first;
286 hcalFraction = itFractionMap->second.second;
287
288 fTrackECalFractions.push_back(ecalFraction);
289 fTrackHCalFractions.push_back(hcalFraction);
290
291 // find eta bin [1, fEtaBins.size - 1]
292 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
293 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
294 etaBin = distance(fEtaBins.begin(), itEtaBin);
295
296 // phi bins for given eta bin
297 phiBins = fPhiBins[etaBin];
298
299 // find phi bin [1, phiBins.size - 1]
300 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
301 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
302 phiBin = distance(phiBins->begin(), itPhiBin);
303
304 flags = 1;
305
306 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
307 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
308
309 fTowerHits.push_back(towerHit);
310 }
311
312 // all hits are sorted first by eta bin number, then by phi bin number,
313 // then by flags and then by particle or track number
314 sort(fTowerHits.begin(), fTowerHits.end());
315
316 // loop over all hits
317 towerEtaPhi = 0;
318 fTower = 0;
319 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
320 {
321 towerHit = (*itTowerHits);
322 flags = (towerHit >> 24) & 0x00000000000000FFLL;
323 number = (towerHit) & 0x0000000000FFFFFFLL;
324 hitEtaPhi = towerHit >> 32;
325
326 if(towerEtaPhi != hitEtaPhi)
327 {
328 // switch to next tower
329 towerEtaPhi = hitEtaPhi;
330
331 // finalize previous tower
332 FinalizeTower();
333
334 // create new tower
335 fTower = factory->NewCandidate();
336
337 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
338 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
339
340 // phi bins for given eta bin
341 phiBins = fPhiBins[etaBin];
342
343 // calculate eta and phi of the tower's center
344 fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
345 fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
346
347 fTowerEdges[0] = fEtaBins[etaBin - 1];
348 fTowerEdges[1] = fEtaBins[etaBin];
349 fTowerEdges[2] = (*phiBins)[phiBin - 1];
350 fTowerEdges[3] = (*phiBins)[phiBin];
351
352 fTowerECalEnergy = 0.0;
353 fTowerHCalEnergy = 0.0;
354
355 fTrackECalEnergy = 0.0;
356 fTrackHCalEnergy = 0.0;
357
358 fTowerECalTime = 0.0;
359 fTowerHCalTime = 0.0;
360
361 fTrackECalTime = 0.0;
362 fTrackHCalTime = 0.0;
363
364 fTowerECalTimeWeight = 0.0;
365 fTowerHCalTimeWeight = 0.0;
366
367 fTowerTrackHits = 0;
368 fTowerPhotonHits = 0;
369
370 fTowerTrackArray->Clear();
371 }
372
373 // check for track hits
374 if(flags & 1)
375 {
376 ++fTowerTrackHits;
377
378 track = static_cast<Candidate*>(fTrackInputArray->At(number));
379 momentum = track->Momentum;
380 position = track->Position;
381
382
383 ecalEnergy = momentum.E() * fTrackECalFractions[number];
384 hcalEnergy = momentum.E() * fTrackHCalFractions[number];
385
386 fTrackECalEnergy += ecalEnergy;
387 fTrackHCalEnergy += hcalEnergy;
388
389 fTrackECalTime += TMath::Sqrt(ecalEnergy)*position.T();
390 fTrackHCalTime += TMath::Sqrt(hcalEnergy)*position.T();
391
392 fTrackECalTimeWeight += TMath::Sqrt(ecalEnergy);
393 fTrackHCalTimeWeight += TMath::Sqrt(hcalEnergy);
394
395 fTowerTrackArray->Add(track);
396
397 continue;
398 }
399
400 // check for photon and electron hits in current tower
401 if(flags & 2) ++fTowerPhotonHits;
402
403 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
404 momentum = particle->Momentum;
405 position = particle->Position;
406
407 // fill current tower
408 ecalEnergy = momentum.E() * fTowerECalFractions[number];
409 hcalEnergy = momentum.E() * fTowerHCalFractions[number];
410
411 fTowerECalEnergy += ecalEnergy;
412 fTowerHCalEnergy += hcalEnergy;
413
414 fTowerECalTime += TMath::Sqrt(ecalEnergy)*position.T();
415 fTowerHCalTime += TMath::Sqrt(hcalEnergy)*position.T();
416
417 fTowerECalTimeWeight += TMath::Sqrt(ecalEnergy);
418 fTowerHCalTimeWeight += TMath::Sqrt(hcalEnergy);
419
420
421 fTower->AddCandidate(particle);
422 }
423
424 // finalize last tower
425 FinalizeTower();
426}
427
428//------------------------------------------------------------------------------
429
430void Calorimeter::FinalizeTower()
431{
432 Candidate *track, *tower;
433 Double_t energy, pt, eta, phi;
434 Double_t ecalEnergy, hcalEnergy;
435 Double_t ecalSigma, hcalSigma;
436 Double_t ecalTime, hcalTime, time;
437
438 if(!fTower) return;
439
440 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerECalEnergy);
441 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerHCalEnergy);
442
443 ecalEnergy = LogNormal(fTowerECalEnergy, ecalSigma);
444 hcalEnergy = LogNormal(fTowerHCalEnergy, hcalSigma);
445
446 ecalTime = (fTowerECalTimeWeight < 1.0E-09 ) ? 0.0 : fTowerECalTime/fTowerECalTimeWeight;
447 hcalTime = (fTowerHCalTimeWeight < 1.0E-09 ) ? 0.0 : fTowerHCalTime/fTowerHCalTimeWeight;
448
449 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
450 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
451
452 if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
453 if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
454
455 energy = ecalEnergy + hcalEnergy;
456 time = (TMath::Sqrt(ecalEnergy)*ecalTime + TMath::Sqrt(hcalEnergy)*hcalTime)/(TMath::Sqrt(ecalEnergy) + TMath::Sqrt(hcalEnergy));
457
458 if(fSmearTowerCenter)
459 {
460 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
461 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
462 }
463 else
464 {
465 eta = fTowerEta;
466 phi = fTowerPhi;
467 }
468
469 pt = energy / TMath::CosH(eta);
470
471 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, time);
472 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
473 fTower->Eem = ecalEnergy;
474 fTower->Ehad = hcalEnergy;
475
476 fTower->Edges[0] = fTowerEdges[0];
477 fTower->Edges[1] = fTowerEdges[1];
478 fTower->Edges[2] = fTowerEdges[2];
479 fTower->Edges[3] = fTowerEdges[3];
480
481 if(energy > 0.0)
482 {
483 if(fTowerPhotonHits > 0 && fTowerTrackHits == 0)
484 {
485 fPhotonOutputArray->Add(fTower);
486 }
487
488 fTowerOutputArray->Add(fTower);
489 }
490
491 // fill energy flow candidates
492
493 // save all the tracks as energy flow tracks
494 fItTowerTrackArray->Reset();
495 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
496 {
497 fEFlowTrackOutputArray->Add(track);
498 }
499
500 ecalEnergy -= fTrackECalEnergy;
501 hcalEnergy -= fTrackHCalEnergy;
502
503 ecalSigma = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, ecalEnergy);
504 hcalSigma = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, hcalEnergy);
505
506 if(ecalEnergy < fECalEnergyMin || ecalEnergy < fECalEnergySignificanceMin*ecalSigma) ecalEnergy = 0.0;
507 if(hcalEnergy < fHCalEnergyMin || hcalEnergy < fHCalEnergySignificanceMin*hcalSigma) hcalEnergy = 0.0;
508
509 energy = ecalEnergy + hcalEnergy;
510
511 if(ecalEnergy > 0.0)
512 {
513 // create new photon tower
514 tower = static_cast<Candidate*>(fTower->Clone());
515
516 pt = ecalEnergy / TMath::CosH(eta);
517
518 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalEnergy);
519 tower->Eem = ecalEnergy;
520 tower->Ehad = 0.0;
521
522 fEFlowPhotonOutputArray->Add(tower);
523 }
524 if(hcalEnergy > 0.0)
525 {
526 // create new neutral hadron tower
527 tower = static_cast<Candidate*>(fTower->Clone());
528
529 pt = hcalEnergy / TMath::CosH(eta);
530
531 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, hcalEnergy);
532 tower->Eem = 0.0;
533 tower->Ehad = hcalEnergy;
534
535 fEFlowNeutralHadronOutputArray->Add(tower);
536 }
537}
538
539//------------------------------------------------------------------------------
540
541Double_t Calorimeter::LogNormal(Double_t mean, Double_t sigma)
542{
543 Double_t a, b;
544
545 if(mean > 0.0)
546 {
547 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
548 a = TMath::Log(mean) - 0.5*b*b;
549
550 return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
551 }
552 else
553 {
554 return 0.0;
555 }
556}
Note: See TracBrowser for help on using the repository browser.