8 #include "G4ParticleTable.hh" 9 #include "G4DecayTable.hh" 10 #include "G4PhaseSpaceDecayChannel.hh" 11 #include "G4ProcessManager.hh" 12 #include "G4ParticleDefinition.hh" 13 #include "G4SystemOfUnits.hh" 22 #ifdef G4MULTITHREADED 23 G4Mutex CustomParticleFactory::customParticleFactoryMutex = G4MUTEX_INITIALIZER;
34 #ifdef G4MULTITHREADED 35 G4MUTEXLOCK(&customParticleFactoryMutex);
43 std::ifstream configFile(filePath.c_str());
46 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Reading Custom Particle and G4DecayTable from \n" 49 G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
50 while (getline(configFile, line)) {
51 line.erase(0, line.find_first_not_of(
" \t"));
52 if (line.length() == 0 || line.at(0) ==
'#') {
56 if (
ToLower(line).find(
"block") < line.npos &&
ToLower(line).find(
"mass") < line.npos) {
57 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Retrieving mass table.";
60 if (line.find(
"DECAY") < line.npos) {
64 std::stringstream lineStream(line);
65 lineStream >> tmpString >> pdgId >>
width;
68 <<
"CustomParticleFactory: entry to G4DecayTable: pdgID, width " << pdgId <<
", " <<
width;
69 G4ParticleDefinition *aParticle = theParticleTable->FindParticle(pdgId);
70 if (!aParticle || width == 0.0) {
73 G4DecayTable *aDecayTable =
getDecayTable(&configFile, pdgId);
74 aParticle->SetDecayTable(aDecayTable);
75 aParticle->SetPDGStable(
false);
77 G4ParticleDefinition *aAntiParticle = theParticleTable->FindAntiParticle(pdgId);
78 if (aAntiParticle && aAntiParticle->GetPDGEncoding() !=
pdgId) {
80 aAntiParticle->SetPDGStable(
false);
81 aAntiParticle->SetPDGLifeTime(1.0 / (width * CLHEP::GeV) * 6.582122
e-22 *
CLHEP::MeV * CLHEP::s);
85 #ifdef G4MULTITHREADED 86 G4MUTEXUNLOCK(&customParticleFactoryMutex);
93 <<
"Pdg code too low " << pdgCode <<
" " <<
std::abs(pdgCode) / 1000000;
98 G4String pType =
"custom";
99 G4String pSubType =
"";
100 G4double spectatormass = 0.0;
101 G4ParticleDefinition *spectator =
nullptr;
119 if (name.compare(0, 4,
"~HIP") == 0) {
120 if ((name.compare(0, 7,
"~HIPbar") == 0)) {
122 charge = CLHEP::eplus * atoi(str.c_str()) / 3.;
125 charge = -CLHEP::eplus * atoi(str.c_str()) / 3.;
128 if (name.compare(0, 9,
"anti_~HIP") == 0) {
129 if ((name.compare(0, 12,
"anti_~HIPbar") == 0)) {
131 charge = -CLHEP::eplus * atoi(str.c_str()) / 3.;
134 charge = CLHEP::eplus * atoi(str.c_str()) / 3.;
146 double lifetime = -1;
162 G4DecayTable *decaytable =
nullptr;
163 G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
183 if (pType ==
"rhadron" && name !=
"~g") {
184 G4String cloudname = name +
"cloud";
185 G4String cloudtype = pType +
"cloud";
186 spectator = theParticleTable->FindParticle(1000021);
187 spectatormass = spectator->GetPDGMass();
188 G4double cloudmass = massGeV - spectatormass;
190 new CustomParticle(cloudname, cloudmass, 0.0, 0, 0, +1, 0, 0, 0, 0, cloudtype, 0, +1, 0,
true, -1.0,
nullptr);
195 <<
"CustomParticleFactory: " << name <<
" being assigned spectator" << spectator->GetParticleName()
196 <<
" and cloud " << cloudname <<
"\n" 197 <<
" Masses: " << mass <<
" Gev, " << spectatormass / CLHEP::GeV <<
" GeV and " 198 << cloudmass / CLHEP::GeV <<
" GeV.";
199 }
else if (pType ==
"mesonino" || pType ==
"sbaryon") {
204 G4String cloudname = name +
"cloud";
205 G4String cloudtype = pType +
"cloud";
207 spectator = theParticleTable->FindParticle(1000006 * sign);
210 spectator = theParticleTable->FindParticle(1000005 * sign);
213 edm::LogError(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Cannot find spectator parton";
217 spectatormass = spectator->GetPDGMass();
219 G4double cloudmass = massGeV - spectatormass;
221 new CustomParticle(cloudname, cloudmass, 0.0, 0, 0, +1, 0, 0, 0, 0, cloudtype, 0, +1, 0,
true, -1.0,
nullptr);
227 <<
"CustomParticleFactory: " << name <<
" being assigned spectator" << spectator->GetParticleName()
228 <<
" and cloud " << cloudname <<
"\n" 229 <<
" Masses: " << mass <<
" Gev, " << spectatormass / CLHEP::GeV <<
" GeV and " 230 << cloudmass / CLHEP::GeV <<
" GeV.";
243 G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
246 while (getline(*configFile, line)) {
247 line.erase(0, line.find_first_not_of(
" \t"));
248 if (line.length() == 0 || line.at(0) ==
'#')
250 if (
ToLower(line).find(
"block") < line.npos) {
251 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Finished the Mass Table ";
254 std::stringstream sstr(line);
255 sstr >> pdgId >> mass >> tmp >>
name;
258 if (theParticleTable->FindParticle(pdgId)) {
263 <<
"CustomParticleFactory: Calling addCustomParticle for pdgId: " << pdgId <<
", mass " << mass <<
" GeV " 269 int pdgIdPartner = pdgId % 100;
270 G4ParticleDefinition *aParticle = theParticleTable->FindParticle(pdgIdPartner);
274 <<
"CustomParticleFactory: Found partner particle for " 275 <<
" pdgId= " << pdgId <<
", pdgIdPartner= " << pdgIdPartner <<
" " << aParticle->GetParticleName();
279 pdgId != 1000006 && pdgId != -1000006 && pdgId != 25 && pdgId != 35 && pdgId != 36 && pdgId != 37) {
280 int sign = aParticle->GetAntiPDGEncoding() / pdgIdPartner;
282 <<
"CustomParticleFactory: For " << aParticle->GetParticleName() <<
" pdg= " << pdgIdPartner
283 <<
", GetAntiPDGEncoding() " << aParticle->GetAntiPDGEncoding() <<
" sign= " <<
sign;
286 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: sign= " << sign
287 <<
" a: " << aParticle->GetAntiPDGEncoding() <<
" b: " << pdgIdPartner;
288 aParticle->DumpTable();
290 if (sign == -1 && pdgId != 25 && pdgId != 35 && pdgId != 36 && pdgId != 37 && pdgId != 1000039) {
291 tmp =
"anti_" +
name;
292 if (!theParticleTable->FindParticle(-pdgId)) {
294 <<
"CustomParticleFactory: Calling addCustomParticle for antiparticle with pdgId: " << -pdgId <<
", mass " 295 << mass <<
" GeV, name " <<
tmp;
297 theParticleTable->FindParticle(pdgId)->SetAntiPDGEncoding(-pdgId);
300 theParticleTable->FindParticle(pdgId)->SetAntiPDGEncoding(pdgId);
303 if (pdgId == 1000039)
304 theParticleTable->FindParticle(pdgId)->SetAntiPDGEncoding(pdgId);
305 if (pdgId == 1000024 || pdgId == 1000037 || pdgId == 37) {
306 tmp =
"anti_" +
name;
307 if (!theParticleTable->FindParticle(-pdgId)) {
309 <<
"CustomParticleFactory: Calling addCustomParticle for antiparticle (2) with pdgId: " << -pdgId
310 <<
", mass " << mass <<
" GeV, name " <<
tmp;
312 theParticleTable->FindParticle(pdgId)->SetAntiPDGEncoding(-pdgId);
324 G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
326 std::string parentName = theParticleTable->FindParticle(pdgId)->GetParticleName();
327 G4DecayTable *decaytable =
new G4DecayTable();
329 while (getline(*configFile, line)) {
330 line.erase(0, line.find_first_not_of(
" \t"));
331 if (line.length() == 0)
333 if (line.at(0) ==
'#' &&
ToLower(line).find(
"br") < line.npos &&
ToLower(line).find(
"nda") < line.npos)
335 if (line.at(0) ==
'#' ||
ToLower(line).find(
"block") < line.npos) {
336 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Finished the Decay Table ";
340 std::stringstream sstr(line);
341 sstr >> br >> nDaughters;
342 LogDebug(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Branching Ratio: " << br
343 <<
", Number of Daughters: " << nDaughters;
344 if (nDaughters > 4) {
346 <<
"CustomParticleFactory: Number of daughters is too large (max = 4): " << nDaughters
347 <<
" for pdgId: " <<
pdgId;
352 <<
"CustomParticleFactory: Branching ratio is " << br <<
" for pdgId: " <<
pdgId;
356 for (
int i = 0;
i < nDaughters; ++
i) {
358 LogDebug(
"SimG4CoreCustomPhysics") <<
"CustomParticleFactory: Daughter ID " << pdg[
i];
359 const G4ParticleDefinition *
part = theParticleTable->FindParticle(pdg[
i]);
362 <<
"CustomParticleFactory: particle with PDG code" << pdg[
i] <<
" not found!";
365 name[
i] = part->GetParticleName();
368 G4PhaseSpaceDecayChannel *aDecayChannel =
369 new G4PhaseSpaceDecayChannel(parentName, br, nDaughters, name[0], name[1], name[2], name[3]);
370 decaytable->Insert(aDecayChannel);
372 <<
"CustomParticleFactory: inserted decay channel " 373 <<
" for pdgID= " << pdgId <<
" " << parentName <<
" BR= " << br <<
" Daughters: " << name[0] <<
" " 374 << name[1] <<
" " << name[2] <<
" " << name[3];
381 G4ParticleTable *theParticleTable = G4ParticleTable::GetParticleTable();
383 std::string parentName = theParticleTable->FindParticle(-pdgId)->GetParticleName();
384 G4DecayTable *decaytable =
new G4DecayTable();
386 for (
int i = 0;
i < theDecayTable->entries(); ++
i) {
387 G4VDecayChannel *theDecayChannel = theDecayTable->GetDecayChannel(
i);
388 int nd =
std::min(4, theDecayChannel->GetNumberOfDaughters());
389 for (
int j = 0; j < nd; ++j) {
390 int id = theDecayChannel->GetDaughter(j)->GetAntiPDGEncoding();
391 const G4ParticleDefinition *
part = theParticleTable->FindParticle(
id);
394 <<
"CustomParticleFactory: antiparticle with PDG code" <<
id <<
" not found!";
397 name[j] = part->GetParticleName();
399 G4PhaseSpaceDecayChannel *aDecayChannel =
400 new G4PhaseSpaceDecayChannel(parentName, theDecayChannel->GetBR(), nd, name[0], name[1], name[2], name[3]);
401 decaytable->Insert(aDecayChannel);
403 <<
"CustomParticleFactory: inserted decay channel " 404 <<
" for pdgID= " << -pdgId <<
" " << parentName <<
" BR= " << theDecayChannel->GetBR()
405 <<
" Daughters: " << name[0] <<
" " << name[1] <<
" " << name[2] <<
" " << name[3];
413 str.at(
i) = std::tolower(str.at(
i),
loc);
static bool s_isDphoton(int pdg)
static double s_charge(int pdg)
static bool s_isSbaryon(int pdg)
void SetCloud(G4ParticleDefinition *theCloud)
void loadCustomParticles(const std::string &filePath)
static std::vector< G4ParticleDefinition * > m_particles
static bool s_isMesonino(int pdg)
static bool s_isSLepton(int pdg)
Abs< T >::type abs(const T &t)
std::string ToLower(std::string str)
void SetSpectator(G4ParticleDefinition *theSpectator)
static double s_spin(int pdg)
void getMassTable(std::ifstream *configFile)
G4DecayTable * getDecayTable(std::ifstream *configFile, int pdgId)
void addCustomParticle(int pdgCode, double mass, const std::string &name)
G4DecayTable * getAntiDecayTable(int pdgId, G4DecayTable *theDecayTable)
static bool s_isRHadron(int pdg)
std::vector< std::vector< double > > tmp
const std::vector< G4ParticleDefinition * > & GetCustomParticles()
static bool s_issbottomHadron(int pdg)
static bool s_isstopHadron(int pdg)