Fork me on GitHub

source: git/modules/SimpleCalorimeter.cc@ f6b6ee7

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

fix EOL characters in GPLv3 header

  • Property mode set to 100644
File size: 13.1 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 SimpleCalorimeter
21 *
22 * Fills SimpleCalorimeter towers, performs SimpleCalorimeter resolution smearing,
23 * and creates energy flow objects (tracks, photons, and neutral hadrons).
24 *
25 * $Date: 2014-04-16 15:29:31 +0200 (Wed, 16 Apr 2014) $
26 * $Revision: 1364 $
27 *
28 *
29 * \author P. Demin - UCL, Louvain-la-Neuve
30 *
31 */
32
33#include "modules/SimpleCalorimeter.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
60SimpleCalorimeter::SimpleCalorimeter() :
61 fResolutionFormula(0),
62 fItParticleInputArray(0), fItTrackInputArray(0),
63 fTowerTrackArray(0), fItTowerTrackArray(0)
64{
65 fResolutionFormula = new DelphesFormula;
66
67 fTowerTrackArray = new TObjArray;
68 fItTowerTrackArray = fTowerTrackArray->MakeIterator();
69}
70
71//------------------------------------------------------------------------------
72
73SimpleCalorimeter::~SimpleCalorimeter()
74{
75 if(fResolutionFormula) delete fResolutionFormula;
76
77 if(fTowerTrackArray) delete fTowerTrackArray;
78 if(fItTowerTrackArray) delete fItTowerTrackArray;
79}
80
81//------------------------------------------------------------------------------
82
83void SimpleCalorimeter::Init()
84{
85 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
86 Long_t i, j, k, size, sizeEtaBins, sizePhiBins;
87 Double_t fraction;
88 TBinMap::iterator itEtaBin;
89 set< Double_t >::iterator itPhiBin;
90 vector< Double_t > *phiBins;
91
92 // read eta and phi bins
93 param = GetParam("EtaPhiBins");
94 size = param.GetSize();
95 fBinMap.clear();
96 fEtaBins.clear();
97 fPhiBins.clear();
98 for(i = 0; i < size/2; ++i)
99 {
100 paramEtaBins = param[i*2];
101 sizeEtaBins = paramEtaBins.GetSize();
102 paramPhiBins = param[i*2 + 1];
103 sizePhiBins = paramPhiBins.GetSize();
104
105 for(j = 0; j < sizeEtaBins; ++j)
106 {
107 for(k = 0; k < sizePhiBins; ++k)
108 {
109 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
110 }
111 }
112 }
113
114 // for better performance we transform map of sets to parallel vectors:
115 // vector< double > and vector< vector< double >* >
116 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
117 {
118 fEtaBins.push_back(itEtaBin->first);
119 phiBins = new vector< double >(itEtaBin->second.size());
120 fPhiBins.push_back(phiBins);
121 phiBins->clear();
122 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
123 {
124 phiBins->push_back(*itPhiBin);
125 }
126 }
127
128 // read energy fractions for different particles
129 param = GetParam("EnergyFraction");
130 size = param.GetSize();
131
132 // set default energy fractions values
133 fFractionMap.clear();
134 fFractionMap[0] = 1.0;
135
136 for(i = 0; i < size/2; ++i)
137 {
138 paramFractions = param[i*2 + 1];
139 fraction = paramFractions[0].GetDouble();
140 fFractionMap[param[i*2].GetInt()] = fraction;
141 }
142/*
143 TFractionMap::iterator itFractionMap;
144 for(itFractionMap = fFractionMap.begin(); itFractionMap != fFractionMap.end(); ++itFractionMap)
145 {
146 cout << itFractionMap->first << " " << itFractionMap->second.first << " " << itFractionMap->second.second << endl;
147 }
148*/
149 // read resolution formulas
150 fResolutionFormula->Compile(GetString("ResolutionFormula", "0"));
151
152 // import array with output from other modules
153 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
154 fItParticleInputArray = fParticleInputArray->MakeIterator();
155
156 fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
157 fItTrackInputArray = fTrackInputArray->MakeIterator();
158
159 // create output arrays
160 fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
161 fEFlowTowerOutputArray = ExportArray(GetString("EFlowTowerOutputArray", "eflowTowers"));
162
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 Int_t pdgCode;
190
191 TFractionMap::iterator itFractionMap;
192
193 vector< Double_t >::iterator itEtaBin;
194 vector< Double_t >::iterator itPhiBin;
195 vector< Double_t > *phiBins;
196
197 vector< Long64_t >::iterator itTowerHits;
198
199 DelphesFactory *factory = GetFactory();
200 fTowerHits.clear();
201 fTowerFractions.clear();
202 fTrackFractions.clear();
203
204 // loop over all particles
205 fItParticleInputArray->Reset();
206 number = -1;
207 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
208 {
209 const TLorentzVector &particlePosition = particle->Position;
210 ++number;
211
212 pdgCode = TMath::Abs(particle->PID);
213
214 itFractionMap = fFractionMap.find(pdgCode);
215 if(itFractionMap == fFractionMap.end())
216 {
217 itFractionMap = fFractionMap.find(0);
218 }
219
220 fraction = itFractionMap->second;
221 fTowerFractions.push_back(fraction);
222
223 if(fraction < 1.0E-9) continue;
224
225 // find eta bin [1, fEtaBins.size - 1]
226 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), particlePosition.Eta());
227 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
228 etaBin = distance(fEtaBins.begin(), itEtaBin);
229
230 // phi bins for given eta bin
231 phiBins = fPhiBins[etaBin];
232
233 // find phi bin [1, phiBins.size - 1]
234 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), particlePosition.Phi());
235 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
236 phiBin = distance(phiBins->begin(), itPhiBin);
237
238 flags = 0;
239 flags |= (pdgCode == 11 || pdgCode == 22) << 1;
240
241 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
242 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
243
244 fTowerHits.push_back(towerHit);
245 }
246
247 // loop over all tracks
248 fItTrackInputArray->Reset();
249 number = -1;
250 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
251 {
252 const TLorentzVector &trackPosition = track->Position;
253 ++number;
254
255 pdgCode = TMath::Abs(track->PID);
256
257 itFractionMap = fFractionMap.find(pdgCode);
258 if(itFractionMap == fFractionMap.end())
259 {
260 itFractionMap = fFractionMap.find(0);
261 }
262
263 fraction = itFractionMap->second;
264
265 fTrackFractions.push_back(fraction);
266
267 // find eta bin [1, fEtaBins.size - 1]
268 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
269 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
270 etaBin = distance(fEtaBins.begin(), itEtaBin);
271
272 // phi bins for given eta bin
273 phiBins = fPhiBins[etaBin];
274
275 // find phi bin [1, phiBins.size - 1]
276 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
277 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
278 phiBin = distance(phiBins->begin(), itPhiBin);
279
280 flags = 1;
281
282 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for track number}
283 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
284
285 fTowerHits.push_back(towerHit);
286 }
287
288 // all hits are sorted first by eta bin number, then by phi bin number,
289 // then by flags and then by particle or track number
290 sort(fTowerHits.begin(), fTowerHits.end());
291
292 // loop over all hits
293 towerEtaPhi = 0;
294 fTower = 0;
295 for(itTowerHits = fTowerHits.begin(); itTowerHits != fTowerHits.end(); ++itTowerHits)
296 {
297 towerHit = (*itTowerHits);
298 flags = (towerHit >> 24) & 0x00000000000000FFLL;
299 number = (towerHit) & 0x0000000000FFFFFFLL;
300 hitEtaPhi = towerHit >> 32;
301
302 if(towerEtaPhi != hitEtaPhi)
303 {
304 // switch to next tower
305 towerEtaPhi = hitEtaPhi;
306
307 // finalize previous tower
308 FinalizeTower();
309
310 // create new tower
311 fTower = factory->NewCandidate();
312
313 phiBin = (towerHit >> 32) & 0x000000000000FFFFLL;
314 etaBin = (towerHit >> 48) & 0x000000000000FFFFLL;
315
316 // phi bins for given eta bin
317 phiBins = fPhiBins[etaBin];
318
319 // calculate eta and phi of the tower's center
320 fTowerEta = 0.5*(fEtaBins[etaBin - 1] + fEtaBins[etaBin]);
321 fTowerPhi = 0.5*((*phiBins)[phiBin - 1] + (*phiBins)[phiBin]);
322
323 fTowerEdges[0] = fEtaBins[etaBin - 1];
324 fTowerEdges[1] = fEtaBins[etaBin];
325 fTowerEdges[2] = (*phiBins)[phiBin - 1];
326 fTowerEdges[3] = (*phiBins)[phiBin];
327
328 fTowerEnergy = 0.0;
329 fTrackEnergy = 0.0;
330
331 fTowerTime = 0.0;
332 fTrackTime = 0.0;
333
334 fTowerWeightTime = 0.0;
335
336 fTowerTrackHits = 0;
337 fTowerPhotonHits = 0;
338
339 fTowerTrackArray->Clear();
340 }
341
342 // check for track hits
343 if(flags & 1)
344 {
345 ++fTowerTrackHits;
346
347 track = static_cast<Candidate*>(fTrackInputArray->At(number));
348 momentum = track->Momentum;
349 position = track->Position;
350
351 energy = momentum.E() * fTrackFractions[number];
352
353 fTrackEnergy += energy;
354
355 fTrackTime += TMath::Sqrt(energy)*position.T();
356 fTrackWeightTime += TMath::Sqrt(energy);
357
358 fTowerTrackArray->Add(track);
359
360 continue;
361 }
362
363 // check for photon and electron hits in current tower
364 if(flags & 2) ++fTowerPhotonHits;
365
366 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
367 momentum = particle->Momentum;
368 position = particle->Position;
369
370 // fill current tower
371 energy = momentum.E() * fTowerFractions[number];
372
373 fTowerEnergy += energy;
374
375 fTowerTime += TMath::Sqrt(energy)*position.T();
376 fTowerWeightTime += TMath::Sqrt(energy);
377
378 fTower->AddCandidate(particle);
379 }
380
381 // finalize last tower
382 FinalizeTower();
383}
384
385//------------------------------------------------------------------------------
386
387void SimpleCalorimeter::FinalizeTower()
388{
389 Candidate *tower;
390 Double_t energy, pt, eta, phi;
391 Double_t sigma;
392 Double_t time;
393
394 if(!fTower) return;
395
396 sigma = fResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerEnergy);
397
398// energy = gRandom->Gaus(fTowerEnergy, sigma);
399// if(energy < 0.0) energy = 0.0;
400
401 energy = LogNormal(fTowerEnergy, sigma);
402 time = (fTowerWeightTime < 1.0E-09 ) ? 0 : fTowerTime/fTowerWeightTime;
403
404 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
405 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
406
407 pt = energy / TMath::CosH(eta);
408
409 // fTower->Position.SetXYZT(-time, 0.0, 0.0, time);
410 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, time);
411 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
412
413 fTower->Edges[0] = fTowerEdges[0];
414 fTower->Edges[1] = fTowerEdges[1];
415 fTower->Edges[2] = fTowerEdges[2];
416 fTower->Edges[3] = fTowerEdges[3];
417
418
419 // fill SimpleCalorimeter towers
420 if(energy > 0.0) fTowerOutputArray->Add(fTower);
421
422
423 // fill energy flow candidates
424 energy -= fTrackEnergy;
425 if(energy < 0.0) energy = 0.0;
426
427 // save energy excess as an energy flow tower
428 if(energy > 0.0)
429 {
430 // create new photon tower
431 tower = static_cast<Candidate*>(fTower->Clone());
432 pt = energy / TMath::CosH(eta);
433
434 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
435 fEFlowTowerOutputArray->Add(tower);
436 }
437
438}
439
440//------------------------------------------------------------------------------
441
442Double_t SimpleCalorimeter::LogNormal(Double_t mean, Double_t sigma)
443{
444 Double_t a, b;
445
446 if(mean > 0.0)
447 {
448 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
449 a = TMath::Log(mean) - 0.5*b*b;
450
451 return TMath::Exp(a + b*gRandom->Gaus(0, 1));
452 }
453 else
454 {
455 return 0.0;
456 }
457}
Note: See TracBrowser for help on using the repository browser.