ATLAS Offline Software
Loading...
Searching...
No Matches
PunchThroughG4Tool.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// PunchThrough headers
9
10// PathResolver
12
13// standard C++ libraries
14#include <iostream>
15#include <sstream>
16#include <string>
17#include <algorithm>
18#include <vector>
19#include <numeric>
20#include <string_view>
21#include <charconv>
22
23// standard C libraries
24#include <cmath>
25
26// CLHEP
27#include "CLHEP/Random/RandFlat.h"
28
29// ROOT
30#include "TFile.h"
31#include "TH2F.h"
32#include "TH1.h"
33#include "TMath.h"
34#include "TKey.h"
35
36// namespace Amg functions
39
40using namespace CxxUtils;
41
42PunchThroughG4Tool::PunchThroughG4Tool(const std::string& type, const std::string& name, const IInterface* parent)
43 : base_class(type,name,parent)
44{
45}
46
47// Athena method, used to get out the G4 geometry and set up the Fast Simulation Models
49 ATH_MSG_DEBUG("[PunchThroughG4Tool] ==>" << name() << "::initialize()");
50
51 //retrieve inverse CDF config file
53 {
54 ATH_MSG_WARNING("[PunchThroughG4Tool] unable to open or read the inverse CDF config");
55 return StatusCode::FAILURE;
56 }
57
58 //retrieve inverse PCA config file
60 {
61 ATH_MSG_WARNING("[PunchThroughG4Tool] unable to open or read the inverse PCA config");
62 return StatusCode::FAILURE;
63 }
64
65 //check first the size of infoMap for both PCA and CDF, they should be equal
66 if (!(m_xml_info_pca.size() == m_xml_info_cdf.size()))
67 {
68 ATH_MSG_WARNING("[PunchThroughG4Tool] size of infoMap for PCA and CDF differs! Something is wrong with input xml files.");
69 return StatusCode::FAILURE;
70 }
71
72 ATH_MSG_INFO( "[PunchThroughG4Tool] initialization is successful" );
73 return StatusCode::SUCCESS;
74}
75
77 G4ParticleTable *ptable = G4ParticleTable::GetParticleTable();
78 for ( unsigned int num = 0; num < m_punchThroughParticles.size(); num++ )
79 {
80 const int pdg = m_punchThroughParticles[num];
81 // if no information is given on the creation of anti-particles -> do not simulate anti-particles
82 const bool doAnti = ( num < m_doAntiParticles.size() ) ? m_doAntiParticles[num] : false;
83 // if no information is given on the minimum energy -> take 50. MeV as default
84 const double minEnergy = ( num < m_minEnergy.size() ) ? m_minEnergy[num] : 50.;
85 // if no information is given on the maximum number of punch-through particles -> take -1 as default
86 const int maxNum = ( num < m_minEnergy.size() ) ? m_maxNumParticles[num] : -1;
87 // if no information is given on the scale factor for the number of particles -> take 1. as defaulft
88 const double numFactor = ( num < m_numParticlesFactor.size() ) ? m_numParticlesFactor[num] : 1.;
89 // if no information is given on the position angle factor -> take 1.
90 const double posAngleFactor = ( num < m_posAngleFactor.size() ) ? m_posAngleFactor[num] : 1.;
91 // if no information is given on the momentum angle factor -> take 1.
92 const double momAngleFactor = ( num < m_momAngleFactor.size() ) ? m_momAngleFactor[num] : 1.;
93 // if no information is given on the scale factor for the energy -> take 1. as default
94 const double energyFactor = ( num < m_energyFactor.size() ) ? m_energyFactor[num] : 1.;
95
96 // register the particle
97 ATH_MSG_VERBOSE("VERBOSE: [PunchThroughG4Tool] registering punch-through particle type with pdg = " << pdg );
98 if (registerPunchThroughParticle( *ptable, pdg, doAnti, minEnergy, maxNum, numFactor, energyFactor, posAngleFactor, momAngleFactor ) != StatusCode::SUCCESS)
99 {
100 ATH_MSG_ERROR("[PunchThroughG4Tool] unable to register punch-through particle type with pdg = " << pdg);
101 return StatusCode::FAILURE;
102 }
103 }
104 // if all goes well
105 return StatusCode::SUCCESS;
106}
107
109 // resolving lookuptable file
110
111 ATH_MSG_INFO("[PunchThroughG4Tool] PunchThroughG4Tool::initializePhysics() called");
112 std::string resolvedFileName = PathResolverFindCalibFile (m_filenameLookupTable);
113 if (resolvedFileName.empty()) {
114 ATH_MSG_ERROR( "[PunchThroughG4Tool] Parametrisation file '" << m_filenameLookupTable << "' not found" );
115 return StatusCode::FAILURE;
116 }
117 ATH_MSG_DEBUG( "[PunchThroughG4Tool] Parametrisation file found: " << resolvedFileName );
118 // open the LookupTable file
119 m_fileLookupTable = new TFile( resolvedFileName.c_str(), "READ");
120 if (!m_fileLookupTable) {
121 ATH_MSG_WARNING("[PunchThroughG4Tool] unable to open the lookup-table for the punch-through simulation (file does not exist)");
122 return StatusCode::FAILURE;
123 }
124
125 if (!m_fileLookupTable->IsOpen()) {
126 ATH_MSG_WARNING("[PunchThroughG4Tool] unable to open the lookup-table for the punch-through simulation (wrong or empty file?)");
127 return StatusCode::FAILURE;
128 }
129
130 //--------------------------------------------------------------------------------
131 // register all the punch-through particles which will be simulated
133
134 //--------------------------------------------------------------------------------
135 // register all the punch-through correlations between particles to be simulated
137
138 //--------------------------------------------------------------------------------
139 // close the file with the lookuptable
140 m_fileLookupTable->Close();
141
142 // Geometry identifier service
143 if ( !m_geoIDSvc.empty() && m_geoIDSvc.retrieve().isFailure())
144 {
145 ATH_MSG_FATAL ( "[PunchThroughG4Tool] Could not retrieve GeometryIdentifier Service. Abort");
146 return StatusCode::FAILURE;
147 }
148
149 //envelope definition service
150 if (m_envDefSvc.retrieve().isFailure() )
151 {
152 ATH_MSG_ERROR( "[PunchThroughG4Tool] Could not retrieve " << m_envDefSvc );
153 return StatusCode::FAILURE;
154 }
155
156 // direct-initialization of the Muon and Calo boundaries (giving R, Z coordinates)
157 // can also do like below (RZPairVector is vector<pair<double, double>>)
158 // const RZPairVector* rzMS = &(m_envDefSvc->getMuonRZBoundary());
159 const std::vector<std::pair<double, double>>* rzMS = &(m_envDefSvc->getMuonRZBoundary());
160 const std::vector<std::pair<double, double>>* rzCalo = &(m_envDefSvc->getCaloRZBoundary());
161 ATH_CHECK(checkCaloMSBoundaries(rzMS, rzCalo));
162 return StatusCode::SUCCESS;
163}
164
166 //--------------------------------------------------------------------------------
167 // TODO: implement punch-through parameters for different m_pdgInitiators
168 // currently m_pdgInitiators is only used to filter out particles,
169 // inside computePunchThroughParticles
170
171 // check if more correlations were given than particle types were registered
172 unsigned int numCorrelations = m_correlatedParticle.size();
173 if ( numCorrelations > m_punchThroughParticles.size() )
174 {
175 ATH_MSG_WARNING("[PunchThroughG4Tool] more punch-through particle correlations are given, than punch-through particle types are registered (skipping the last ones)");
176 numCorrelations = m_punchThroughParticles.size();
177 }
178
179 // now register correlation between particles
180 for ( unsigned int num = 0; num < numCorrelations; num++ )
181 {
182 const int pdg1 = m_punchThroughParticles[num];
183 const int pdg2 = m_correlatedParticle[num];
184 const double fullCorrEnergy = ( num < m_fullCorrEnergy.size() ) ? m_fullCorrEnergy[num] : 0.;
185 const double minCorrEnergy = ( num < m_minCorrEnergy.size() ) ? m_minCorrEnergy[num] : 0.;
186
187 // if correlatedParticle==0 is given -> no correlation
188 if ( ! pdg2) continue;
189 // register it
190 if ( registerCorrelation(pdg1, pdg2, minCorrEnergy, fullCorrEnergy) != StatusCode::SUCCESS )
191 {
192 ATH_MSG_ERROR("[PunchThroughG4Tool] unable to register punch-through particle correlation for pdg1=" << pdg1 << " pdg2=" << pdg2 );
193 return StatusCode::FAILURE;
194 }
195 }
196 // if all goes well
197 return StatusCode::SUCCESS;
198}
199
201 //m_R1, m_R2, m_z1, m_z2
202 std::vector<double> caloMSVarsVect = {m_R1, m_R2, m_z1, m_z2};
203 ATH_MSG_DEBUG("[PunchThroughG4Tool] getCaloMSVars()==> m_R1 = "<< caloMSVarsVect[0] << ", m_R2 = " << caloMSVarsVect[1] << ", m_z1 = " << caloMSVarsVect[2] << ", m_z2 = " << caloMSVarsVect[3]);
204
205 return caloMSVarsVect;
206}
207
208StatusCode PunchThroughG4Tool::checkCaloMSBoundaries(const std::vector<std::pair<double, double>>* rzMS,
209 const std::vector<std::pair<double, double>>* rzCalo)
210{
211 bool found1, found2;
212 found1=false; found2=false;
213
214 // procedure below is to find intersection of calo and MS boundary (right at the edge of calo and into MS)
215 // at least 2 distinct calo-MS boundary points must be found to proceed
216 for ( const auto & [r_tempCalo , z_tempCalo] : *rzCalo )
217 {
218 // calo point must be more than beampipe radius
219 if (r_tempCalo> m_beamPipe)
220 {
221 for ( const auto & [r_tempMS , z_tempMS] : *rzMS )
222 {
223 // first point, if R and Z of calo matches R and Z of MS
224 if (r_tempCalo==r_tempMS && z_tempCalo==z_tempMS && found1==false )
225 {
226 m_R1=r_tempMS;
227 m_z1=std::fabs(z_tempMS);
228 found1=true;
229 continue;
230 }
231 // second point, if R and Z of calo matches R and Z of MS
232 else if (r_tempCalo==r_tempMS && z_tempCalo==z_tempMS && r_tempCalo!=m_R1 && std::fabs(z_tempCalo)!=m_z1)
233 {
234 m_R2=r_tempMS;
235 m_z2=std::fabs(z_tempMS);
236 found2=true;
237 }
238 }
239
240 if (found1==true && found2==true) break;
241 }
242 }
243
244 //in case geometry description changes
245 if (found1 == false){
246 ATH_MSG_ERROR ("[PunchThroughG4Tool] first coordinate of calo-MS border not found");
247 return StatusCode::FAILURE;
248 }
249 if (found2 == false){
250 ATH_MSG_ERROR ("[PunchThroughG4Tool] second coordinate of calo-MS border not found; first one is: R1 ="<<m_R1<<" z1 ="<<m_z1);
251 return StatusCode::FAILURE;
252 }
253
254 //now order the found values
255 if (m_R1>m_R2) { std::swap(m_R1, m_R2); } //m_R1 - smaller one
256 if (m_z1<m_z2) { std::swap(m_z1, m_z2); } //m_z1 - bigger one
257
258 if (m_R1==m_R2 || m_z1==m_z2){
259 ATH_MSG_ERROR ("[PunchThroughG4Tool] Bug in propagation calculation! R1="<<m_R1<<" R2 = "<<m_R2<<" z1="<<m_z1<<" z2= "<<m_z2 );
260 return StatusCode::FAILURE;
261 }
262
263 ATH_MSG_DEBUG ("[PunchThroughG4Tool] calo-MS boundary coordinates: R1="<<m_R1<<" R2 = "<<m_R2<<" z1="<<m_z1<<" z2= "<<m_z2);
264
265 return StatusCode::SUCCESS;
266}
267
269{
270 ATH_MSG_DEBUG( "[PunchThroughG4Tool] finalize() starting" );
271 for(auto & each : m_particles) {
272 delete each.second;
273 }
274
275 ATH_MSG_DEBUG( "[PunchThroughG4Tool] finalize() successful" );
276
277 return StatusCode::SUCCESS;
278}
279
280void PunchThroughG4Tool::checkParticleTable(G4ParticleTable &ptable, int secondarySignedPDG){
281 // Get Geant4 particle definition
282 G4ParticleDefinition* secG4Particle = ptable.FindParticle(secondarySignedPDG);
283 double mass = secG4Particle->GetPDGMass();
284 ATH_MSG_DEBUG("[PunchThroughG4Tool::checkParticleTable] secondarySignedPDG = " << secondarySignedPDG << ", mass = "<< mass);
285}
286
287std::vector<std::map<std::string, double>> PunchThroughG4Tool::computePunchThroughParticles(const G4FastTrack& fastTrack, CLHEP::HepRandomEngine* rndmEngine, double punchThroughProbability, double punchThroughClassifierRand)
288{
289 ATH_MSG_DEBUG( "[PunchThroughG4Tool] starting punch-through simulation");
290
291 // holder to all secondary particle kinematics
292 std::vector<std::map<std::string, double>> secKinematicsMapVect;
293
294 // Get primary Geant4 track (main particle producing punchthrough)
295 const G4Track * g4PrimaryTrack = fastTrack.GetPrimaryTrack();
296
297 // Get Geant4 particle definition
298 const G4ParticleDefinition * mainG4Particle = g4PrimaryTrack->GetDefinition();
299
300 // Printing for debugging
301 ATH_MSG_DEBUG("PunchThroughG4Tool::computePunchThroughParticles Debug:"
302 << " Incoming particle:"
303 << " E=" << g4PrimaryTrack->GetKineticEnergy()/CLHEP::GeV << " GeV"
304 << " theta=" << g4PrimaryTrack->GetMomentumDirection().theta()
305 << " phi=" << g4PrimaryTrack->GetMomentumDirection().phi()
306 << " PDG=" << (mainG4Particle ? mainG4Particle->GetPDGEncoding() : -999)
307 << " Name=" << (mainG4Particle ? mainG4Particle->GetParticleName() : "UNKNOWN"));
308
309 // Get primary Geant4 particle properties
310 int pdgID = mainG4Particle->GetPDGEncoding();
311 float mainPartMass = mainG4Particle->GetPDGMass();
312 float mainMomMag2 = g4PrimaryTrack->GetMomentum().mag2();
313 float mainPartEta = g4PrimaryTrack->GetPosition().eta();
314 // get primary track's current position and momentum direction
315 const G4ThreeVector mainMomentumDir = g4PrimaryTrack->GetMomentumDirection();
316 const G4ThreeVector mainMomentum = g4PrimaryTrack->GetMomentum();
317 const G4ThreeVector mainPosition = g4PrimaryTrack->GetPosition();
318
319 // calculate incoming energy and eta
320 const double initEnergy = std::sqrt( mainMomMag2 + mainPartMass*mainPartMass );
321 const double initEta = mainPartEta;
322
323 ATH_MSG_DEBUG("[PunchThroughG4Tool] position of the input particle: r"<<mainPosition.perp()<<" z= "<<mainPosition.z() );
324
325 // check if it points to the calorimeter - if not, don't simulate
326 // first convert to Amg::Vector3D
327 const Amg::Vector3D amg3DPosVec(mainPosition.x(), mainPosition.y(), mainPosition.z());
328 const Amg::Vector3D amg3DMomDirVec(mainMomentumDir.x(), mainMomentumDir.y(), mainMomentumDir.z());
329 const Amg::Vector3D amg3DMomVec(mainMomentum.x(), mainMomentum.y(), mainMomentum.z());
330 // find nextGeoID to get AtlasDetDescr::AtlasRegion
331 AtlasDetDescr::AtlasRegion nextGeoID = m_geoIDSvc->identifyNextGeoID(amg3DPosVec,amg3DMomVec);
332 if ( nextGeoID != AtlasDetDescr::fAtlasCalo)
333 {
334 ATH_MSG_DEBUG("[PunchThroughG4Tool](GeoIDSvc) input particle doesn't point to calorimeter"<< "Next GeoID: "<< nextGeoID );
335 //return nullptr;
336 return {};
337 }
338
339 // check if the particle's pdg is registered as a punch-through-causing type
340 {
341 std::vector<int>::const_iterator pdgIt = m_pdgInitiators.begin();
342 std::vector<int>::const_iterator pdgItEnd = m_pdgInitiators.end();
343
344 std::vector<int>::const_iterator minEnergyIt = m_initiatorsMinEnergy.begin();
345 // loop over all known punch-through initiators
346 for ( ; pdgIt != pdgItEnd; ++pdgIt, ++minEnergyIt)
347 {
348 if (std::abs(pdgID) == *pdgIt){
349 if(initEnergy < *minEnergyIt){
350 ATH_MSG_DEBUG("[PunchThroughG4Tool] particle does not meet initiator min energy requirement. Dropping it in the calo.");
351 //return nullptr;
352 return {};
353 }
354 break;
355 }
356 }
357
358 // particle will not cause punch-through -> bail out
359 if (pdgIt == pdgItEnd)
360 {
361 ATH_MSG_DEBUG("[PunchThroughG4Tool] particle is not registered as punch-through initiator. Dropping it in the calo.");
362 //return nullptr;
363 return {};
364 }
365 }
366
367 if(initEta < m_initiatorsEtaRange.value().at(0) || initEta > m_initiatorsEtaRange.value().at(1) ){
368 ATH_MSG_DEBUG("[PunchThroughG4Tool] particle does not meet initiator eta range requirement. Dropping it in the calo.");
369 //return nullptr;
370 return {};
371 }
372
373 ATH_MSG_DEBUG("[PunchThroughG4Tool] punchThroughProbability output: " << punchThroughProbability << " RandFlat: " << punchThroughClassifierRand );
374
375 //Safety check: If probability < random number then don't simulate punch through
376 if( punchThroughClassifierRand > punchThroughProbability){
377 ATH_MSG_DEBUG("[PunchThroughG4Tool] particle not classified to create punch through. Dropping it in the calo.");
378 //return nullptr;
379 return {};
380 }
381
382 //if initial particle is on ID surface, points to the calorimeter, is a punch-through initiator, meets initiator min enery and eta range
383
384 // this is the place where the magic is done:
385 // test for each registered punch-through pdg if a punch-through
386 // occures and create these particles
387 // -> therefore loop over all registered pdg ids
388 // to keep track of the correlated particles which were already simulated:
389 // first int is pdg, second int is number of particles created
390
391 // interpolate energy and eta
392 const double interpEnergy = interpolateEnergy(initEnergy, rndmEngine);
393 const double interpEta = interpolateEta(initEta, rndmEngine);
394
395 std::map<int, int> corrPdgNumDone;
396
397 int maxTries = 10;
398 int nTries = 0;
399
400 // loop over all particle pdgs
401 ATH_MSG_DEBUG("[DEBUG] Enter computePunchThroughParticles: starting while (retry) loop");
402 while(secKinematicsMapVect.empty() && nTries < maxTries) { //ensure we always create at least one punch through particle, maxTries to catch very rare cases
403 // loop over all particle pdgs
404 ATH_MSG_DEBUG("[DEBUG] Retry " << nTries << ", secKinematicsMapVect.size=" << secKinematicsMapVect.size());
405 for (const auto& currentParticle : m_particles)
406 {
407 // the pdg that is currently treated
408 int doPdg = currentParticle.first;
409 // get the current particle's correlated pdg
410 int corrPdg = currentParticle.second->getCorrelatedPdg();
411
412 // if there is a correlated particle type to this one
413 if (corrPdg)
414 {
415 // find out if the current pdg was already simulated
416 std::map<int,int>::iterator itPdgPos = corrPdgNumDone.find(doPdg);
417 // if the current pdg was not simulated yet, find out if
418 // it's correlated one was simulated
419 if ( itPdgPos == corrPdgNumDone.end() ) itPdgPos = corrPdgNumDone.find(corrPdg);
420
421 // neither this nor the correlated particle type was simulated
422 // so far:
423 if ( itPdgPos == corrPdgNumDone.end() )
424 {
425 // -> roll a dice if we create this particle or its correlated one
426 if ( CLHEP::RandFlat::shoot(rndmEngine) > 0.5 ) doPdg = corrPdg;
427 // now create the particles with the given pdg and note how many
428 // particles of this pdg are created
429 corrPdgNumDone[doPdg] = getAllParticles(*g4PrimaryTrack, secKinematicsMapVect, rndmEngine, doPdg, interpEnergy, interpEta);
430 }
431
432 // one of the two correlated particle types was already simulated
433 // 'itPdgPos' points to the already simulated one
434 else
435 {
436 // get the pdg of the already simulated particle and the number
437 // of these particles that were created
438 const int donePdg = itPdgPos->first;
439 const int doneNumPart = itPdgPos->second;
440 // set the pdg of the particle type that will be done
441 if (donePdg == doPdg) doPdg = corrPdg;
442
443 // now create the correlated particles
444 getCorrelatedParticles(*g4PrimaryTrack, secKinematicsMapVect, doPdg, doneNumPart, rndmEngine, interpEnergy, interpEta);
445 // note: no need to take note, that this particle type is now simulated,
446 // since this is the second of two correlated particles, which is
447 // simulated and we do not have correlations of more than two particles.
448 }
449
450 // if no correlation for this particle
451 // -> directly create all particles with the current pdg
452 }
453 else{
454 getAllParticles(*g4PrimaryTrack, secKinematicsMapVect, rndmEngine, doPdg, interpEnergy, interpEta);
455 }
456 } // for-loop over all particle pdgs
457
458 // Debugging: Stored properties secondary particles to be simulated
459 if (msgLvl(MSG::DEBUG) && !secKinematicsMapVect.empty()) {
460 for (size_t i = 0; i < secKinematicsMapVect.size(); ++i) {
461 ATH_MSG_DEBUG("Secondary[" << i << "] kinematics:");
462 for (const auto& kv : secKinematicsMapVect[i]) {
463 ATH_MSG_DEBUG(" " << kv.first << " = " << kv.second);
464 }
465 }
466 }
467
468 // increment so that we can escape the while loop at maxTries
469 nTries++;
470 }
471
472 return secKinematicsMapVect;
473}
474
475int PunchThroughG4Tool::getAllParticles(const G4Track& g4PrimaryTrack, std::vector<std::map<std::string, double>> &secKinematicsMapVect, CLHEP::HepRandomEngine* rndmEngine, int pdg, double interpEnergy, double interpEta, int numParticles)
476{
477 // Get initial Geant4 position theta and phi
478 float initParticleTheta = g4PrimaryTrack.GetPosition().theta();
479 float initParticlePhi = g4PrimaryTrack.GetPosition().phi();
480
481 // get the current particle
483 double minAllowedEnergy = p->getMinEnergy();
484
485 // if no number of particles (=-1) was handed over as an argument
486 // -> get the number of particles from the pdf
487 if ( numParticles < 0 )
488 {
489 // prepare the function arguments for the PunchThroughPDFCreator class
490 std::vector<int> parameters;
491 parameters.push_back( std::round(interpEnergy) );
492 parameters.push_back( std::round(interpEta*100) );
493 // the maximum number of particles which should be produced
494 // if no maximum number is given, this is -1
495 int maxParticles = p->getMaxNumParticles();
496
497 // get the right number of punch-through particles
498 // and ensure that we do not create too many particles
499 do
500 {
501 numParticles = int( p->getNumParticlesPDF()->getRand(rndmEngine, parameters) );
502
503 // scale the number of particles if requested
504 numParticles = lround( numParticles *= p->getNumParticlesFactor() );
505 }
506 while ( (maxParticles >= 0.) && (numParticles > maxParticles) );
507 }
508
509 ATH_MSG_DEBUG("[PunchThroughG4Tool] adding " << numParticles << " punch-through particles with pdg code: " << pdg);
510
511 // get how many secondary particles to be created
512 for ( int numCreated = 0; numCreated < numParticles; numCreated++ )
513 {
514 // create one particle which fullfills the right energy distribution
515 std::map<std::string, double> secondaryKinematicsMap = getOneParticleKinematics(rndmEngine, pdg, initParticleTheta, initParticlePhi, interpEnergy, interpEta);
516
517 //safety to cater minimum energy as rest mass
518 double energy = secondaryKinematicsMap.at("energy");
519 if(energy > minAllowedEnergy){
520 // push to the vector of maps
521 secKinematicsMapVect.push_back(secondaryKinematicsMap);
522 }
523 }
524
525 // Get the size of the vector
526 std::size_t numSecondaries = secKinematicsMapVect.size();
527
528 // the actual (pre) number of particles which was created is numSecondaries, which is incremented in the loop
529 return (numSecondaries);
530}
531
532int PunchThroughG4Tool::getCorrelatedParticles(const G4Track& g4PrimaryTrack, std::vector<std::map<std::string, double>> &secKinematicsMapVect, int pdg, int corrParticles, CLHEP::HepRandomEngine* rndmEngine, double interpEnergy, double interpEta)
533{
534 // Get Geant4 particle definition
535 const G4ParticleDefinition * mainG4Particle = g4PrimaryTrack.GetDefinition();
536 // Get Geant4 particle pdgID
537 float mainMomMag2 = g4PrimaryTrack.GetMomentum().mag2();
538 float mainPartMass = mainG4Particle->GetPDGMass();
539
540 // get the PunchThroughParticle class
542
543 const double initEnergy = std::sqrt( mainMomMag2 + mainPartMass*mainPartMass );
544
545 // (1.) decide if we do correlation or not
546 double rand = CLHEP::RandFlat::shoot(rndmEngine)
547 *(p->getFullCorrelationEnergy() - p->getMinCorrelationEnergy())
548 + p->getMinCorrelationEnergy();
549
550 // if initEnergy less than random corr energy (meaning threshold for min energy for doing correlation not passed)
551 if ( initEnergy < rand )
552 {
553 // here we do not do correlation
554 return getAllParticles(g4PrimaryTrack, secKinematicsMapVect, rndmEngine, pdg, interpEnergy, interpEta);
555 }
556
557 // (2.) if this point is reached, we do correlation
558 // decide which 2d correlation histogram to use
559 double *histDomains = p->getCorrelationHistDomains();
560 TH2F *hist2d = nullptr;
561 // compute the center values of the lowE and highE
562 // correlation histogram domains
563 if ( initEnergy < histDomains[1])
564 {
565 // initial energy lower than border between lowEnergy and highEnergy histogram domain
566 // --> choose lowEnergy correlation histogram
567 hist2d = p->getCorrelationLowEHist();
568 }
569 else
570 {
571 double rand = CLHEP::RandFlat::shoot(rndmEngine)*(histDomains[2]-histDomains[1]) + histDomains[1];
572 hist2d = ( initEnergy < rand) ? p->getCorrelationLowEHist() : p->getCorrelationHighEHist();
573 }
574
575 // get the correlation 2d histogram
576 // now find out where on the x-axis the the bin for number of
577 // correlated particles is
578 Int_t xbin = hist2d->GetXaxis()->FindFixBin(corrParticles);
579 int numParticles = 0;
580 int maxParticles = p->getMaxNumParticles();
581 // now the distribution along the y-axis is a PDF for the number
582 // of 'pdg' particles
583 do
584 {
585 double rand = CLHEP::RandFlat::shoot(rndmEngine);
586 double sum = 0.;
587 for ( int ybin = 1; ybin <= hist2d->GetNbinsY(); ybin++ )
588 {
589 sum += hist2d->GetBinContent(xbin, ybin);
590 // check if we choose the current bin or not
591 if ( sum >= rand )
592 {
593 numParticles = ybin - 1;
594 break;
595 }
596 }
597 // scale the number of particles is requested
598 numParticles = lround( numParticles * p->getNumParticlesFactor() );
599 }
600 while ( (maxParticles >= 0.) && (numParticles > maxParticles) );
601
602 // finally create this exact number of particles
603 return getAllParticles(g4PrimaryTrack, secKinematicsMapVect, rndmEngine, pdg, interpEnergy, interpEta, numParticles);
604}
605
606std::map<std::string, double> PunchThroughG4Tool::getOneParticleKinematics(CLHEP::HepRandomEngine* rndmEngine, int secondaryPDG, float initParticleTheta, float initParticlePhi, double interpEnergy, double interpEta) const
607{
608 // get a local copy of the needed punch-through particle class
609 PunchThroughParticle *p = m_particles.at(secondaryPDG);
610
611 // (0.) get the pca / cdf group based on pdgId and eta, eta times 100, e.g eta -4 to 4 is from eta -400 to 400
612 int pcaCdfIterator = passedParamIterator(secondaryPDG, interpEta*100, m_xml_info_pca); //pca and cdf info should be of same size
613
614 ATH_MSG_DEBUG("[PunchThroughG4Tool] passedPCAIterator ==> pcaCdfIterator = "<< pcaCdfIterator <<" , pdg = "<< secondaryPDG <<" , interpEnergy = "<< interpEnergy <<" MeV, interpEta(*100) = "<< interpEta*100);
615
616 // (1.) decide if we create a particle or an anti-particle,
617 // this is from m_doAntiParticles properties
618 int anti = 1;
619 if ( p->getdoAnti() )
620 {
621 // get a random-value
622 double rand = CLHEP::RandFlat::shoot(rndmEngine);
623 // 50/50 chance to be a particle or its anti-particle
624 if (rand > 0.5) anti = -1;
625 }
626
627 // (2.) get the right punch-through distributions
628 // prepare the function arguments for the PunchThroughPDFCreator class
629 std::vector<int> parInitEnergyEta;
630 parInitEnergyEta.push_back( std::round(interpEnergy) );
631 parInitEnergyEta.push_back( std::round(interpEta*100) );
632
633 //initialise variables to store punch through particle kinematics
634 double energy = 0.;
635 double deltaTheta = 0.;
636 double deltaPhi = 0.;
637 double momDeltaTheta = 0.;
638 double momDeltaPhi = 0.;
639
640 double principal_component_0 = 0.;
641 double principal_component_1 = 0.;
642 double principal_component_2 = 0.;
643 double principal_component_3 = 0.;
644 double principal_component_4 = 0.;
645 std::vector<double> transformed_variables;
646
647 principal_component_0 = p->getPCA0PDF()->getRand(rndmEngine, parInitEnergyEta);
648 principal_component_1 = p->getPCA1PDF()->getRand(rndmEngine, parInitEnergyEta);
649 principal_component_2 = p->getPCA2PDF()->getRand(rndmEngine, parInitEnergyEta);
650 principal_component_3 = p->getPCA3PDF()->getRand(rndmEngine, parInitEnergyEta);
651 principal_component_4 = p->getPCA4PDF()->getRand(rndmEngine, parInitEnergyEta);
652
653 ATH_MSG_DEBUG("[PunchThroughG4Tool] Drawn punch through kinematics PCA components: PCA0 = "<< principal_component_0 <<" PCA1 = "<< principal_component_1 <<" PCA2 = "<< principal_component_2 <<" PCA3 = "<< principal_component_3 <<" PCA4 = "<< principal_component_4 );
654
655 std::vector<double> principal_components {
656 principal_component_0,
657 principal_component_1,
658 principal_component_2,
659 principal_component_3,
660 principal_component_4
661 };
662
663 transformed_variables = inversePCA(pcaCdfIterator,principal_components);
664
665 energy = inverseCdfTransform(transformed_variables.at(0), m_variable0_inverse_cdf[pcaCdfIterator]);
666 deltaTheta = inverseCdfTransform(transformed_variables.at(1), m_variable1_inverse_cdf[pcaCdfIterator]);
667 deltaPhi = inverseCdfTransform(transformed_variables.at(2), m_variable2_inverse_cdf[pcaCdfIterator]);
668 momDeltaTheta = inverseCdfTransform(transformed_variables.at(3), m_variable3_inverse_cdf[pcaCdfIterator]);
669 momDeltaPhi = inverseCdfTransform(transformed_variables.at(4), m_variable4_inverse_cdf[pcaCdfIterator]);
670
671 ATH_MSG_DEBUG("[PunchThroughG4Tool] Transformed punch through kinematics: energy = "<< energy <<" MeV deltaTheta = "<< deltaTheta <<" deltaPhi = "<< deltaPhi <<" momDeltaTheta = "<< momDeltaTheta <<" momDeltaPhi = "<< momDeltaPhi );
672
673 energy *= p->getEnergyFactor(); // scale the energy if requested
674
675 // safety: if energy less than minimum energy, add a small energy (10MeV to it to continue)
676 if (energy < p->getMinEnergy()) {
677 energy = p->getMinEnergy() + 10;
678 }
679
680 // (2.2) get the particles DeltaTheta relative to the incoming particle
681 double theta = 0;
682 // loop to keep theta within range [0,PI]
683 do
684 {
685 // decide if delta positive/negative
686 deltaTheta *= ( CLHEP::RandFlat::shoot(rndmEngine) > 0.5 ) ? 1. : -1.;
687
688 // calculate the exact theta value of the later created
689 // punch-through particle
690 theta = initParticleTheta + deltaTheta*p->getPosAngleFactor();
691 }
692 while ( (theta > M_PI) || (theta < 0.) );
693
694 // (2.3) get the particle's delta phi relative to the incoming particle
695 deltaPhi *= ( CLHEP::RandFlat::shoot(rndmEngine) > 0.5 ) ? 1. : -1.;
696
697 // keep phi within range [-PI,PI]
698 double phi = initParticlePhi + deltaPhi*p->getPosAngleFactor();
699 while ( std::fabs(phi) > 2*M_PI) phi /= 2.;
700 while (phi > M_PI) phi -= 2*M_PI;
701 while (phi < -M_PI) phi += 2*M_PI;
702
703 // (2.4) get the particle momentum delta theta, relative to its position
704 //
705 // loop to keep momTheta within range [0,PI]
706 double momTheta = 0.;
707 do
708 {
709 // decide if delta positive/negative
710 momDeltaTheta *= ( CLHEP::RandFlat::shoot(rndmEngine) > 0.5 ) ? 1. : -1.;
711 // calculate the exact momentum theta value of the later created
712 // punch-through particle
713 momTheta = theta + momDeltaTheta*p->getMomAngleFactor();
714 }
715 while ( (momTheta > M_PI) || (momTheta < 0.) );
716
717 // (2.5) get the particle momentum delta phi, relative to its position
718 momDeltaPhi *= ( CLHEP::RandFlat::shoot(rndmEngine) > 0.5 ) ? 1. : -1.;
719
720 double momPhi = phi + momDeltaPhi*p->getMomAngleFactor();
721 // keep momPhi within range [-PI,PI]
722 while ( std::fabs(momPhi) > 2*M_PI) momPhi /= 2.;
723 while (momPhi > M_PI) momPhi -= 2*M_PI;
724 while (momPhi < -M_PI) momPhi += 2*M_PI;
725
726 // store the kinematics in a map and return it
727 std::map<std::string, double> secondaryKinematicsMap;
728 secondaryKinematicsMap.insert({ "anti" , anti });
729 secondaryKinematicsMap.insert({ "secondaryPDG" , secondaryPDG });
730 secondaryKinematicsMap.insert({ "energy" , energy });
731 secondaryKinematicsMap.insert({ "theta" , theta });
732 secondaryKinematicsMap.insert({ "phi" , phi });
733 secondaryKinematicsMap.insert({ "momTheta" , momTheta });
734 secondaryKinematicsMap.insert({ "momPhi" , momPhi });
735
736 return secondaryKinematicsMap;
737}
738
739std::vector<std::map<std::string, double>> PunchThroughG4Tool::checkEnergySumFromSecondaries(double mainEnergyInit, std::vector<std::map<std::string, double>> &secKinematicsMapVect){
740 //initialize variables
741 double energy;
742 double totEnergySecondaries = 0;
743
744 // Check energy conservation
745 // Step 1: sum all of the energies from the vector of secondaries (map) together, to initialize procedure later
746 for (std::size_t i = 0; i < secKinematicsMapVect.size(); i++) {
747 energy = secKinematicsMapVect[i].at("energy");
748 totEnergySecondaries += energy;
749 }
750
751 // Step 2: Sort the vector based on "energy" in descending order using lambda expression
752 std::sort(
753 secKinematicsMapVect.begin(),
754 secKinematicsMapVect.end(),
755 [](const std::map<std::string, double>& lhs, const std::map<std::string, double>& rhs) {
756 return lhs.at("energy") > rhs.at("energy");
757 }
758 );
759
760 // Just for printing here
761 if(totEnergySecondaries > mainEnergyInit){
762 ATH_MSG_DEBUG("[PunchThroughG4Tool::checkEnergySumFromSecondaries] Case where energy of created secondaries more than parent track identified! ");
763 ATH_MSG_DEBUG("[PunchThroughG4Tool::checkEnergySumFromSecondaries] ==> TotalSecondariesEnergy = " << totEnergySecondaries << ", ParentEnergy = "<< mainEnergyInit);
764 }
765
766 // Step 3: Adjust the energy and remove secondary with lowest energy to ensure total energy of secondaries <= mainEnergyInit
767 while (totEnergySecondaries > mainEnergyInit && !secKinematicsMapVect.empty()) {
768 // Get the last map in the vector (which has the lowest energy due to sorting)
769 std::map<std::string, double> lastMap = secKinematicsMapVect.back();
770
771 // Get the energy value of the last map
772 double lastEnergy = lastMap.at("energy");
773
774 // Remove the last map from the vector
775 secKinematicsMapVect.pop_back();
776
777 // Subtract the removed energy from totEnergySecondaries
778 totEnergySecondaries -= lastEnergy;
779 }
780
781 //finally return the secKinematicsMapVect after the checks
782 return secKinematicsMapVect;
783}
784
785void PunchThroughG4Tool::createAllSecondaryTracks(G4ParticleTable &ptable, G4FastStep& fastStep, const G4Track& g4PrimaryTrack, std::vector<std::map<std::string, double>> &secKinematicsMapVect, G4TrackVector& secTrackCont, const std::vector<double> &caloMSVars){
786 // Get Geant4 particle definition
787 const G4ParticleDefinition* mainG4Particle = g4PrimaryTrack.GetDefinition();
788 float mainMomMag2 = g4PrimaryTrack.GetMomentum().mag2();
789 float mainPartMass = mainG4Particle->GetPDGMass();
790 double mainEnergyInit = std::sqrt(mainMomMag2 + mainPartMass * mainPartMass);
791
792 // Energy conservation: call checkEnergySumFromSecondaries to check/update secKinematicsMapVect
793 secKinematicsMapVect = checkEnergySumFromSecondaries(mainEnergyInit, secKinematicsMapVect);
794
795 // set number of secondary tracks to be produced
796 int numSecondaries = secKinematicsMapVect.size();
797
798 if(numSecondaries>0){ //safety
799 // set number of secondary tracks from numSecondaries above
800 fastStep.SetNumberOfSecondaryTracks(numSecondaries);
801
802 // Get current (global) time from original main G4Track
803 double currentTime = g4PrimaryTrack.GetGlobalTime();
804
805 // create the secondaries given the number of secondaries to be created
806 for(int i=0;i < numSecondaries; i++){
807 int anti = (int)(secKinematicsMapVect[i].at("anti"));
808 int secondaryPDG = (int)(secKinematicsMapVect[i].at("secondaryPDG"));
809 int signedPDG = secondaryPDG*anti;
810 double energy = secKinematicsMapVect[i].at("energy");
811 double theta = secKinematicsMapVect[i].at("theta");
812 double phi = secKinematicsMapVect[i].at("phi");
813 double momTheta = secKinematicsMapVect[i].at("momTheta");
814 double momPhi = secKinematicsMapVect[i].at("momPhi");
815
816 // (**) finally create the punch-through particle as a G4Track (secondary particle track)
817 ATH_MSG_DEBUG("[PunchThroughG4Tool::createAllSecondaryTracks] createSecondaryTrack input parameters: currentTime = " << currentTime << " anti? = "<< anti <<" signedPDG = "<< signedPDG <<" energy = "<< energy <<" theta = "<< theta <<" phi = "<< phi <<" momTheta = "<< momTheta << " momPhi " << momPhi);
818 G4Track* newSecTrack = createSecondaryTrack( ptable, fastStep, currentTime, signedPDG, energy, theta, phi, momTheta, momPhi, caloMSVars);
819
820 // if something went wrong and the secondary track is simply not created
821 if (!newSecTrack)
822 {
823 G4Exception("[PunchThroughG4Tool::createAllSecondaryTracks]", "ExceptionError", FatalException, "something went wrong while creating punch-through particle tracks");
824 return;
825 }
826
827 // add this secondary track to the g4trackvector
828 secTrackCont.push_back( newSecTrack );
829 }
830 }
831
832 //printout: fastStep.DumpInfo()
833}
834
835G4Track* PunchThroughG4Tool::createSecondaryTrack( G4ParticleTable &ptable, G4FastStep& fastStep, double currentTime, int secondarySignedPDG,
836 double energy, double theta, double phi,double momTheta, double momPhi, const std::vector<double> &caloMSVars)
837{
838 //initialize to nullptr
839 G4Track *newSecTrack = nullptr;
840
841 G4ParticleDefinition* secG4Particle = ptable.FindParticle(secondarySignedPDG);
842 double mass = secG4Particle->GetPDGMass();
843
844 // the intersection point with Calo-MS surface
845 G4ThreeVector posVec = punchTroughPosPropagator(theta, phi, caloMSVars[0], caloMSVars[1], caloMSVars[2], caloMSVars[3]);
846
847 // set up the real punch-through particle at this position
848 // set up the momentum vector of this particle as a GlobalMomentum
849 // by using the given energy and mass of the particle and also using
850 // the given theta and phi
851 G4ThreeVector momVec;
852 double momMag = std::sqrt(energy*energy - mass*mass);
853 momVec.setRThetaPhi(momMag,momTheta,momPhi);
854 ATH_MSG_DEBUG("[PunchThroughG4Tool]::createSecondaryTrack] setRThetaPhi pre input parameters: energy = " << energy <<" mass = "<< mass);
855 ATH_MSG_DEBUG("[PunchThroughG4Tool]::createSecondaryTrack] setRThetaPhi input parameters: std::sqrt(energy*energy - mass*mass) = " << std::sqrt(energy*energy - mass*mass) << " momTheta = "<< momTheta <<" momPhi = "<< momPhi);
856
857 // create dynamic particle
858 G4DynamicParticle dynParticle(secG4Particle,momVec);
859 newSecTrack = fastStep.CreateSecondaryTrack(dynParticle, posVec, currentTime, false); // use position in global coordinates instead
860
861 return newSecTrack;
862}
863
865 return 0.5 * TMath::Erfc(-x * M_SQRT1_2);
866}
867
868std::vector<double> PunchThroughG4Tool::dotProduct(const std::vector<std::vector<double>> &m, const std::vector<double> &v)
869{
870 std::vector<double> result;
871 result.reserve(m.size());
872 for (const auto& r : m){
873 result.push_back(std::inner_product(v.begin(), v.end(), r.begin(), 0.0));
874 }
875
876 return result;
877}
878
879template<typename T>
880std::vector<T> str_to_list(const std::string_view str)
881{
882 constexpr char delimiters = ',';
883 std::vector<T> tokens;
884 // Skip delimiters at beginning.
885 std::string_view::size_type lastPos = str.find_first_not_of(delimiters, 0);
886 // Find first "non-delimiter".
887 std::string_view::size_type pos = str.find_first_of(delimiters, lastPos);
888
889 while (std::string_view::npos != pos || std::string_view::npos != lastPos) {
890 // Found a token, add it to the vector.
891 std::string_view numbStr = str.substr(lastPos, pos - lastPos);
892 T num = -9999;
893 std::from_chars(numbStr.data(), numbStr.data() + numbStr.size(), num);
894 tokens.push_back(num);
895 // Skip delimiters. Note the "not_of"
896 lastPos = str.find_first_not_of(delimiters, pos);
897 // Find next "non-delimiter"
898 pos = str.find_first_of(delimiters, lastPos);
899 }
900 return tokens;
901}
902
903int PunchThroughG4Tool::passedParamIterator(int pid, double eta, const std::vector<std::map<std::string, std::string>> &mapvect) const
904{
905 //convert the pid to absolute value and string for query
906 int pidStrSingle = std::abs(pid);
907 //STEP 1
908 //filter items matching pid first
909
910 for (unsigned int i = 0; i < mapvect.size(); i++){
911 const std::string &pidStr = mapvect[i].at("pidStr");
912 auto v = str_to_list<int>(pidStr);
913 if(std::find(v.begin(), v.end(),pidStrSingle)==v.end()) continue;
914 const std::string &etaMinsStr = mapvect[i].at("etaMins");
915 const std::string &etaMaxsStr = mapvect[i].at("etaMaxs");
916 std::vector<double> etaMinsVect = str_to_list<double>(etaMinsStr);
917 std::vector<double> etaMaxsVect = str_to_list<double>(etaMaxsStr);
918 assert(etaMaxsVect.size() == etaMinsVect.size());
919 for (unsigned int j = 0; j < etaMinsVect.size(); j++){ // assume size etaMinsVect == etaMaxsVect
920 double etaMinToCompare = etaMinsVect[j];
921 double etaMaxToCompare = etaMaxsVect[j];
922 if((eta >= etaMinToCompare) && (eta < etaMaxToCompare)){
923 //PASS CONDITION
924 //then choose the passing one and note it's iterator
925 return (i); //in case more than 1 match (ambiguous case)
926 }
927 }
928 }
929 return 0;
930}
931
932std::vector<std::map<std::string, std::string>> PunchThroughG4Tool::getInfoMap(const std::string& mainNode, const std::string &xmlFilePath){
933 // Initialize pointers
934 xmlDocPtr doc;
935
936 std::vector<std::map<std::string, std::string>> xml_info;
937 doc = xmlParseFile( xmlFilePath.c_str() );
938
939 //check info first
940 for( xmlNodePtr nodeRoot = doc->children; nodeRoot != nullptr; nodeRoot = nodeRoot->next) {
941 if (xmlStrEqual( nodeRoot->name, BAD_CAST mainNode.c_str() )) {
942 for( xmlNodePtr nodeRootChild = nodeRoot->children; nodeRootChild != nullptr; nodeRootChild = nodeRootChild->next ) {
943 if (xmlStrEqual( nodeRootChild->name, BAD_CAST "info" )) {
944 if (nodeRootChild->children != NULL) {
945 for( xmlNodePtr infoNode = nodeRootChild->children; infoNode != nullptr; infoNode = infoNode->next) {
946 if(xmlStrEqual( infoNode->name, BAD_CAST "item" )){
947 std::map<std::string, std::string> xml_info_item;
952 xml_info.emplace_back(std::move(xml_info_item));
953 }
954 }
955 }
956 }
957 }
958 }
959 }
960
961 // free memory when done
962 xmlFreeDoc(doc);
963
964 return xml_info;
965}
966
967std::vector<double> PunchThroughG4Tool::inversePCA(int pcaCdfIterator, std::vector<double> &variables) const
968{
969 std::vector<double> transformed_variables = dotProduct(m_inverse_PCA_matrix[pcaCdfIterator], variables);
970
971 std::transform (transformed_variables.begin(), transformed_variables.end(), m_PCA_means[pcaCdfIterator].begin(), transformed_variables.begin(), std::plus<double>()); // + means
972
973 return transformed_variables;
974}
975
976StatusCode PunchThroughG4Tool::initializeInversePCA(const std::string & inversePCAConfigFile){
977 // Initialize pointers
978 xmlDocPtr doc;
979
980 doc = xmlParseFile( inversePCAConfigFile.c_str() );
981
982 ATH_MSG_DEBUG( "[PunchThroughG4Tool] Loading inversePCA: " << inversePCAConfigFile);
983
984 //check info first
985 m_xml_info_pca = getInfoMap("PCAinverse",inversePCAConfigFile);
986
987 //do the saving
988 for (unsigned int i = 0; i < m_xml_info_pca.size(); i++) {
989 std::vector<std::vector<double>> PCA_matrix;
990 ATH_MSG_DEBUG( "[PunchThroughG4Tool] m_xml_info_pca[" << i << "].at('name') = " << m_xml_info_pca[i].at("name"));
991
992 for( xmlNodePtr nodeRoot = doc->children; nodeRoot != nullptr; nodeRoot = nodeRoot->next) {
993 if (xmlStrEqual( nodeRoot->name, BAD_CAST "PCAinverse" )) {
994 for( xmlNodePtr nodePCAinverse = nodeRoot->children; nodePCAinverse != nullptr; nodePCAinverse = nodePCAinverse->next ) {
995
996 if (xmlStrEqual( nodePCAinverse->name, BAD_CAST m_xml_info_pca[i].at("name").c_str() )) {
997 if (nodePCAinverse->children != NULL) {
998 for( xmlNodePtr pcaNode = nodePCAinverse->children; pcaNode != nullptr; pcaNode = pcaNode->next) {
999
1000 if (xmlStrEqual( pcaNode->name, BAD_CAST "PCAmatrix" )) {
1001 std::vector<double> PCA_matrix_row;
1002 for (int i = 0; i <= 4; ++i) {
1003 std::string propName = "comp_" + std::to_string(i); // Dynamically create property name
1004 AddXmlToCollection<double, std::vector<double> >(pcaNode, propName.c_str(), PCA_matrix_row);
1005 }
1006 PCA_matrix.push_back(std::move(PCA_matrix_row));
1007 }
1008 else if (xmlStrEqual( pcaNode->name, BAD_CAST "PCAmeans" )) {
1009 std::vector<double> PCA_means_row;
1010 for (int i = 0; i <= 4; ++i) {
1011 std::string propName = "mean_" + std::to_string(i); // Dynamically create property name
1012 AddXmlToCollection<double, std::vector<double> >(pcaNode, propName.c_str(), PCA_means_row);
1013 }
1014 m_PCA_means.push_back(std::move(PCA_means_row));
1015 }
1016
1017 }
1018
1019 }
1020 }
1021
1022 }
1023 }
1024 }
1025 m_inverse_PCA_matrix.push_back(PCA_matrix);
1026 }
1027
1028 // free memory when done
1029 xmlFreeDoc(doc);
1030
1031 return StatusCode::SUCCESS;
1032}
1033
1034StatusCode PunchThroughG4Tool::initializeInverseCDF(const std::string & inverseCdfConfigFile){
1035 std::map<double, double> variable0_inverse_cdf_row;
1036 std::map<double, double> variable1_inverse_cdf_row;
1037 std::map<double, double> variable2_inverse_cdf_row;
1038 std::map<double, double> variable3_inverse_cdf_row;
1039 std::map<double, double> variable4_inverse_cdf_row;
1040
1041 // Initialize pointers
1042 xmlDocPtr doc;
1043
1044 //parse xml that contains config for inverse CDF for each of punch through particle kinematics
1045
1046 doc = xmlParseFile( inverseCdfConfigFile.c_str() );
1047
1048 ATH_MSG_DEBUG( "[PunchThroughG4Tool] Loading inverse CDF: " << inverseCdfConfigFile);
1049
1050 //check info first
1051 m_xml_info_cdf = getInfoMap("CDFMappings",inverseCdfConfigFile);
1052
1053 //do the saving
1054 for (unsigned int i = 0; i < m_xml_info_cdf.size(); i++) {
1055 ATH_MSG_DEBUG( "[PunchThroughG4Tool] m_xml_info_cdf[" << i << "].at('name') = " << m_xml_info_cdf[i].at("name"));
1056
1057 for( xmlNodePtr nodeRoot = doc->children; nodeRoot != nullptr; nodeRoot = nodeRoot->next) {
1058 if (xmlStrEqual( nodeRoot->name, BAD_CAST "CDFMappings" )) {
1059 for( xmlNodePtr typeMappings = nodeRoot->children; typeMappings != nullptr; typeMappings = typeMappings->next ) {
1060 if (xmlStrEqual( typeMappings->name, BAD_CAST m_xml_info_cdf[i].at("name").c_str() )) {
1061 if (typeMappings->children != NULL) {
1062 for( xmlNodePtr nodeMappings = typeMappings->children; nodeMappings != nullptr; nodeMappings = nodeMappings->next) {
1063
1064 if (xmlStrEqual( nodeMappings->name, BAD_CAST "variable0" )) {
1065 variable0_inverse_cdf_row = getVariableCDFmappings(nodeMappings);
1066 }
1067 else if (xmlStrEqual( nodeMappings->name, BAD_CAST "variable1" )) {
1068 variable1_inverse_cdf_row = getVariableCDFmappings(nodeMappings);
1069 }
1070 else if (xmlStrEqual( nodeMappings->name, BAD_CAST "variable2" )) {
1071 variable2_inverse_cdf_row = getVariableCDFmappings(nodeMappings);
1072 }
1073 else if (xmlStrEqual( nodeMappings->name, BAD_CAST "variable3" )) {
1074 variable3_inverse_cdf_row = getVariableCDFmappings(nodeMappings);
1075 }
1076 else if (xmlStrEqual( nodeMappings->name, BAD_CAST "variable4" )) {
1077 variable4_inverse_cdf_row = getVariableCDFmappings(nodeMappings);
1078 }
1079 }
1080
1081 }
1082 }
1083 }
1084 }
1085 }
1086 m_variable0_inverse_cdf.push_back(variable0_inverse_cdf_row);
1087 m_variable1_inverse_cdf.push_back(variable1_inverse_cdf_row);
1088 m_variable2_inverse_cdf.push_back(variable2_inverse_cdf_row);
1089 m_variable3_inverse_cdf.push_back(variable3_inverse_cdf_row);
1090 m_variable4_inverse_cdf.push_back(variable4_inverse_cdf_row);
1091 }
1092
1093 // free memory when done
1094 xmlFreeDoc(doc);
1095
1096 return StatusCode::SUCCESS;
1097}
1098
1099std::map<double, double> PunchThroughG4Tool::getVariableCDFmappings(xmlNodePtr& nodeParent){
1100 std::map<double, double> mappings;
1101 double ref = -1;
1102 double quant = -1;
1103 for( xmlNodePtr node = nodeParent->children; node != nullptr; node = node->next ) {
1104 //Get min and max values that we normalise values to
1105 if (xmlStrEqual( node->name, BAD_CAST "CDFmap" )) {
1106 ref = GetXmlAttr(node, "ref", -1);
1107 quant = GetXmlAttr(node, "quant", -1);
1108 mappings.emplace(ref, quant );
1109 }
1110 }
1111
1112 return mappings;
1113}
1114
1115double PunchThroughG4Tool::inverseCdfTransform(double variable, const std::map<double, double> &inverse_cdf_map) {
1116
1117 double norm_cdf = normal_cdf(variable);
1118
1119 auto upper = inverse_cdf_map.upper_bound(norm_cdf);
1120 auto lower = upper--;
1121
1122 double m = (upper->second - lower->second)/(upper->first - lower->first);
1123 double c = lower->second - m * lower->first;
1124 double transformed = m * norm_cdf + c;
1125
1126 return transformed;
1127
1128}
1129
1130double PunchThroughG4Tool::interpolateEnergy(const double &energy, CLHEP::HepRandomEngine* rndmEngine) const{
1131
1132 ATH_MSG_DEBUG("[PunchThroughG4Tool] interpolating incoming energy: " << energy);
1133
1134 std::string energyPointsString;
1135 for (auto element:m_energyPoints){
1136 energyPointsString += std::to_string(element) + " ";
1137 }
1138
1139 ATH_MSG_DEBUG("[PunchThroughG4Tool] available energy points: " << energyPointsString);
1140
1141 auto const upperEnergy = std::upper_bound(m_energyPoints.begin(), m_energyPoints.end(), energy);
1142
1143 if(upperEnergy == m_etaPoints.end()){ //if no energy greater than input energy, choose greatest energy
1144 ATH_MSG_DEBUG("[PunchThroughG4Tool] incoming energy > largest energy point, returning greatest energy point: " << m_energyPoints.back());
1145 return m_energyPoints.back();
1146 }
1147 else if(upperEnergy == m_etaPoints.begin()){ //if smallest energy greater than input energy, choose smallest energy
1148 ATH_MSG_DEBUG("[PunchThroughG4Tool] incoming energy < smallest energy point, returning smallest energy point: " << *upperEnergy);
1149 return *upperEnergy;
1150 }
1151
1152 ATH_MSG_DEBUG("[PunchThroughG4Tool] energy points upper_bound: "<< *upperEnergy);
1153
1154 double randomShoot = CLHEP::RandFlat::shoot(rndmEngine);
1155
1156 ATH_MSG_DEBUG("[PunchThroughG4Tool] Shooting random number: "<< randomShoot);
1157
1158 double midPoint = *std::prev(upperEnergy)*M_SQRT2;
1159
1160 if(energy < midPoint){ //if energy smaller than mid point in log(energy)
1161
1162 double distance = std::abs(energy - *std::prev(upperEnergy))/((midPoint) - *std::prev(upperEnergy));
1163
1164 ATH_MSG_DEBUG( "[PunchThroughG4Tool] incoming energy is closest to prev(upper_bound) in log(energy), distance: " << distance );
1165
1166 if(randomShoot < distance){
1167 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot < distance, returning upper_bound " << *upperEnergy );
1168 return *upperEnergy;
1169 }
1170 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot > distance, returning prev(upper_bound) " << *std::prev(upperEnergy) );
1171
1172 return *std::prev(upperEnergy);
1173 }
1174 else if(energy > midPoint){ //if energy greater than mid point in log(energy)
1175
1176 double distance = std::abs(energy - *upperEnergy)/((*upperEnergy - midPoint));
1177
1178 ATH_MSG_DEBUG( "[PunchThroughG4Tool] incoming energy is closest to upper_bound in log(energy), distance: " << distance );
1179
1180 if(randomShoot < distance){
1181 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot < distance, returning prev(upper_bound) " << *std::prev(upperEnergy) );
1182 return *std::prev(upperEnergy);
1183 }
1184 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot > distance, returning upper_bound " << *upperEnergy );
1185 return *upperEnergy;
1186 }
1187
1188 return *upperEnergy;
1189}
1190
1191double PunchThroughG4Tool::interpolateEta(const double &eta, CLHEP::HepRandomEngine* rndmEngine) const{
1192
1193 double absEta = std::abs(eta);
1194
1195 ATH_MSG_DEBUG("[PunchThroughG4Tool] interpolating incoming abs(eta): " << absEta);
1196
1197 std::string etaPointsString;
1198 for (auto element:m_etaPoints){
1199 etaPointsString += std::to_string(element) + " ";
1200 }
1201
1202 ATH_MSG_DEBUG("[PunchThroughG4Tool] available eta points: " << etaPointsString);
1203
1204 auto const upperEta = std::upper_bound(m_etaPoints.begin(), m_etaPoints.end(), absEta);
1205
1206 if(upperEta == m_etaPoints.end()){
1207 ATH_MSG_DEBUG("[PunchThroughG4Tool] incoming abs(eta) > largest eta point, returning greatest eta point: " << m_etaPoints.back());
1208 return m_etaPoints.back();
1209 }
1210
1211
1212 ATH_MSG_DEBUG("[PunchThroughG4Tool] eta points upper_bound: "<< *upperEta);
1213
1214 double randomShoot = CLHEP::RandFlat::shoot(rndmEngine);
1215
1216 ATH_MSG_DEBUG("[PunchThroughG4Tool] Shooting random number: "<< randomShoot);
1217
1218 if(std::abs(absEta - *upperEta) < std::abs(absEta - *std::prev(upperEta))){
1219
1220 double distance = std::abs(absEta - *upperEta)/((*upperEta - *std::prev(upperEta))/2);
1221
1222 ATH_MSG_DEBUG( "[PunchThroughG4Tool] abs(eta) is closer to eta points upper_bound, distance: " << distance );
1223
1224 if(randomShoot > distance){
1225 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot > distance, returning upper_bound " << *upperEta );
1226 return *upperEta;
1227 }
1228
1229 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot < distance, returning prev(upper_bound) " << *std::prev(upperEta) );
1230
1231 return *std::prev(upperEta);
1232 }
1233 else if(std::abs(absEta - *std::prev(upperEta)) < std::abs(absEta - *upperEta)){
1234
1235 if(std::prev(upperEta) == m_etaPoints.begin()){
1236 ATH_MSG_DEBUG( "[PunchThroughG4Tool] prev of upper bound is begin, returning that: " << *std::prev(upperEta) );
1237 return *std::prev(upperEta);
1238 }
1239
1240 double distance = std::abs(absEta - *std::prev(upperEta))/((*std::prev(upperEta) - *std::prev(std::prev(upperEta)))/2);
1241
1242 ATH_MSG_DEBUG( "[PunchThroughG4Tool] abs(eta) is closer to eta points prev(upper_bound), distance: " << distance );
1243
1244 if(randomShoot > distance){
1245 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot > distance, returning prev(prev(upper_bound)) " << *std::prev(std::prev(upperEta)) );
1246
1247 return *std::prev(std::prev(upperEta));
1248 }
1249 ATH_MSG_DEBUG( "[PunchThroughG4Tool] randomShoot < distance, returning prev(upper_bound) " << *std::prev(upperEta) );
1250
1251 return *std::prev(upperEta);
1252 }
1253
1254 return *std::prev(upperEta);
1255}
1256
1257StatusCode
1258PunchThroughG4Tool::registerPunchThroughParticle(G4ParticleTable &ptable, int pdg, bool doAntiparticle,
1259 double minEnergy, int maxNumParticles, double numParticlesFactor,
1260 double energyFactor, double posAngleFactor, double momAngleFactor)
1261{
1262 G4ParticleDefinition* secG4Particle = ptable.FindParticle(std::abs(pdg));
1263 double restMass = secG4Particle->GetPDGMass();
1264
1265 // read in the data needed to construct the distributions for the number of punch-through particles
1266 // (1.) get the distribution function for the number of punch-through particles
1267 std::unique_ptr<PunchThroughPDFCreator> pdf_num(readLookuptablePDF(pdg, m_fileLookupTable, "FREQ_PDG"));
1268 if (!pdf_num ) return StatusCode::FAILURE; // return error if something went wrong
1269
1270 // (2.) get the PDF for the punch-through energy
1271 std::unique_ptr<PunchThroughPDFCreator> pdf_pca0 (readLookuptablePDF(pdg, m_fileLookupTable, "PCA0_PDG"));
1272 if (!pdf_pca0)
1273 {
1274 return StatusCode::FAILURE; // return error if something went wrong
1275 }
1276
1277 // (3.) get the PDF for the punch-through particles difference in
1278 // theta compared to the incoming particle
1279 std::unique_ptr<PunchThroughPDFCreator> pdf_pca1 (readLookuptablePDF(pdg, m_fileLookupTable, "PCA1_PDG"));
1280 if (!pdf_pca1)
1281 {
1282 return StatusCode::FAILURE;
1283 }
1284
1285 // (4.) get the PDF for the punch-through particles difference in
1286 // phi compared to the incoming particle
1287 std::unique_ptr<PunchThroughPDFCreator> pdf_pca2 (readLookuptablePDF(pdg, m_fileLookupTable, "PCA2_PDG"));
1288 if (!pdf_pca2)
1289 {
1290 return StatusCode::FAILURE;
1291 }
1292
1293 // (5.) get the PDF for the punch-through particle momentum delta theta angle
1294 std::unique_ptr<PunchThroughPDFCreator> pdf_pca3 (readLookuptablePDF(pdg, m_fileLookupTable, "PCA3_PDG"));
1295 if (!pdf_pca3)
1296 {
1297 return StatusCode::FAILURE;
1298 }
1299
1300 // (6.) get the PDF for the punch-through particle momentum delta phi angle
1301 std::unique_ptr<PunchThroughPDFCreator> pdf_pca4 (readLookuptablePDF(pdg, m_fileLookupTable, "PCA4_PDG"));
1302 if (!pdf_pca4)
1303 {
1304 return StatusCode::FAILURE;
1305 }
1306
1307 // (7.) now finally store all this in the right std::map
1308 PunchThroughParticle *particle = new PunchThroughParticle(pdg, doAntiparticle);
1309 particle->setNumParticlesPDF(std::move(pdf_num));
1310 particle->setPCA0PDF(std::move(pdf_pca0));
1311 particle->setPCA1PDF(std::move(pdf_pca1));
1312 particle->setPCA2PDF(std::move(pdf_pca2));
1313 particle->setPCA3PDF(std::move(pdf_pca3));
1314 particle->setPCA4PDF(std::move(pdf_pca4));
1315
1316 // (8.) set some additional particle and simulation properties
1317 minEnergy = ( minEnergy > restMass ) ? minEnergy : restMass;
1318 particle->setMinEnergy(minEnergy);
1319 particle->setMaxNumParticles(maxNumParticles);
1320 particle->setNumParticlesFactor(numParticlesFactor);
1321 particle->setEnergyFactor(energyFactor);
1322 particle->setPosAngleFactor(posAngleFactor);
1323 particle->setMomAngleFactor(momAngleFactor);
1324
1325 // (9.) insert this PunchThroughParticle instance into the std::map class member
1326 m_particles[pdg] = particle;
1327
1328 return StatusCode::SUCCESS;
1329}
1330
1331StatusCode PunchThroughG4Tool::registerCorrelation(int pdgID1, int pdgID2,
1332 double minCorrEnergy, double fullCorrEnergy)
1333{
1334 // find the given pdgs in the registered particle ids
1335 std::map<int, PunchThroughParticle*>::iterator location1 = m_particles.find(pdgID1);
1336 std::map<int, PunchThroughParticle*>::iterator location2 = m_particles.find(pdgID2);
1337
1338 // if at least one of the given pdgs was not registered yet -> return an error
1339 if ( (location1 == m_particles.end()) || (location2 == m_particles.end()) )
1340 return StatusCode::FAILURE;
1341
1342 // now look for the correlation histograms
1343 std::stringstream name;
1344 name << "NumExitCorrelations/x_PDG" << std::abs(pdgID1) << "__y_PDG" << std::abs(pdgID2) << "__lowE";
1345 TH2F *hist1_2_lowE = (TH2F*)m_fileLookupTable->Get(name.str().c_str());
1346 name.str("");
1347 name << "NumExitCorrelations/x_PDG" << std::abs(pdgID1) << "__y_PDG" << std::abs(pdgID2) << "__highE";
1348 TH2F *hist1_2_highE = (TH2F*)m_fileLookupTable->Get(name.str().c_str());
1349 name.str("");
1350 name << "NumExitCorrelations/x_PDG" << std::abs(pdgID2) << "__y_PDG" << std::abs(pdgID1) << "__lowE";
1351 TH2F *hist2_1_lowE = (TH2F*)m_fileLookupTable->Get(name.str().c_str());
1352 name.str("");
1353 name << "NumExitCorrelations/x_PDG" << std::abs(pdgID2) << "__y_PDG" << std::abs(pdgID1) << "__highE";
1354 TH2F *hist2_1_highE = (TH2F*)m_fileLookupTable->Get(name.str().c_str());
1355 // check if the histograms exist
1356 if ( (!hist1_2_lowE) || (!hist2_1_lowE) || (!hist1_2_highE) || (!hist2_1_highE) )
1357 {
1358 ATH_MSG_ERROR("[PunchThroughG4Tool] unable to retrieve the correlation data for PDG IDs " << pdgID1 << " and " << pdgID2);
1359 return StatusCode::FAILURE;
1360 }
1361
1362 // TODO: if only one of the two histograms exists, create the other one
1363 // by mirroring the data
1364 const double lowE = getFloatAfterPatternInStr( hist1_2_lowE->GetTitle(), "elow_");
1365 const double midE = getFloatAfterPatternInStr( hist1_2_lowE->GetTitle(), "ehigh_");
1366 //TODO: check if the same:
1367 // double midE = getFloatAfterPatternInStr( hist1_2_lowE->GetTitle(), "elow_");
1368 const double upperE = getFloatAfterPatternInStr( hist1_2_highE->GetTitle(), "ehigh_");
1369 // now store the correlation either way id1->id2 and id2->id1
1370 m_particles[pdgID1]->setCorrelation(pdgID2, hist2_1_lowE, hist2_1_highE,
1371 minCorrEnergy, fullCorrEnergy,
1372 lowE, midE, upperE);
1373 m_particles[pdgID2]->setCorrelation(pdgID1, hist1_2_lowE, hist1_2_highE,
1374 minCorrEnergy, fullCorrEnergy,
1375 lowE, midE, upperE);
1376 return StatusCode::SUCCESS;
1377}
1378
1379std::unique_ptr<PunchThroughPDFCreator> PunchThroughG4Tool::readLookuptablePDF(int pdg, TFile* fileLookupTable, const std::string& folderName)
1380{
1381
1382 // will hold the PDFcreator class which will be returned at the end
1383 // this will store the distributions for the punch through particles
1384 // (as map of energy & eta of the incoming particle)
1385 //PunchThroughPDFCreator *pdf = new PunchThroughPDFCreator();
1386 std::unique_ptr<PunchThroughPDFCreator> pdf = std::make_unique<PunchThroughPDFCreator>();
1387
1388 //Get directory object
1389 std::stringstream dirName;
1390 dirName << folderName << pdg;
1391 pdf->setName(dirName.str());
1392
1393 TDirectory * dir = (TDirectory*)fileLookupTable->Get(dirName.str().c_str());
1394 if(! dir)
1395 {
1396 ATH_MSG_ERROR( "[PunchThroughG4Tool] unable to retrieve directory object ("<< folderName << pdg << ")" );
1397 return nullptr;
1398 }
1399
1400 //Get list of all objects in directory
1401 TIter keyList(dir->GetListOfKeys());
1402 TKey *key;
1403
1404 while ((key = (TKey*)keyList())) {
1405
1406 //Get histogram object from key and its name
1407 TH1* hist = nullptr;
1408
1409 std::string histName;
1410 if(strcmp(key->GetClassName(), "TH1F") == 0){
1411 hist = (TH1*)key->ReadObj();
1412 histName = hist->GetName();
1413 }
1414
1415 //extract energy and eta from hist name 6 and 1 to position delimeters correctly
1416 std::string strEnergy = histName.substr( histName.find_first_of('E') + 1, histName.find_first_of('_')-histName.find_first_of('E') - 1 );
1417 histName.erase(0, histName.find_first_of('_') + 1);
1418 std::string strEtaMin = histName.substr( histName.find("etaMin") + 6, histName.find_first_of('_') - histName.find("etaMin") - 6 );
1419 histName.erase(0, histName.find('_') + 1);
1420 std::string strEtaMax = histName.substr( histName.find("etaMax") + 6, histName.length());
1421
1422 //create integers to store in map
1423 const int energy = std::stoi(strEnergy);
1424 const int etaMin = std::stoi(strEtaMin);
1425
1426 //Add entry to pdf map
1427 pdf->addToEnergyEtaHist1DMap(energy, etaMin, hist);
1428
1429 //create doubles to store energy and eta points for interpolation
1430 const double energyDbl = static_cast<double>(energy);
1431 const double etaDbl = static_cast<double>(etaMin)/100.;
1432
1433 //create vectors to store the eta and energy points, this allows us to interpolate
1434 if (std::find(m_energyPoints.begin(), m_energyPoints.end(), energyDbl) == m_energyPoints.end()) {
1435 m_energyPoints.push_back(energyDbl);
1436 }
1437 if (std::find(m_etaPoints.begin(), m_etaPoints.end(), etaDbl) == m_etaPoints.end()) {
1438 m_etaPoints.push_back(etaDbl);
1439 }
1440
1441 }
1442
1443 return pdf;
1444}
1445
1446double PunchThroughG4Tool::getFloatAfterPatternInStr(const char *cstr, const char *cpattern)
1447{
1448 double num = 0.;
1449
1450 const std::string_view str( cstr);
1451 const std::string_view pattern( cpattern);
1452 const size_t pos = str.find(pattern);
1453
1454 if ( pos == std::string::npos)
1455 {
1456 ATH_MSG_WARNING("[PunchThroughG4Tool] unable to retrieve floating point number from string");
1457 return -999999999999.;
1458 }
1459 const std::string_view substring = str.substr(pos+pattern.length());
1460 std::from_chars(substring.data(), substring.data() + substring.size(), num);
1461 return num;
1462}
1463
1464
1465G4ThreeVector PunchThroughG4Tool::punchTroughPosPropagator(double theta, double phi, double R1, double R2, double z1, double z2) const
1466{
1467 // phi, theta - direction of the punch-through particle coming into calo
1468 // particle propagates inside the calorimeter along the straight line
1469 // coordinates of this particles when exiting the calo (on calo-MS boundary)
1470
1471 double x, y, z, r;
1472
1473 // cylinders border angles
1474 const double theta1 = atan (R1 / z1);
1475 const double theta2 = atan (R1 / z2);
1476 const double theta3 = atan (R2 / z2);
1477 //where is the particle
1478
1479 if (theta >= 0 && theta < theta1)
1480 {
1481 z = z1;
1482 r = std::fabs (z1 * tan(theta));
1483 }
1484 else if (theta >= theta1 && theta < theta2)
1485 {
1486 z = R1 / tan(theta);
1487 r = R1;
1488 }
1489 else if (theta >= theta2 && theta < theta3)
1490 {
1491 z = z2;
1492 r = std::fabs(z2 * tan(theta));;
1493 }
1494 else if (theta >= theta3 && theta < (TMath::Pi() - theta3) )
1495 {
1496 z = R2 / tan(theta);
1497 r = R2;
1498 }
1499 else if (theta >= (TMath::Pi() - theta3) && theta < (TMath::Pi() - theta2) )
1500 {
1501 z = -z2;
1502 r = std::fabs(z2 * tan(theta));
1503 }
1504 else if (theta >= (TMath::Pi() - theta2) && theta < (TMath::Pi() - theta1) )
1505 {
1506 z = R1 / tan(theta);
1507 r = R1;
1508 }
1509 else if (theta >= (TMath::Pi() - theta1) && theta <= TMath::Pi() )
1510 {
1511 z = -z1;
1512 r = std::fabs(z1 * tan(theta));
1513 }
1514
1515 //parallel universe
1516 else
1517 {
1518 ATH_MSG_WARNING("[PunchThroughG4Tool::punchTroughPosPropagator] Given theta angle is incorrect, setting particle position to (0, 0, 0)");
1519 x = 0.0; y = 0.0; z = 0.0; r = 0.0;
1520 }
1521
1522 x = r * cos(phi);
1523 y = r * sin(phi);
1524 G4ThreeVector posVec(x, y, z);
1525
1526 ATH_MSG_DEBUG("[PunchThroughG4Tool::punchTroughPosPropagator] position of produced punch-through particle: x = " << x <<" y = "<< y <<" z = "<< z<<" r = "<< posVec.perp() <<"std::sqrt(x^2+y^2) = "<< std::sqrt(x * x + y * y));
1527
1528 return posVec;
1529}
const boost::regex ref(r_ef)
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static int quant(double min, double max, unsigned nSteps, double val)
int upper(int c)
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
std::vector< T > str_to_list(const std::string_view str)
#define y
#define x
#define z
std::vector< std::map< std::string, std::string > > m_xml_info_cdf
double getFloatAfterPatternInStr(const char *str, const char *pattern)
get the floating point number in a string, after the given pattern
static double normal_cdf(double x)
std::vector< std::map< std::string, std::string > > getInfoMap(const std::string &mainNode, const std::string &xmlFilePath)
virtual StatusCode initialize() override
AlgTool initialize method.
IntegerArrayProperty m_pdgInitiators
DoubleArrayProperty m_minCorrEnergy
static std::vector< double > dotProduct(const std::vector< std::vector< double > > &m, const std::vector< double > &v)
ServiceHandle< IEnvelopeDefSvc > m_envDefSvc
DoubleArrayProperty m_numParticlesFactor
DoubleArrayProperty m_initiatorsEtaRange
StatusCode registerCorrelation(int pdgID1, int pdgID2, double minCorrEnergy=0., double fullCorrEnergy=0.)
register a correlation for the two given types of punch-through particles with a given energy thresho...
std::vector< std::map< double, double > > m_variable2_inverse_cdf
TFile * m_fileLookupTable
ROOT objects.
DoubleArrayProperty m_posAngleFactor
IntegerArrayProperty m_initiatorsMinEnergy
StatusCode initializeRegisterCorrelations()
initialize register all correlations between particles
std::map< int, PunchThroughParticle * > m_particles
needed to initially create punch-through particles with the right distributions
DoubleArrayProperty m_minEnergy
std::vector< std::vector< std::vector< double > > > m_inverse_PCA_matrix
pca vectors
virtual void createAllSecondaryTracks(G4ParticleTable &ptable, G4FastStep &fastStep, const G4Track &g4PrimaryTrack, std::vector< std::map< std::string, double > > &secKinematicsMapVect, G4TrackVector &secTrackCont, const std::vector< double > &caloMSVars) override
create all secondary tracks from kinematics map
static double inverseCdfTransform(double variable, const std::map< double, double > &inverse_cdf_map)
G4Track * createSecondaryTrack(G4ParticleTable &ptable, G4FastStep &fastStep, double currentTime, int secondarySignedPDG, double energy, double theta, double phi, double momTheta, double momPhi, const std::vector< double > &caloMSVars)
create secondary track for each given the kinematics
DoubleArrayProperty m_fullCorrEnergy
std::vector< std::map< std::string, double > > checkEnergySumFromSecondaries(double mainEnergyInit, std::vector< std::map< std::string, double > > &secKinematicsMapVect)
check the energies satisfying energy condition
std::unique_ptr< PunchThroughPDFCreator > readLookuptablePDF(int pdgID, TFile *fileLookupTable, const std::string &folderName)
reads out the lookuptable for the given type of particle
IntegerArrayProperty m_correlatedParticle
StatusCode initializeInverseCDF(const std::string &quantileTransformerConfigFile)
StatusCode initializeInversePCA(const std::string &inversePCAConfigFile)
void checkParticleTable(G4ParticleTable &ptable, int secondarySignedPDG)
int getCorrelatedParticles(const G4Track &g4PrimaryTrack, std::vector< std::map< std::string, double > > &secKinematicsMapVect, int pdg, int corrParticles, CLHEP::HepRandomEngine *rndmEngine, double interpEnergy, double interpEta)
get the right number of particles for the given pdg while considering the correlation to an other par...
static std::map< double, double > getVariableCDFmappings(xmlNodePtr &nodeParent)
int passedParamIterator(int pid, double eta, const std::vector< std::map< std::string, std::string > > &mapvect) const
std::vector< std::map< double, double > > m_variable0_inverse_cdf
(vector of map) for CDF mappings
DoubleProperty m_beamPipe
beam pipe radius
std::vector< std::vector< double > > m_PCA_means
int getAllParticles(const G4Track &g4PrimaryTrack, std::vector< std::map< std::string, double > > &secKinematicsMapVect, CLHEP::HepRandomEngine *rndmEngine, int pdg, double interpEnergy, double interpEta, int numParticles=-1)
create the right number of punch-through particles for the given pdg and return the number of particl...
virtual std::vector< double > getCaloMSVars() override
double m_R1
calo-MS borders
StatusCode initializePhysics() override
double interpolateEnergy(const double &energy, CLHEP::HepRandomEngine *rndmEngine) const
StringProperty m_filenameInversePCA
StatusCode checkCaloMSBoundaries(const std::vector< std::pair< double, double > > *rzMS, const std::vector< std::pair< double, double > > *rzCalo)
Check calo-MS boundaries.
IntegerArrayProperty m_maxNumParticles
std::vector< std::map< std::string, std::string > > m_xml_info_pca
infoMaps
ServiceHandle< ISF::IGeoIDSvc > m_geoIDSvc
BooleanArrayProperty m_doAntiParticles
IntegerArrayProperty m_punchThroughParticles
virtual StatusCode finalize() override
AlgTool finalize method.
DoubleArrayProperty m_energyFactor
double interpolateEta(const double &eta, CLHEP::HepRandomEngine *rndmEngine) const
StringProperty m_filenameInverseCDF
std::map< std::string, double > getOneParticleKinematics(CLHEP::HepRandomEngine *rndmEngine, int secondaryPDG, float initParticleTheta, float initParticlePhi, double interpEnergy, double interpEta) const
create exactly one punch-through particle with the given pdg and the given max energy
std::vector< std::map< double, double > > m_variable4_inverse_cdf
StatusCode initializeRegisterPunchThroughParticles()
initialize register all the punch-through particles which will be simulated
virtual std::vector< std::map< std::string, double > > computePunchThroughParticles(const G4FastTrack &fastTrack, CLHEP::HepRandomEngine *rndmEngine, double punchThroughProbability, double punchThroughClassifierRand) override
interface function: fill a vector with the punch-through particles
DoubleArrayProperty m_momAngleFactor
std::vector< double > m_energyPoints
energy and eta points in param
StringProperty m_filenameLookupTable
std::vector< double > inversePCA(int pcaCdfIterator, std::vector< double > &variables) const
std::vector< double > m_etaPoints
PunchThroughG4Tool(const std::string &, const std::string &, const IInterface *)
StatusCode registerPunchThroughParticle(G4ParticleTable &ptable, int pdg, bool doAntiparticle=false, double minEnergy=0., int maxNumParticles=-1, double numParticlesFactor=1., double energyFactor=1., double posAngleFactor=1., double momAngleFactor=1.)
registers a type of punch-through particles which will be simulated
std::vector< std::map< double, double > > m_variable3_inverse_cdf
G4ThreeVector punchTroughPosPropagator(double theta, double phi, double R1, double R2, double z1, double z2) const
get particle through the calorimeter
std::vector< std::map< double, double > > m_variable1_inverse_cdf
This class holds information for different properties of a punch-through particle (energy,...
Definition node.h:24
void name(const std::string &n)
Definition node.h:40
int r
Definition globals.cxx:22
Eigen::Matrix< double, 3, 1 > Vector3D
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
void AddXmlToCollectionMap(xmlNodePtr node, const char *attrName, COL &collection)
void AddXmlToCollection(xmlNodePtr node, const char *attrName, COL &collection)
T deltaPhi(T phiA, T phiB)
Return difference phiA - phiB in range [-pi, pi].
Definition phihelper.h:42
T GetXmlAttr(xmlNodePtr node, const char *attrName, const T &defaultValue=T{}) noexcept(std::is_arithmetic_v< T >)
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void swap(ElementLinkVector< DOBJ > &lhs, ElementLinkVector< DOBJ > &rhs)