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