Fork me on GitHub

source: svn/trunk/modules/Calorimeter.cc@ 1142

Last change on this file since 1142 was 1142, checked in by Pavel Demin, 11 years ago

use log normal for calorimeter resolution

  • Property svn:keywords set to Id Revision Date
File size: 14.8 KB
Line 
1
2/** \class Calorimeter
3 *
4 * Fills calorimeter towers, performs calorimeter resolution smearing,
5 * preselects towers hit by photons and creates energy flow objects.
6 *
7 * $Date: 2013-06-25 12:43:27 +0000 (Tue, 25 Jun 2013) $
8 * $Revision: 1142 $
9 *
10 *
11 * \author P. Demin - UCL, Louvain-la-Neuve
12 *
13 */
14
15#include "modules/Calorimeter.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
42Calorimeter::Calorimeter() :
43 fECalResolutionFormula(0), fHCalResolutionFormula(0),
44 fItParticleInputArray(0), fItTrackInputArray(0),
45 fTowerTrackArray(0), fItTowerTrackArray(0),
46 fTowerPhotonArray(0), fItTowerPhotonArray(0)
47{
48 fECalResolutionFormula = new DelphesFormula;
49 fHCalResolutionFormula = new DelphesFormula;
50 fTowerTrackArray = new TObjArray;
51 fItTowerTrackArray = fTowerTrackArray->MakeIterator();
52 fTowerPhotonArray = new TObjArray;
53 fItTowerPhotonArray = fTowerPhotonArray->MakeIterator();
54}
55
56//------------------------------------------------------------------------------
57
58Calorimeter::~Calorimeter()
59{
60 if(fECalResolutionFormula) delete fECalResolutionFormula;
61 if(fHCalResolutionFormula) delete fHCalResolutionFormula;
62 if(fTowerTrackArray) delete fTowerTrackArray;
63 if(fItTowerTrackArray) delete fItTowerTrackArray;
64 if(fTowerPhotonArray) delete fTowerPhotonArray;
65 if(fItTowerPhotonArray) delete fItTowerPhotonArray;
66}
67
68//------------------------------------------------------------------------------
69
70void Calorimeter::Init()
71{
72 ExRootConfParam param, paramEtaBins, paramPhiBins, paramFractions;
73 Long_t i, j, k, size, sizeEtaBins, sizePhiBins, sizeFractions;
74 Double_t ecalFraction, hcalFraction;
75 TBinMap::iterator itEtaBin;
76 set< Double_t >::iterator itPhiBin;
77 vector< Double_t > *phiBins;
78
79 // read eta and phi bins
80 param = GetParam("EtaPhiBins");
81 size = param.GetSize();
82 fBinMap.clear();
83 fEtaBins.clear();
84 fPhiBins.clear();
85 for(i = 0; i < size/2; ++i)
86 {
87 paramEtaBins = param[i*2];
88 sizeEtaBins = paramEtaBins.GetSize();
89 paramPhiBins = param[i*2 + 1];
90 sizePhiBins = paramPhiBins.GetSize();
91
92 for(j = 0; j < sizeEtaBins; ++j)
93 {
94 for(k = 0; k < sizePhiBins; ++k)
95 {
96 fBinMap[paramEtaBins[j].GetDouble()].insert(paramPhiBins[k].GetDouble());
97 }
98 }
99 }
100
101 // for better performance we transform map of sets to parallel vectors:
102 // vector< double > and vector< vector< double >* >
103 for(itEtaBin = fBinMap.begin(); itEtaBin != fBinMap.end(); ++itEtaBin)
104 {
105 fEtaBins.push_back(itEtaBin->first);
106 phiBins = new vector< double >(itEtaBin->second.size());
107 fPhiBins.push_back(phiBins);
108 phiBins->clear();
109 for(itPhiBin = itEtaBin->second.begin(); itPhiBin != itEtaBin->second.end(); ++itPhiBin)
110 {
111 phiBins->push_back(*itPhiBin);
112 }
113 }
114
115 // read energy fractions for different particles
116 param = GetParam("EnergyFraction");
117 size = param.GetSize();
118
119 // set default energy fractions values
120 fFractionMap.clear();
121 fFractionMap[0] = make_pair(0.0, 1.0);
122
123 for(i = 0; i < size/2; ++i)
124 {
125 paramFractions = param[i*2 + 1];
126 sizeFractions = paramFractions.GetSize();
127
128 ecalFraction = paramFractions[0].GetDouble();
129 hcalFraction = paramFractions[1].GetDouble();
130
131 fFractionMap[param[i*2].GetInt()] = make_pair(ecalFraction, hcalFraction);
132 }
133/*
134 TFractionMap::iterator itFractionMap;
135 for(itFractionMap = fFractionMap.begin(); itFractionMap != fFractionMap.end(); ++itFractionMap)
136 {
137 cout << itFractionMap->first << " " << itFractionMap->second.first << " " << itFractionMap->second.second << endl;
138 }
139*/
140 // read resolution formulas
141 fECalResolutionFormula->Compile(GetString("ECalResolutionFormula", "0"));
142 fHCalResolutionFormula->Compile(GetString("HCalResolutionFormula", "0"));
143
144 // import array with output from other modules
145 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "ParticlePropagator/particles"));
146 fItParticleInputArray = fParticleInputArray->MakeIterator();
147
148 fTrackInputArray = ImportArray(GetString("TrackInputArray", "ParticlePropagator/tracks"));
149 fItTrackInputArray = fTrackInputArray->MakeIterator();
150
151 // create output arrays
152 fTowerOutputArray = ExportArray(GetString("TowerOutputArray", "towers"));
153 fPhotonOutputArray = ExportArray(GetString("PhotonOutputArray", "photons"));
154
155 fEFlowTrackOutputArray = ExportArray(GetString("EFlowTrackOutputArray", "eflowTracks"));
156 fEFlowTowerOutputArray = ExportArray(GetString("EFlowTowerOutputArray", "eflowTowers"));
157
158
159
160}
161
162//------------------------------------------------------------------------------
163
164void Calorimeter::Finish()
165{
166 vector< vector< Double_t>* >::iterator itPhiBin;
167 if(fItParticleInputArray) delete fItParticleInputArray;
168 if(fItTrackInputArray) delete fItTrackInputArray;
169 for(itPhiBin = fPhiBins.begin(); itPhiBin != fPhiBins.end(); ++itPhiBin)
170 {
171 delete *itPhiBin;
172 }
173}
174
175//------------------------------------------------------------------------------
176
177void Calorimeter::Process()
178{
179 Candidate *particle, *track;
180 TLorentzVector position, momentum;
181 Short_t etaBin, phiBin, flags;
182 Int_t number;
183 Long64_t towerHit, towerEtaPhi, hitEtaPhi;
184 Double_t ecalFraction, hcalFraction;
185 Double_t ecalEnergy, hcalEnergy;
186 Int_t pdgCode;
187
188 TFractionMap::iterator itFractionMap;
189
190 vector< Double_t >::iterator itEtaBin;
191 vector< Double_t >::iterator itPhiBin;
192 vector< Double_t > *phiBins;
193
194 vector< Long64_t >::iterator itTowerHits;
195
196 DelphesFactory *factory = GetFactory();
197 fTowerHits.clear();
198 fECalFractions.clear();
199 fHCalFractions.clear();
200
201 // loop over all particles
202 fItParticleInputArray->Reset();
203 number = -1;
204 while((particle = static_cast<Candidate*>(fItParticleInputArray->Next())))
205 {
206 const TLorentzVector &particlePosition = particle->Position;
207 ++number;
208
209 pdgCode = TMath::Abs(particle->PID);
210
211 itFractionMap = fFractionMap.find(pdgCode);
212 if(itFractionMap == fFractionMap.end())
213 {
214 itFractionMap = fFractionMap.find(0);
215 }
216
217 ecalFraction = itFractionMap->second.first;
218 hcalFraction = itFractionMap->second.second;
219
220 fECalFractions.push_back(ecalFraction);
221 fHCalFractions.push_back(hcalFraction);
222
223 if(ecalFraction < 1.0E-9 && hcalFraction < 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 = (particle->Charge == 0);
239 flags |= (pdgCode == 22) << 1;
240 flags |= (pdgCode == 11) << 2;
241
242 // make tower hit {16-bits for eta bin number, 16-bits for phi bin number, 8-bits for flags, 24-bits for particle number}
243 towerHit = (Long64_t(etaBin) << 48) | (Long64_t(phiBin) << 32) | (Long64_t(flags) << 24) | Long64_t(number);
244
245 fTowerHits.push_back(towerHit);
246 }
247
248 // loop over all tracks
249 fItTrackInputArray->Reset();
250 number = -1;
251 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
252 {
253 const TLorentzVector &trackPosition = track->Position;
254 ++number;
255
256 // find eta bin [1, fEtaBins.size - 1]
257 itEtaBin = lower_bound(fEtaBins.begin(), fEtaBins.end(), trackPosition.Eta());
258 if(itEtaBin == fEtaBins.begin() || itEtaBin == fEtaBins.end()) continue;
259 etaBin = distance(fEtaBins.begin(), itEtaBin);
260
261 // phi bins for given eta bin
262 phiBins = fPhiBins[etaBin];
263
264 // find phi bin [1, phiBins.size - 1]
265 itPhiBin = lower_bound(phiBins->begin(), phiBins->end(), trackPosition.Phi());
266 if(itPhiBin == phiBins->begin() || itPhiBin == phiBins->end()) continue;
267 phiBin = distance(phiBins->begin(), itPhiBin);
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(1) << 27) | 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 fTowerECalEnergy = 0.0;
316 fTowerHCalEnergy = 0.0;
317
318 fTowerECalNeutralEnergy = 0.0;
319 fTowerHCalNeutralEnergy = 0.0;
320
321 fTowerNeutralHits = 0;
322 fTowerPhotonHits = 0;
323 fTowerElectronHits = 0;
324 fTowerTrackHits = 0;
325 fTowerAllHits = 0;
326
327 fTowerTrackArray->Clear();
328 fTowerPhotonArray->Clear();
329 }
330
331 // check for track hits
332 if(flags & 8)
333 {
334 ++fTowerTrackHits;
335 track = static_cast<Candidate*>(fTrackInputArray->At(number));
336 fTowerTrackArray->Add(track);
337 continue;
338 }
339
340 particle = static_cast<Candidate*>(fParticleInputArray->At(number));
341 momentum = particle->Momentum;
342
343 // fill current tower
344 ecalEnergy = momentum.E() * fECalFractions[number];
345 hcalEnergy = momentum.E() * fHCalFractions[number];
346
347 fTowerECalEnergy += ecalEnergy;
348 fTowerHCalEnergy += hcalEnergy;
349
350 ++fTowerAllHits;
351 fTower->AddCandidate(particle);
352
353 // check for neutral hits in current tower
354 if(flags & 1) ++fTowerNeutralHits;
355
356 // check for photon hits in current tower
357 if(flags & 2)
358 {
359 ++fTowerPhotonHits;
360 fTowerPhotonArray->Add(particle);
361 }
362
363 // check for electron hits in current tower
364 if(flags & 4) ++fTowerElectronHits;
365 }
366
367 // finalize last tower
368 FinalizeTower();
369}
370
371//------------------------------------------------------------------------------
372
373void Calorimeter::FinalizeTower()
374{
375 Candidate *particle, *track, *tower;
376 Double_t energy, pt, eta, phi;
377 Double_t ecalEnergy, hcalEnergy;
378 Double_t sE,mE,sH,mH,SE=0,ME=0,SH=0,MH=0;
379
380 if(!fTower) return;
381
382 // cout<<"new tower -----------------------"<<endl;
383
384
385 sE = fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerECalEnergy);
386 mE = fTowerECalEnergy;
387
388 sH = fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerHCalEnergy);
389 mH = fTowerHCalEnergy;
390
391 // cout<<sE<<","<<mE<<endl;
392 // cout<<sH<<","<<mH<<endl;
393 // cout<<"---"<<endl;
394
395 if(mE>0)SE = TMath::Sqrt(TMath::Log((1+(sE*sE)/(mE*mE))));
396 if(mE>0)ME = TMath::Log(mE)-0.5*SE*SE;
397
398 if(mH>0)SH = TMath::Sqrt(TMath::Log((1+(sH*sH)/(mH*mH))));
399 if(mH>0)MH = TMath::Log(mH)-0.5*SH*SH;
400
401 // cout<<SE<<","<<ME<<endl;
402 // cout<<SH<<","<<MH<<endl;
403 // cout<<"---"<<endl;
404
405// ecalEnergy = gRandom->Gaus(fTowerECalEnergy, fECalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerECalEnergy));
406 ecalEnergy = (ME != 0.0 ? Draw_lognormal(ME,SE) : 0.0);
407 // cout<<ecalEnergy<<endl;
408
409 if(ecalEnergy < 0.0) ecalEnergy = 0.0;
410
411 hcalEnergy = (MH != 0.0 ? Draw_lognormal(MH,SH) : 0.0);
412 //hcalEnergy = gRandom->Gaus(fTowerHCalEnergy, fHCalResolutionFormula->Eval(0.0, fTowerEta, 0.0, fTowerHCalEnergy));
413 if(hcalEnergy < 0.0) hcalEnergy = 0.0;
414
415// cout<<hcalEnergy<<endl;
416
417 energy = ecalEnergy + hcalEnergy;
418
419 // eta = fTowerEta;
420 // phi = fTowerPhi;
421
422 eta = gRandom->Uniform(fTowerEdges[0], fTowerEdges[1]);
423 phi = gRandom->Uniform(fTowerEdges[2], fTowerEdges[3]);
424
425 pt = energy / TMath::CosH(eta);
426
427 fTower->Position.SetPtEtaPhiE(1.0, eta, phi, 0.0);
428 fTower->Momentum.SetPtEtaPhiE(pt, eta, phi, energy);
429 fTower->Eem = ecalEnergy;
430 fTower->Ehad = hcalEnergy;
431
432 fTower->Edges[0] = fTowerEdges[0];
433 fTower->Edges[1] = fTowerEdges[1];
434 fTower->Edges[2] = fTowerEdges[2];
435 fTower->Edges[3] = fTowerEdges[3];
436
437 // fill calorimeter towers and photon candidates
438 if(energy > 0.0)
439 {
440 if((fTowerPhotonHits > 0 || fTowerElectronHits > 0) &&
441 fTowerTrackHits == 0)
442 {
443 fPhotonOutputArray->Add(fTower);
444 }
445
446 fTowerOutputArray->Add(fTower);
447 }
448
449 // fill energy flow candidates
450 if(fTowerTrackHits == fTowerAllHits)
451 {
452 fItTowerTrackArray->Reset();
453 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
454 {
455 fEFlowTrackOutputArray->Add(track);
456 }
457 }
458 else if(fTowerTrackHits > 0 &&
459 fTowerElectronHits == 0 &&
460 fTowerPhotonHits + fTowerTrackHits == fTowerAllHits)
461 {
462 fItTowerTrackArray->Reset();
463 while((track = static_cast<Candidate*>(fItTowerTrackArray->Next())))
464 {
465 fEFlowTrackOutputArray->Add(track);
466 }
467
468 if(ecalEnergy > 0.0)
469 {
470 DelphesFactory *factory = GetFactory();
471
472 // create new tower
473 tower = factory->NewCandidate();
474
475 fItTowerPhotonArray->Reset();
476 while((particle = static_cast<Candidate*>(fItTowerPhotonArray->Next())))
477 {
478 tower->AddCandidate(particle);
479 }
480
481 pt = ecalEnergy / TMath::CosH(eta);
482
483 tower->Position.SetPtEtaPhiE(1.0, eta, phi, 0.0);
484 tower->Momentum.SetPtEtaPhiE(pt, eta, phi, ecalEnergy);
485 tower->Eem = ecalEnergy;
486 tower->Ehad = 0.0;
487
488 tower->Edges[0] = fTowerEdges[0];
489 tower->Edges[1] = fTowerEdges[1];
490 tower->Edges[2] = fTowerEdges[2];
491 tower->Edges[3] = fTowerEdges[3];
492
493 fEFlowTowerOutputArray->Add(tower);
494 }
495 }
496 else if(energy > 0.0)
497 {
498 fEFlowTowerOutputArray->Add(fTower);
499 }
500}
501
502//------------------------------------------------------------------------------
503
504Double_t Calorimeter::Draw_lognormal(Double_t mu, Double_t sigma)
505 {
506 Double_t g = gRandom->Gaus(0, 1);
507 return TMath::Exp(mu + sigma * g);
508 }
Note: See TracBrowser for help on using the repository browser.