Fork me on GitHub

source: git/modules/Calorimeter.cc@ f3fcb6e

ImprovedOutputFile Timing dual_readout llp
Last change on this file since f3fcb6e was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

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