ATLAS Offline Software
G4ProcessHelper.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #include "G4ProcessHelper.hh"
6 #include "CustomParticle.h"
8 #include "G4ParticleTable.hh"
9 #include "G4DecayTable.hh"
10 #include "CLHEP/Random/RandFlat.h"
11 #include "CLHEP/Units/PhysicalConstants.h"
12 #include <iostream>
13 #include <fstream>
14 #include <stdexcept>
15 
17 
18 G4ProcessHelper::G4ProcessHelper()
19  : theRmesoncloud(0)
20  , theRbaryoncloud(0)
21 {
22  G4cout << "G4ProcessHelper constructor: start" << G4endl;
23  particleTable = G4ParticleTable::GetParticleTable();
24  theProton = particleTable->FindParticle("proton");
25  theNeutron = particleTable->FindParticle("neutron");
26 
27 
28  // I opted for string based read-in, as it would make physics debugging easier later on
29 
30  std::ifstream process_stream ("ProcessList.txt");
31  G4String line;
32  while(getline(process_stream,line)){
33  std::vector<G4String> tokens;
34 
35  //Getting a line
36  ReadAndParse(line,tokens,"#");
37 
38  //Important info
39  G4String incident = tokens[0];
40  // G4cout << "Incident particle: " << incident << G4endl;
41  G4ParticleDefinition* incidentDef = particleTable->FindParticle(incident);
42  G4int incidentPDG = incidentDef->GetPDGEncoding();
43  known_particles[incidentDef]=true;
44 
45  G4String target = tokens[1];
46  // G4cout << "Target particle: " << target << G4endl;
47 
48  // Making a ReactionProduct
49  ReactionProduct prod;
50  for (unsigned int i = 2; i != tokens.size();i++){
51  G4String part = tokens[i];
52  if (particleTable->contains(part))
53  {
54  prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
55  } else {
56  G4cout<<"Particle: "<<part<<" is unknown."<<G4endl;
57  G4Exception("G4ProcessHelper", "UnkownParticle", FatalException,
58  "Initialization: The reaction product list contained an unknown particle");
59  }
60  }
61  if (target == "proton"){
62  pReactionMap[incidentPDG].push_back(prod);
63  } else if (target == "neutron") {
64  nReactionMap[incidentPDG].push_back(prod);
65  } else {
66  G4Exception("G4ProcessHelper", "IllegalTarget", FatalException,
67  "Initialization: The reaction product list contained an illegal target particle");
68  }
69  }
70 
71  process_stream.close();
72  G4cout << "Found " << pReactionMap.size() << " proton interactions and " << nReactionMap.size() << " neutron interactions in ProcessList.txt." << G4endl;
73 
74  std::map<G4String,G4double> parameters;
75  ReadInPhysicsParameters(parameters);
76 
77  resonant = false;
78  reggemodel = false;
79  if (parameters["Resonant"]!=0.) resonant=true;
80  ek_0 = parameters["ResonanceEnergy"]*CLHEP::GeV;
81  gamma = parameters["Gamma"]*CLHEP::GeV;
82  amplitude = parameters["Amplitude"]*CLHEP::millibarn;
83  xsecmultiplier = parameters["XsecMultiplier"];
84  suppressionfactor = parameters["ReggeSuppression"];
85  hadronlifetime = parameters["HadronLifeTime"];
86  mixing = parameters["Mixing"];
87  if(parameters["ReggeModel"]!=0.) reggemodel=true;
88  doDecays=parameters["DoDecays"];
89 
90  G4cout<<"Read in physics parameters:"<<G4endl;
91  G4cout<<"Resonant = "<< resonant <<G4endl;
92  G4cout<<"ResonanceEnergy = "<<ek_0/CLHEP::GeV<<" GeV"<<G4endl;
93  G4cout<<"XsecMultiplier = "<<xsecmultiplier<<G4endl;
94  G4cout<<"Gamma = "<<gamma/CLHEP::GeV<<" GeV"<<G4endl;
95  G4cout<<"Amplitude = "<<amplitude/CLHEP::millibarn<<" millibarn"<<G4endl;
96  G4cout<<"ReggeSuppression = "<<100*suppressionfactor<<" %"<<G4endl;
97  G4cout<<"HadronLifeTime = "<<hadronlifetime;
98  if (doDecays) G4cout<<" ns"<<G4endl;
99  else G4cout<<" s"<<G4endl;
100  G4cout<<"ReggeModel = "<< reggemodel <<G4endl;
101  G4cout<<"Mixing = "<< mixing*100 <<" %"<<G4endl;
102  G4cout<<"DoDecays = "<< doDecays << G4endl;
103 
104  if ((!doDecays && hadronlifetime>0.) ||
105  (doDecays && hadronlifetime<=0.) ){
106  G4cout << "WARNING: Inconsistent treatment of R-Hadron properties! Lifetime of " << hadronlifetime
107  << " and doDecays= " << doDecays << G4endl;
108  }
109 
110  G4ParticleTable::G4PTblDicIterator* theParticleIterator;
111  theParticleIterator = particleTable->GetIterator();
112 
113  theParticleIterator->reset();
114  while( (*theParticleIterator)() ){
115  CustomParticle* particle = dynamic_cast<CustomParticle*>(theParticleIterator->value());
116  std::string name = theParticleIterator->value()->GetParticleName();
117  G4DecayTable* table = theParticleIterator->value()->GetDecayTable();
118  if(particle!=0&&table!=0&&name.find("cloud")>name.size()) {
119  particle->SetPDGLifeTime(hadronlifetime*CLHEP::s);
120  particle->SetPDGStable(false);
121  G4cout<<"Lifetime of: "<<name<<" set to: "<<particle->GetPDGLifeTime()/CLHEP::s<<" s."<<G4endl;
122  G4cout<<"Stable: "<<particle->GetPDGStable()<<G4endl;
123  }
124  if (particle && doDecays==1){
125  // Make them decay immediately!!
126  particle->SetPDGStable(false);
127  particle->SetPDGLifeTime(hadronlifetime*CLHEP::ns);
128  G4cout<<"Forcing a decay for "<<name<<G4endl;
129  G4cout<<"Lifetime of: "<<name<<" set to: "<<particle->GetPDGLifeTime()/CLHEP::ns<<" ns."<<G4endl;
130  G4cout<<"Stable: "<<particle->GetPDGStable()<<G4endl;
131  }
132 
133  G4cout << "Done with particle " << name << G4endl;
134  }
135  theParticleIterator->reset();
136  G4cout << "G4ProcessHelper constructor: end" << G4endl;
137  return;
138 }
139 
140 
141 void G4ProcessHelper::ReadInPhysicsParameters(std::map<G4String,G4double>& parameters) const
142  {
143  parameters["Resonant"]=0.;
144  parameters["ResonanceEnergy"]=0.;
145  parameters["XsecMultiplier"]=1.;
146  parameters["Gamma"]=0.;
147  parameters["Amplitude"]=0.;
148  parameters["ReggeSuppression"]=0.;
149  parameters["HadronLifeTime"]=0.;
150  parameters["ReggeModel"]=0.;
151  parameters["Mixing"]=0.;
152  parameters["DoDecays"]=0;
153 
154  std::ifstream physics_stream ("PhysicsConfiguration.txt");
155  G4String line;
156  char** endptr=0;
157  while (getline(physics_stream,line)) {
158  std::vector<G4String> tokens;
159  //Getting a line
160  ReadAndParse(line,tokens,"=");
161  G4String key = tokens[0];
162  G4double val = strtod(tokens[1],endptr);
163  parameters[key]=val;
164  }
165  physics_stream.close();
166  }
167 
168 
169 const G4ProcessHelper* G4ProcessHelper::Instance()
170 {
171  static const G4ProcessHelper instance;
172  return &instance;
173 }
174 
175 
176 G4bool G4ProcessHelper::ApplicabilityTester(const G4ParticleDefinition& aPart) const {
177  try {
178  return known_particles.at(&aPart);
179  }
180  catch (const std::out_of_range& e) {
181  return false;
182  }
183 }
184 
185 G4double G4ProcessHelper::GetInclusiveCrossSection(const G4DynamicParticle *aParticle,
186  const G4Element *anElement) const {
187  //We really do need a dedicated class to handle the cross sections. They might not always be constant
188 
189  //Disassemble the PDG-code
190  G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
191  double boost = (aParticle->GetKineticEnergy()+aParticle->GetMass())/aParticle->GetMass();
192  G4double theXsec = 0;
193  G4String name = aParticle->GetDefinition()->GetParticleName();
194 
195  if(!reggemodel){
196  //Flat cross section
197  if(MC::SUSY::isRGlueball(thePDGCode)) {
198  theXsec = 24 * CLHEP::millibarn;
199  } else {
200  std::vector<G4int> nq=MC::SUSY::containedQuarks(thePDGCode);
201  for (std::vector<G4int>::iterator it = nq.begin();
202  it != nq.end();
203  ++it)
204  {
205  // 12 mb taken from asymptotic pion-nucleon scattering cross sections
206  if (*it == 1 || *it == 2) theXsec += 12 * CLHEP::millibarn;
207  // 6 mb taken from asymptotic kaon-nucleon scattering cross sections
208  // No data for D or B, so setting to behave like a kaon
209  if (*it == 3 || *it == 4 || *it == 5) theXsec += 6 * CLHEP::millibarn;
210  }
211  }
212  } else {
213  // From Eur. Phys. J. C (2010) 66: 493-501
214  // DOI 10.1140/epjc/s10052-010-1262-1
215  double R = Regge(boost);
216  double P = Pom(boost);
217  if(thePDGCode>0)
218  {
219  if(MC::SUSY::isSMeson(thePDGCode)) theXsec=(P+R)*CLHEP::millibarn;
220  if(MC::SUSY::isSBaryon(thePDGCode)) theXsec=2*P*CLHEP::millibarn;
221  if(MC::SUSY::isRMeson(thePDGCode)||MC::SUSY::isRGlueball(thePDGCode)) theXsec=(R+2*P)*CLHEP::millibarn;
222  if(MC::SUSY::isRBaryon(thePDGCode)) theXsec=3*P*CLHEP::millibarn;
223  }
224  else
225  {
226  if(MC::SUSY::isSMeson(thePDGCode)) theXsec=P*CLHEP::millibarn;
227  if(MC::SUSY::isSBaryon(thePDGCode)) theXsec=(2*(P+R)+30/sqrt(boost))*CLHEP::millibarn;
228  if(MC::SUSY::isRMeson(thePDGCode)||MC::SUSY::isRGlueball(thePDGCode)) theXsec=(R+2*P)*CLHEP::millibarn;
229  if(MC::SUSY::isRBaryon(thePDGCode)) theXsec=3*P*CLHEP::millibarn;
230  }
231  }
232 
233 
234 
235  //Adding resonance
236 
237  if(resonant)
238  {
239  // Described in Section 5.1 of http://r-hadrons.web.cern.ch/r-hadrons/download/mackeprang_thesis.pdf
240  // mentioned but dismissed in Section 3.3 of https://arxiv.org/pdf/hep-ex/0404001.pdf
241  double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass(); //Now total energy
242 
243  e_0 = sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
244  + theProton->GetPDGMass()*theProton->GetPDGMass()
245  + 2.*e_0*theProton->GetPDGMass());
246  double sqrts=sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
247  + theProton->GetPDGMass()*theProton->GetPDGMass() + 2*aParticle->GetTotalEnergy()*theProton->GetPDGMass());
248 
249  double res_result = amplitude*(gamma*gamma/4.)/((sqrts-e_0)*(sqrts-e_0)+(gamma*gamma/4.));//Non-relativistic Breit Wigner
250 
251  theXsec += res_result;
252  }
253 
254 
255  return theXsec * pow(anElement->GetN(),0.7)*1.25 * xsecmultiplier;// * 0.523598775598299;
256 
257 }
258 
259 ReactionProduct G4ProcessHelper::GetFinalState(const G4Track& aTrack, G4ParticleDefinition*& aTarget) const {
260  return GetFinalStateInternal(aTrack,aTarget,false);
261 }
262 
263 // Version where we know if we baryonize already
264 ReactionProduct G4ProcessHelper::GetFinalStateInternal(const G4Track& aTrack,G4ParticleDefinition*& aTarget, const bool baryonize_failed) const {
265 
266  const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
267 
268  //-----------------------------------------------
269  // Choose n / p as target
270  // and get ReactionProductList pointer
271  //-----------------------------------------------
272 
273  G4Material* aMaterial = aTrack.GetMaterial();
274  const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
275  const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
276 
277  G4double NumberOfProtons=0;
278  G4double NumberOfNucleons=0;
279 
280  for ( size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
281  {
282  //Summing number of protons per unit volume
283  NumberOfProtons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetZ();
284  //Summing nucleons (not neutrons)
285  NumberOfNucleons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetN();
286  }
287 
288  const ReactionMap* reactionMap;
289  if(CLHEP::RandFlat::shoot()<NumberOfProtons/NumberOfNucleons){
290  reactionMap = &pReactionMap;
291  aTarget = theProton;
292  } else {
293  reactionMap = &nReactionMap;
294  aTarget = theNeutron;
295  }
296 
297  G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
298 
299  if(reggemodel
300  &&MC::SUSY::isSMeson(theIncidentPDG)
301  &&CLHEP::RandFlat::shoot()*mixing>0.5
302  &&aDynamicParticle->GetDefinition()->GetPDGCharge()==0.
303  )
304  {
305  // G4cout<<"Oscillating..."<<G4endl;
306  theIncidentPDG *= -1;
307  }
308 
309 
310  bool baryonise=false;
311 
312  if(!baryonize_failed
313  && reggemodel
314  && CLHEP::RandFlat::shoot()>0.9
315  && (
316  (MC::SUSY::isSMeson(theIncidentPDG)&&theIncidentPDG>0)
317  ||
318  MC::SUSY::isRMeson(theIncidentPDG)
319  )
320  ){
321  baryonise=true;
322  }
323 
324  // Reference directly to the ReactionProductList we are looking at. Makes life easier :-)
325  const ReactionProductList& aReactionProductList = reactionMap->at(theIncidentPDG);
326 
327  //-----------------------------------------------
328  // Count processes
329  // kinematic check
330  // compute number of 2 -> 2 and 2 -> 3 processes
331  //-----------------------------------------------
332 
333  G4int N22 = 0; //Number of 2 -> 2 processes
334  G4int N23 = 0; //Number of 2 -> 3 processes. Couldn't think of more informative names
335 
336  //This is the list to be populated
337  ReactionProductList theReactionProductList;
338  std::vector<bool> theChargeChangeList;
339 
340  for (const ReactionProduct& prod : aReactionProductList) {
341  G4int secondaries = prod.size();
342  // If the reaction is not possible we will not consider it
343  /* if(ReactionIsPossible(*prod_it,aDynamicParticle)
344  &&(
345  !baryonise||(baryonise&&ReactionGivesBaryon(*prod_it))
346  ))*/
347  if(ReactionIsPossible(prod,*aTarget,aDynamicParticle)
348  &&(
349  (baryonise&&ReactionGivesBaryon(prod))
350  ||
351  (!baryonise&&!ReactionGivesBaryon(prod))
352  ||
353  (MC::SUSY::isSBaryon(theIncidentPDG))
354  ||
355  (MC::SUSY::isRBaryon(theIncidentPDG))
356  ||!reggemodel
357  )
358  )
359  {
360  // The reaction is possible. Let's store and count it
361  theReactionProductList.push_back(prod);
362  if (secondaries == 2){
363  N22++;
364  } else if (secondaries ==3) {
365  N23++;
366  } else {
367  G4cerr << "ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
368  }
369  }
370  }
371 
372  if (theReactionProductList.size()==0 && baryonize_failed){
373  G4Exception("G4ProcessHelper", "NoProcessPossible", FatalException,
374  "GetFinalState: No process could be selected from the given list.");
375  } else if (theReactionProductList.size()==0 && !baryonize_failed) {
376  // Baryonization had not yet failed -- try again
377  G4cout << "G4ProcessHelper::GetFinalStateInternal WARNING Could not select an appropriate process in first pass" << G4endl;
378  return GetFinalStateInternal(aTrack,aTarget,true);
379  }
380 
381  // For the Regge model no phase space considerations. We pick a process at random
382  if(reggemodel)
383  {
384  int n_rps = theReactionProductList.size();
385  int select = (int)(CLHEP::RandFlat::shoot()*n_rps);
386  // G4cout<<"Possible: "<<n_rps<<", chosen: "<<select<<G4endl;
387  return theReactionProductList[select];
388  }
389 
390  // Fill a probability map. Remember total probability
391  // 2->2 is 0.15*1/n_22 2->3 uses phase space
392  G4double p22 = 0.15;
393  G4double p23 = 1-p22; // :-)
394 
395  std::vector<G4double> Probabilities;
396  std::vector<G4bool> TwotoThreeFlag;
397 
398  G4double CumulatedProbability = 0;
399 
400  // To each ReactionProduct we assign a cumulated probability and a flag
401  // discerning between 2 -> 2 and 2 -> 3
402  for (unsigned int i = 0; i != theReactionProductList.size(); i++){
403  if (theReactionProductList[i].size() == 2) {
404  if(0==N22) { throw std::runtime_error("G4ProcessHelper::GetFinalState: N22 is zero!"); }
405  CumulatedProbability += p22/N22;
406  TwotoThreeFlag.push_back(false);
407  } else {
408  if(0==N23) { throw std::runtime_error("G4ProcessHelper::GetFinalState: N23 is zero!"); }
409  CumulatedProbability += p23/N23;
410  TwotoThreeFlag.push_back(true);
411  }
412  Probabilities.push_back(CumulatedProbability);
413  }
414 
415  //Renormalising probabilities
416  for (std::vector<G4double>::iterator it = Probabilities.begin();
417  it != Probabilities.end();
418  ++it)
419  {
420  *it /= CumulatedProbability;
421  }
422 
423  // Choosing ReactionProduct
424 
425  G4bool selected = false;
426  G4int tries = 0;
427  // ReactionProductList::iterator prod_it;
428 
429  //Keep looping over the list until we have a choice, or until we have tried 100 times
430  unsigned int i=0;
431  while(!selected && tries < 100){
432  i=0;
433  G4double dice = CLHEP::RandFlat::shoot();
434  while( i<theReactionProductList.size() && dice>Probabilities[i] ){
435  i++;
436  }
437 
438  if(!TwotoThreeFlag[i]) {
439  // 2 -> 2 processes are chosen immediately
440  selected = true;
441  } else {
442  // 2 -> 3 processes require a phase space lookup
443  if (PhaseSpace(theReactionProductList[i],*aTarget,aDynamicParticle)>CLHEP::RandFlat::shoot()) selected = true;
444  }
445  // double suppressionfactor=0.5;
446  auto table ATLAS_THREAD_SAFE = particleTable; // safe because table has been loaded by now
447  if(selected && table->FindParticle(theReactionProductList[i][0])->GetPDGCharge()!=aDynamicParticle->GetDefinition()->GetPDGCharge())
448  {
449  /*
450  G4cout<<"Incoming particle "<<aDynamicParticle->GetDefinition()->GetParticleName()
451  <<" has charge "<<aDynamicParticle->GetDefinition()->GetPDGCharge()<<G4endl;
452  G4cout<<"Suggested particle "<<particleTable->FindParticle(theReactionProductList[i][0])->GetParticleName()
453  <<" has charge "<<particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()<<G4endl;
454  */
455  if(CLHEP::RandFlat::shoot()<suppressionfactor) selected = false;
456  }
457  tries++;
458  // G4cout<<"Tries: "<<tries<<G4endl;
459  }
460  if(tries>=100) G4cerr<<"Could not select process!!!!"<<G4endl;
461 
462  //Return the chosen ReactionProduct
463  return theReactionProductList[i];
464 }
465 
466 G4double G4ProcessHelper::ReactionProductMass(const ReactionProduct& aReaction,const G4ParticleDefinition& aTarget,const G4DynamicParticle* aDynamicParticle) const{
467  // Incident energy:
468  G4double E_incident = aDynamicParticle->GetTotalEnergy();
469  //G4cout<<"Total energy: "<<E_incident<<" Kinetic: "<<aDynamicParticle->GetKineticEnergy()<<G4endl;
470  // sqrt(s)= sqrt(m_1^2 + m_2^2 + 2 E_1 m_2)
471  G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
472  G4double m_2 = aTarget.GetPDGMass();
473  //G4cout<<"M_R: "<<m_1/CLHEP::GeV<<" GeV, M_np: "<<m_2/CLHEP::GeV<<" GeV"<<G4endl;
474  G4double sqrts = sqrt(m_1*m_1 + m_2*(m_2 + 2 * E_incident));
475  //G4cout<<"sqrt(s) = "<<sqrts/CLHEP::GeV<<" GeV"<<G4endl;
476  // Sum of rest masses after reaction:
477  G4double M_after = 0;
478  for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it !=aReaction.end(); ++r_it){
479  //G4cout<<"Mass contrib: "<<(particleTable->FindParticle(*r_it)->GetPDGMass())/CLHEP::MeV<<" MeV"<<G4endl;
480  auto table ATLAS_THREAD_SAFE = particleTable; // safe because table has been loaded by now
481  M_after += table->FindParticle(*r_it)->GetPDGMass();
482  }
483  //G4cout<<"Intending to return this ReactionProductMass: " << sqrts << " - " << M_after << " MeV"<<G4endl;
484  return sqrts - M_after;
485 }
486 
487 G4bool G4ProcessHelper::ReactionIsPossible(const ReactionProduct& aReaction,const G4ParticleDefinition& aTarget,const G4DynamicParticle* aDynamicParticle) const{
488  if (ReactionProductMass(aReaction,aTarget,aDynamicParticle)>0) return true;
489  return false;
490 }
491 
492 G4bool G4ProcessHelper::ReactionGivesBaryon(const ReactionProduct& aReaction) const{
493  for (ReactionProduct::const_iterator it = aReaction.begin();it!=aReaction.end();++it)
494  if(MC::SUSY::isSBaryon(*it)||MC::SUSY::isRBaryon(*it)) return true;
495  return false;
496 }
497 
498 G4double G4ProcessHelper::PhaseSpace(const ReactionProduct& aReaction,const G4ParticleDefinition& aTarget,const G4DynamicParticle* aDynamicParticle) const {
499  G4double qValue = ReactionProductMass(aReaction,aTarget,aDynamicParticle);
500  // Eq 4 of https://arxiv.org/pdf/hep-ex/0404001.pdf
501  G4double phi = sqrt(1+qValue/(2*0.139*CLHEP::GeV))*pow(qValue/(1.1*CLHEP::GeV),3./2.);
502  return (phi/(1+phi));
503 }
504 
505 void G4ProcessHelper::ReadAndParse(const G4String& str,
506  std::vector<G4String>& tokens,
507  const G4String& delimiters) const
508 {
509  // Skip delimiters at beginning.
510  G4String::size_type lastPos = str.find_first_not_of(delimiters, 0);
511  if (lastPos==G4String::npos) return;
512 
513  // Find first "non-delimiter".
514  G4String::size_type pos = str.find_first_of(delimiters, lastPos);
515 
516  while (G4String::npos != pos || G4String::npos != lastPos)
517  {
518  //Skipping leading / trailing whitespaces
519  G4String temp = str.substr(lastPos, pos - lastPos);
520  while(temp.c_str()[0] == ' ') temp.erase(0,1);
521  while(temp[temp.size()-1] == ' ') temp.erase(temp.size()-1,1);
522 
523  // Found a token, add it to the vector.
524  tokens.push_back(temp);
525 
526  // Skip delimiters. Note the "not_of"
527  lastPos = str.find_first_not_of(delimiters, pos);
528  if (lastPos==G4String::npos) continue;
529 
530  // Find next "non-delimiter"
531  pos = str.find_first_of(delimiters, lastPos);
532  }
533 }
534 
535 double G4ProcessHelper::Regge(const double boost) const
536 {
537  // https://link.springer.com/content/pdf/10.1140%2Fepjc%2Fs10052-010-1262-1.pdf Eq 1
538  // Originally from https://arxiv.org/pdf/0710.3930.pdf
539  double a=2.165635078566177;
540  double b=0.1467453738547229;
541  double c=-0.9607903711871166;
542  return 1.5*exp(a+b/boost+c*log(boost));
543 }
544 
545 
546 double G4ProcessHelper::Pom(const double boost) const
547 {
548  // https://link.springer.com/content/pdf/10.1140%2Fepjc%2Fs10052-010-1262-1.pdf Eq 2
549  // Originally from https://arxiv.org/pdf/0710.3930.pdf
550  double a=4.138224000651535;
551  double b=1.50377557581421;
552  double c=-0.05449742257808247;
553  double d=0.0008221235048211401;
554  return a + b*sqrt(boost) + c*boost + d*pow(boost,1.5);
555 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
checkFileSG.line
line
Definition: checkFileSG.py:75
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
python.SystemOfUnits.s
int s
Definition: SystemOfUnits.py:131
SUSY::isRMeson
bool isRMeson(const T &p)
Definition: AtlasPID.h:559
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
IDTPM::R
float R(const U &p)
Definition: TrackParametersHelper.h:101
SUSY::isRBaryon
bool isRBaryon(const T &p)
Definition: AtlasPID.h:563
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
hist_file_dump.d
d
Definition: hist_file_dump.py:137
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
skel.it
it
Definition: skel.GENtoEVGEN.py:423
python.SystemOfUnits.millibarn
int millibarn
Definition: SystemOfUnits.py:77
boost
Definition: DVLIterator.h:29
beamspotman.tokens
tokens
Definition: beamspotman.py:1284
drawFromPickle.exp
exp
Definition: drawFromPickle.py:36
instance
std::map< std::string, double > instance
Definition: Run_To_Get_Tags.h:8
python.setupRTTAlg.size
int size
Definition: setupRTTAlg.py:39
TrigVtx::gamma
@ gamma
Definition: TrigParticleTable.h:26
lumiFormat.i
int i
Definition: lumiFormat.py:92
SUSY::containedQuarks
std::vector< int > containedQuarks(int p)
Definition: AtlasPID.h:618
CustomParticle
Definition: CustomParticle.h:17
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
python.ext.table_printer.table
list table
Definition: table_printer.py:81
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
SUSY::isSBaryon
bool isSBaryon(const T &p)
Definition: AtlasPID.h:584
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
a
TList * a
Definition: liststreamerinfos.cxx:10
CustomParticle.h
SUSY::isRGlueball
bool isRGlueball(const T &p)
Definition: AtlasPID.h:545
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
SUSY::isSMeson
bool isSMeson(const T &p)
Definition: AtlasPID.h:589
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
python.SystemOfUnits.ns
int ns
Definition: SystemOfUnits.py:130
str
Definition: BTagTrackIpAccessor.cxx:11
ATLAS_THREAD_SAFE
#define ATLAS_THREAD_SAFE
Definition: checker_macros.h:211
COOLRates.target
target
Definition: COOLRates.py:1106
checker_macros.h
Define macros for attributes used to control the static checker.
DerivationFramework::ClustersInCone::select
void select(const xAOD::IParticle *particle, const float coneSize, const xAOD::CaloClusterContainer *clusters, std::vector< bool > &mask)
Definition: ClustersInCone.cxx:14
GeV
#define GeV
Definition: CaloTransverseBalanceVecMon.cxx:30
python.compressB64.c
def c
Definition: compressB64.py:93
HepMCHelpers.h
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37