ATLAS Offline Software
Loading...
Searching...
No Matches
FCALDistEnergyShowerLib.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6// this header file
8
11
12#include <sstream>
13#include <fstream>
14
15#include <iostream>
16#include <iomanip>
17
18// G4 includes
19#include "G4Track.hh"
20
23
24#include "TTree.h"
25#include "TFile.h"
26#include "TParameter.h"
27
28#define LIB_VERSION 4
29
30namespace ShowerLib {
31
32 inline int calcKey(float dist) {
33 return (int)(dist*1000);
34 }
35
36 float FCALDistEnergyShowerLib::distance (double x, double y) const
37 {
38 double stepy = m_step * sqrt(3.)/2.;
39 double dy = fmod(y-m_yrodcent,stepy);
40 int ny = int( (y-m_yrodcent)/stepy );
41 //additional half shift:
42 double dx0 = (ny % 2 )*m_step/2.;
43 double dx = fmod(x-m_xrodcent-dx0,m_step);
44 double d1 = sqrt(dx*dx+dy*dy);
45 double d2 = sqrt((m_step-dx)*(m_step-dx)+dy*dy);
46 double d3 = sqrt((m_step/2.-dx)*(m_step/2.-dx)+(stepy-dy)*(stepy-dy));
47 float dd = d1;
48 if (dd > d2) dd = d2;
49 if (dd > d3) dd = d3;
50 return dd;
51 }
52
54 {
55 TParameter<int>* ver;
56 ver = (TParameter<int>*)source->Get("version");
57
58 if ((ver == nullptr) || (ver->GetVal() != LIB_VERSION)) return nullptr;
59
60 TTree* TTreeMeta = (TTree*)source->Get("meta");
61 TTree* TTreeLib = (TTree*)source->Get("library");
62
63 if ((TTreeMeta == nullptr) || (TTreeLib == nullptr)) return nullptr;
64
65 std::cout << "FCALDistEnergyShowerLib header found." << std::endl;
66
68
69 if (!(newlib->readMeta(TTreeMeta)) || !(newlib->read(TTreeLib))) {
70 delete newlib;
71 std::cout << "FCALDistEnergyShowerLib read unsuccessful." << std::endl;
72 return nullptr;
73 }
74
75 return newlib;
76
77 }
78
80 {
81 /*
82 * Eta Energy Structure format:
83 *
84 * VER PART DET DIST1 DIST2 ...
85 * XRODCENT YRODCENT STEP
86 * COMMENT
87 *
88 * where
89 *
90 * VER == 4
91 * DIST(N+1) > DIST(N)
92 * DIST(N) >= 0
93 * First DIST is the beginning of the lib. If not zero, zero will be added
94 * Last DIST is the end of the lib. If less then 4.5, 4.5 will be added
95 */
96 std::ifstream filestr(inputFile.c_str(),std::ios::in);
97
98
99 if (!filestr.is_open()) {
100 std::cout << "FCALDistEnergyShowerLib " << inputFile << ": bad file!" << std::endl;
101 return nullptr;
102 }
103
104 std::string instr;
105 std::getline(filestr,instr);
106 std::stringstream ss(instr);
107
108 int ver;
109
110 ss >> ver;
111
112 if (ver != LIB_VERSION) {
113 return nullptr;
114 }
115
116
117 int part;
118 std::string det;
119
120 ss >> part >> det;
121
122 std::vector<float> etalist;
123
124 float eta;
125 float etaold;
126
127 ss >> eta;
128
129 if (ss.fail()) {
130 std::cout << "file reading failed! (not enough data)" << std::endl;
131 return nullptr;
132 }
133
134 if (eta != 0.) {
135 if (eta < 0) {
136 std::cout << "no negative dists allowed (" << eta << ")" << std::endl;
137 return nullptr;
138 }
139 etalist.push_back(0.);
140 }
141
142 etalist.push_back(eta);
143 etaold = eta;
144
145 while (!ss.eof()) {
146 ss >> eta;
147 if ((ss.fail()) || (eta <= etaold)) {
148 std::cout << "screwed dists! (" << eta << "<=" << etaold << ")" << std::endl;
149 return nullptr;
150 }
151 if (eta < 4.5) {
152 etalist.push_back(eta);
153 } else if (eta > 4.5) {
154 std::cout << "dist can't be bigger than 4.5 (" << eta << "). ignored" << std::endl;
155 }
156 etaold = eta;
157 }
158
159 //if (etalist.back() < 4.5) etalist.push_back(4.5);
160
162 if (!newlib->readStructure(etalist)) {
163 std::cout << "this structure is not valid" << std::endl;
164 delete newlib;
165 return nullptr;
166 }
167
168 newlib->m_detector = std::move(det);
169 newlib->m_particle = part;
170
171 newlib->m_filled = false;
172
173 std::getline(filestr,instr);
174 std::stringstream ss1(instr);
175
176 ss1 >> (newlib->m_xrodcent) >> (newlib->m_yrodcent) >> (newlib->m_step);
177
178 std::getline(filestr,instr);
179 newlib->m_comment = std::move(instr);
180
181 return newlib;
182 }
183
184 std::vector<EnergySpot>* FCALDistEnergyShowerLib::getShower(const G4Track* track, ShowerLibStatistics* stats, int randomShift) const
185 {
186 if (!m_filled) {
187 std::cout << "Library is not created for production use" << std::endl;
188 return nullptr;
189 }
190
191 //std::cout << "Starting getShower()" << std::endl;
192
193 double x = track->GetPosition().getX();
194 double y = track->GetPosition().getY();
195
196 float dist = this->distance(x,y);
197
198 G4int particleCode = track->GetDefinition()->GetPDGEncoding();
199 if ( particleCode < 0 ) particleCode = -particleCode; // hack for positrons.
200
201 if ( particleCode != m_particle ) {
202 std::cout << "wrong particle: " << particleCode << "!=" << m_particle << std::endl;
203 return nullptr;
204 }
205
206 //std::cout << "Extracted particle parameters: " << eta << std::endl;
207
208 library::const_iterator libit = m_libData.upper_bound(dist);
209
210 if (libit == m_libData.begin()) {
211 //this is really weird
212 std::cout << "Something is wrong with dist: x=" << x << " y=" << y << " dist=" << dist << std::endl;
213 } else {
214 --libit;
215 }
216
217 //std::cout << "Found eta bin" << std::endl;
218
219 if ((*libit).second.empty()) {
220 std::cout << "The etabin corresponding to the eta is empty" << std::endl;
221 return nullptr;
222 }
223 double trenergy = track->GetKineticEnergy();
224 //std::cout << "Energy: " << trenergy << std::endl;
225 distbin::const_iterator distit = (*libit).second.lower_bound(trenergy);
226 if (distit == (*libit).second.end()) --distit; //if it points to the end, repoint it to the last shower in bin
227 else if (distit != (*libit).second.begin()) { //find the closest energy. it's either found item or previous (if there is a previous)
228 distbin::const_iterator distitch = distit;
229 --distitch;
230 if (((*distit).first - trenergy) > (trenergy - (*distitch).first)) { // the closest is the previous
231 --distit;
232 }
233 }
234 //std::cout << "Found closest energy:" << (*distit).first << std::endl;
235 if (randomShift > 0) {
236 double upperEnergy = (*distit).first * 1.01; //we allow 1% off
237 for (int i = 0; i < randomShift; i++) {
238 ++distit;
239 if (distit == (*libit).second.end()) {
240 --distit;
241 break;
242 }
243 if ((*distit).first > upperEnergy) break; //energy diff too high
244 }
245 }
246 if ((randomShift < 0)&&(distit != (*libit).second.begin())) {
247 double lowerEnergy = (*distit).first * 0.99; //we allow 1% off
248 for (int i = 0; i > randomShift; i--) {
249 --distit;
250 if (distit == (*libit).second.begin()) {
251 //distit++;
252 break;
253 }
254 if (lowerEnergy > (*distit).first) break; //energy diff too high
255 }
256 }
257 //std::cout << "Found out energy" << std::endl;
258 //std::cout << "Shower with num hits:" << (*etait).second.size() << std::endl;
259 std::vector<EnergySpot>* outshower = new std::vector<EnergySpot>();//((*etait).second);
260 Shower::const_iterator iter;
261 //std::cout << "Created out shower" << std::endl;
262
263 float energyScale = track->GetKineticEnergy() / (*distit).first;
264 //std::cout << "Scale: " << energyScale << std::endl;
265
266 for (iter = (*distit).second.begin() /*outshower->begin()*/; iter != (*distit).second.end() /*outshower->end()*/; ++iter) {
267 EnergySpot tmp( (*iter)->GetPosition(), (*iter)->GetEnergy(), (*iter)->GetTime() );
268 tmp.SetEnergy(tmp.GetEnergy() * energyScale);
269 outshower->push_back(tmp);
270 //(*iter).SetEnergy((*iter).GetEnergy() * energyScale);
271 }
272 //std::cout << "Scaled" << std::endl;
273 if (stats != nullptr) {
274 stats->recordShowerLibUse(calcKey((*libit).first));
275 }
276 //std::cout << "Done" << std::endl;
277 return outshower;
278 }
279
280 double FCALDistEnergyShowerLib::getContainmentZ(const G4Track* track) const
281 {
282 if (!m_filled) {
283 std::cout << "Library is not created for production use" << std::endl;
284 return 0;
285 }
286
287 //std::cout << "Starting getShower()" << std::endl;
288
289 double x = track->GetPosition().getX();
290 double y = track->GetPosition().getY();
291
292 float dist = this->distance(x,y);
293
294 G4int particleCode = track->GetDefinition()->GetPDGEncoding();
295 if ( particleCode < 0 ) particleCode = -particleCode; // hack for positrons.
296
297 if ( particleCode != m_particle ) {
298 std::cout << "wrong particle: " << particleCode << "!=" << m_particle << std::endl;
299 return 0;
300 }
301
302 library::const_iterator libit = m_libData.upper_bound(dist);
303
304 if (libit == m_libData.begin()) {
305 //this is really weird
306 std::cout << "Something is wrong with dist: x=" << x << " y=" << y << " dist=" << dist << std::endl;
307 } else {
308 --libit;
309 }
310
311 //std::cout << "Found eta bin" << std::endl;
312
313 if ((*libit).second.empty()) {
314 std::cout << "The etabin corresponding to the eta is empty" << std::endl;
315 return 0;
316 }
317 double trenergy = track->GetKineticEnergy();
318 //std::cout << "Energy: " << trenergy << std::endl;
319 distbin::const_iterator distit = (*libit).second.lower_bound(trenergy);
320 if (distit == (*libit).second.end()) --distit; //if it points to the end, repoint it to the last shower in bin
321 else if (distit != (*libit).second.begin()) { //find the closest energy. it's either found item or previous (if there is a previous)
322 distbin::const_iterator distitch = distit;
323 --distitch;
324 if (((*distit).first - trenergy) > (trenergy - (*distitch).first)) { // the closest is the previous
325 --distit;
326 }
327 }
328 //std::cout << "Found closest energy:" << (*etait).first << std::endl;
329 double rezZ = (*distit).second.getZSize();
330 distbin::const_iterator distiter = distit;
331 int actualNumFS = 1;
332 int spread = 2; //will calculate average Z for 5 showers max ( -2 .. 0 .. +2 )
333 double upperEnergy = (*distit).first * 1.01; //we allow 1% off
334 for (int i = 0; i < spread; i++) {
335 ++distiter;
336 if (distiter == (*libit).second.end()) {
337 break;
338 }
339 if (upperEnergy < (*distiter).first) break; //energy diff too high
340 //the shower is OK, including it to the average
341 rezZ += (*distiter).second.getZSize();
342 actualNumFS++;
343 }
344 distiter = distit;
345 if (distiter != (*libit).second.begin()) {
346 double lowerEnergy = (*distit).first * 0.99; //we allow 1% off
347 for (int i = 0; i < spread; i++) {
348 --distiter;
349 if (lowerEnergy > (*distiter).first) break; //energy diff too high
350 //the shower is OK, including it to the average
351 rezZ += (*distiter).second.getZSize();
352 actualNumFS++;
353 if (distiter == (*libit).second.begin()) {
354 break;
355 }
356 }
357 }
358 return rezZ/actualNumFS; //average Z size
359 }
360
361 double FCALDistEnergyShowerLib::getContainmentR(const G4Track* track) const
362 {
363 if (!m_filled) {
364 std::cout << "Library is not created for production use" << std::endl;
365 return 0;
366 }
367
368 //std::cout << "Starting getShower()" << std::endl;
369
370 double x = track->GetPosition().getX();
371 double y = track->GetPosition().getY();
372
373 float dist = this->distance(x,y);
374
375 G4int particleCode = track->GetDefinition()->GetPDGEncoding();
376 if ( particleCode < 0 ) particleCode = -particleCode; // hack for positrons.
377
378 if ( particleCode != m_particle ) {
379 std::cout << "wrong particle: " << particleCode << "!=" << m_particle << std::endl;
380 return 0;
381 }
382
383 library::const_iterator libit = m_libData.upper_bound(dist);
384
385 if (libit == m_libData.begin()) {
386 //this is really weird
387 std::cout << "Something is wrong with dist: x=" << x << " y=" << y << " dist=" << dist << std::endl;
388 } else {
389 --libit;
390 }
391
392 //std::cout << "Found eta bin" << std::endl;
393
394 if ((*libit).second.empty()) {
395 std::cout << "The etabin corresponding to the eta is empty" << std::endl;
396 return 0;
397 }
398 double trenergy = track->GetKineticEnergy();
399 //std::cout << "Energy: " << trenergy << std::endl;
400 distbin::const_iterator distit = (*libit).second.lower_bound(trenergy);
401 if (distit == (*libit).second.end()) --distit; //if it points to the end, repoint it to the last shower in bin
402 else if (distit != (*libit).second.begin()) { //find the closest energy. it's either found item or previous (if there is a previous)
403 distbin::const_iterator distitch = distit;
404 --distitch;
405 if (((*distit).first - trenergy) > (trenergy - (*distitch).first)) { // the closest is the previous
406 --distit;
407 }
408 }
409 //std::cout << "Found closest energy:" << (*etait).first << std::endl;
410 double rezR = (*distit).second.getRSize();
411 distbin::const_iterator distiter = distit;
412 int actualNumFS = 1;
413 int spread = 2; //will calculate average Z for 5 showers max ( -2 .. 0 .. +2 )
414 double upperEnergy = (*distit).first * 1.01; //we allow 1% off
415 for (int i = 0; i < spread; i++) {
416 ++distiter;
417 if (distiter == (*libit).second.end()) {
418 break;
419 }
420 if (upperEnergy < (*distiter).first) break; //energy diff too high
421 //the shower is OK, including it to the average
422 rezR += (*distiter).second.getRSize();
423 actualNumFS++;
424 }
425 distiter = distit;
426 if (distiter != (*libit).second.begin()) {
427 double lowerEnergy = (*distit).first * 0.99; //we allow 1% off
428 for (int i = 0; i < spread; i++) {
429 --distiter;
430 if (lowerEnergy > (*distiter).first) break; //energy diff too high
431 //the shower is OK, including it to the average
432 rezR += (*distiter).second.getRSize();
433 actualNumFS++;
434 if (distiter == (*libit).second.begin()) {
435 break;
436 }
437 }
438 }
439 return rezR/actualNumFS; //average Z size
440 }
441
443 {
444 if (m_filled) {
445 std::cout << "ERROR: filled" << std::endl;
446 return false;
447 }
448
449 double x = genParticle->production_vertex()->position().x();
450 double y = genParticle->production_vertex()->position().y();
451
452 float dist = this->distance(x,y);
453
454 if ( genParticle->pdg_id() != m_particle ) {
455 std::cout << "ERROR: wrong pdgcode: " << m_particle << " != " << genParticle->pdg_id() << std::endl;
456 return false;
457 }
458
459 library::iterator libit = m_libData.upper_bound(dist);
460 if (libit == m_libData.begin()) {
461 //this is really weird
462 std::cout << "Something is wrong with dist: x=" << x << " y=" << y << " dist=" << dist << std::endl;
463 } else {
464 --libit;
465 }
466 (*libit).second.insert(distbin::value_type(genParticle->momentum().e(),(*shower)));
467 return true;
468 }
469
471 {
472 if (m_libData.empty()) return false;
473 TParameter<int> ver("version",LIB_VERSION);
474
475 dest->WriteObject(&ver,"version");
476
477 TTree TTreeMeta;
478 TTree TTreeLib;
479
480 write(&TTreeLib);
481 writeMeta(&TTreeMeta);
482
483 dest->WriteObject(&TTreeLib,"library");
484 dest->WriteObject(&TTreeMeta,"meta");
485
486 return true;
487 }
488
489
491 {
492 /*
493 * Dist Energy library format:
494 * | x | y | z | e | time | - name of branch in TTree
495 * ------------------------------------------------------------------
496 * | xrod cent | yrod cent |step (roddist)| not | not | - library header
497 * | (parameter) | (parameter) | (parameter) | used | used |
498 * ------------------------------------------------------------------
499 * |num of showers| min dist for | not | not | not | - dist bin header
500 * | in dist bin | cur dist bin | used | used | used |
501 * ------------------------------------------------------------------
502 * | num of hits |shower r-size |shower z-size | truth | not | - shower header
503 * | in shower |for cont.check|for cont.check| energy | used |
504 * ------------------------------------------------------------------
505 * |x-coord of hit|y-coord of hit|z-coord of hit|dep.energy|hit time| - hit
506 */
507 int nentr = source->GetEntriesFast();
508 if (nentr < 3) return false;
509 Float_t x,y,z,e,time;
510 source->SetBranchAddress("x",&x);
511 source->SetBranchAddress("y",&y);
512 source->SetBranchAddress("z",&z);
513 source->SetBranchAddress("e",&e);
514 source->SetBranchAddress("time",&time);
515 int entr = 0;
516
517 source->GetEntry(entr++);
518 m_xrodcent = x;
519 m_yrodcent = y;
520 m_step = z;
521
522 do {
523 //read eta bin header
524 source->GetEntry(entr++); //x - nshowers, y - min dist in the current dist bin
525 int nsh = (int)(x+0.1); // +0.1 just in case - c++ has low round
526 float curDist = y;
527 distbin * curbin = &(m_libData[curDist]); //creating a new dist bin
528 for(int i = 0; i < nsh; i++) {
529 //read shower header
530 source->GetEntry(entr++); //x - nhits, y - r size, z - z size, e - gen energy
531 int nhits = (int)(x+0.1);
532 float curEnergy = e;
533 Shower * shower = &((*curbin)[curEnergy]);
534 shower->setRSize(y);
535 shower->setZSize(z);
536 for(int j = 0; j < nhits; j++) {
537 source->GetEntry(entr++); //variables mean what the name suggests
538 shower->push_back(new ShowerEnergySpot(G4ThreeVector(x,y,z),e,time));
539 }
540 }
541 } while (entr < nentr);
542
543 if (entr != nentr) {
544 return false;
545 }
546
547 m_filled = true;
548 return true;
549 }
550
551 bool FCALDistEnergyShowerLib::write(TTree* dest) const
552 {
553 /*
554 * Dist Energy library format:
555 * | x | y | z | e | time | - name of branch in TTree
556 * ------------------------------------------------------------------
557 * | xrod cent | yrod cent |step (roddist)| not | not | - library header
558 * | (parameter) | (parameter) | (parameter) | used | used |
559 * ------------------------------------------------------------------
560 * |num of showers| min dist for | not | not | not | - dist bin header
561 * | in dist bin | cur dist bin | used | used | used |
562 * ------------------------------------------------------------------
563 * | num of hits |shower r-size |shower z-size | truth | not | - shower header
564 * | in shower |for cont.check|for cont.check| energy | used |
565 * ------------------------------------------------------------------
566 * |x-coord of hit|y-coord of hit|z-coord of hit|dep.energy|hit time| - hit
567 */
568 Float_t x,y,z,e,time;
569 dest->Branch("x",&x);
570 dest->Branch("y",&y);
571 dest->Branch("z",&z);
572 dest->Branch("e",&e);
573 dest->Branch("time",&time);
574
575 x = m_xrodcent;
576 y = m_yrodcent;
577 z = m_step;
578 e = 0;
579 time = 0;
580 dest->Fill();
581
582 library::const_iterator libit;
583 for (libit = m_libData.begin(); libit != m_libData.end(); ++libit) {
584 x = (*libit).second.size();
585 y = (*libit).first;
586 z = 0;
587 e = 0;
588 time = 0;
589 dest->Fill(); //eta bin header
590 distbin::const_iterator distit;
591 for (distit = (*libit).second.begin(); distit != (*libit).second.end(); ++distit) {
592 x = (*distit).second.size();
593 y = (*distit).second.getRSize();
594 z = (*distit).second.getZSize();
595 e = (*distit).first;
596 time = 0;
597 dest->Fill(); //shower header
598 Shower::const_iterator iter;
599 for (iter = (*distit).second.begin(); iter != (*distit).second.end(); ++iter) {
600 x = (*iter)->GetPosition().x();
601 y = (*iter)->GetPosition().y();
602 z = (*iter)->GetPosition().z();
603 e = (*iter)->GetEnergy();
604 time = (*iter)->GetTime();
605 dest->Fill();
606 }
607 }
608 }
609 //dest->Write();
610 return true;
611 }
612
614 {
615 std::stringstream ss;
616 ss << std::fixed << std::setprecision(3);
617 ss << "Distance calculator parameters: xrod_cent=" << m_xrodcent << " yrod_cent=" << m_yrodcent << " step=" << m_step;
618 return ss.str();
619 }
620
621 bool FCALDistEnergyShowerLib::readStructure(std::vector<float>& structure)
622 {
623 std::vector<float>::const_iterator iter;
624
625 for (iter = structure.begin(); iter != structure.end(); ++iter) {
626 m_libData[(*iter)];
627 }
628
629 return true;
630 }
631
633 {
634 std::map<int, std::string> names;
635 std::map<int, int> sizes;
636 for(library::const_iterator it = m_libData.begin(); it != m_libData.end(); ++it) {
637 sizes[calcKey(it->first)]=it->second.size();
638 float distlow = it->first;
639 float disthigh;
640 library::const_iterator it_copy = it;
641 ++it_copy;
642 if (it_copy == m_libData.end()) {
643 disthigh = 4.5;
644 } else {
645 disthigh = it_copy->first;
646 }
647 std::stringstream ss;
648 ss << std::showpos << std::fixed << std::setprecision(2);
649 ss << "DIST: " << distlow << " .. " << disthigh;
650 names[calcKey(it->first)] = ss.str();
651 }
652 return new ShowerLibStatistics(names, sizes);
653 }
654
655} // namespace ShowerLib
Scalar eta() const
pseudorapidity method
#define LIB_VERSION
static Double_t ss
#define y
#define x
#define z
virtual ShowerLibStatistics * createStatistics() const
bool write(TTree *dest) const
write library to given TTree
virtual const std::string printParameters() const
static IShowerLib * readFromROOTFile(TFile *source)
factory method. create a library from root file. returns NULL if file is invalid.
virtual bool storeShower(HepMC::ConstGenParticlePtr genParticle, const Shower *shower)
store shower in the library
bool read(TTree *source)
read library from given TTree
virtual bool writeToROOT(TFile *dest)
write library to ROOT file
static IShowerLib * createEmptyLib(const std::string &inputFile)
factory method. create empty library with the given structure. returns NULL if file is invalid.
float distance(double x, double y) const
virtual std::vector< EnergySpot > * getShower(const G4Track *track, ShowerLibStatistics *stats, int randomShift) const
get shower for given G4 track
bool readStructure(std::vector< float > &structure)
virtual double getContainmentR(const G4Track *track) const
get average lateral spread of the showers for the given energy
virtual double getContainmentZ(const G4Track *track) const
get average length of showers for the given energy
bool readMeta(TTree *source)
read metadata from the given TTree
IShowerLib()
default constructor
Definition IShowerLib.h:98
int m_particle
ID of the generated particles.
Definition IShowerLib.h:107
bool writeMeta(TTree *dest) const
write metadata to the given TTree
std::string m_detector
name of the detector
Definition IShowerLib.h:106
std::string m_comment
comment
Definition IShowerLib.h:112
bool m_filled
is the library read from ROOT or from structure file
Definition IShowerLib.h:114
Class for shower library shower.
Definition Shower.h:36
void setZSize(const float zsize)
Definition Shower.h:64
void setRSize(const float rsize)
Definition Shower.h:65
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
Namespace for the ShowerLib related classes.
Definition StepInfo.h:17
int calcKey(float eta)