Fork me on GitHub

source: git/modules/SimpleCalorimeter.cc@ 4b9a2dc

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 4b9a2dc was 4b9a2dc, checked in by Michele <michele.selvaggi@…>, 10 years ago

added EnergyMin and SigmaMin parameters

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