CMS 3D CMS Logo

G4ProcessHelper.cc
Go to the documentation of this file.
1 #include"G4ParticleTable.hh"
2 #include "Randomize.hh"
3 
4 #include<iostream>
5 #include<fstream>
6 #include <string>
7 
11 
12 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
16 
17 using namespace CLHEP;
18 
19 G4ProcessHelper::G4ProcessHelper(const edm::ParameterSet & p){
20 
21  particleTable = G4ParticleTable::GetParticleTable();
22 
23  theProton = particleTable->FindParticle("proton");
24  theNeutron = particleTable->FindParticle("neutron");
25 
26  G4String line;
27 
28  edm::FileInPath fp = p.getParameter<edm::FileInPath>("processesDef");
29  std::string processDefFilePath = fp.fullPath();
30  std::ifstream process_stream (processDefFilePath.c_str());
31 
32  resonant = p.getParameter<bool>("resonant");
33  ek_0 = p.getParameter<double>("resonanceEnergy")*GeV;
34  gamma = p.getParameter<double>("gamma")*GeV;
35  amplitude = p.getParameter<double>("amplitude")*millibarn;
36  suppressionfactor = p.getParameter<double>("reggeSuppression");
37  hadronlifetime = p.getParameter<double>("hadronLifeTime");
38  reggemodel = p.getParameter<bool>("reggeModel");
39  mixing = p.getParameter<double>("mixing");
40 
41  edm::LogInfo("SimG4CoreCustomPhysics")
42  <<"ProcessHelper: Read in physics parameters:"
43  <<"\n Resonant = "<< resonant
44  <<"\n ResonanceEnergy = "<<ek_0/GeV<<" GeV"
45  <<"\n Gamma = "<<gamma/GeV<<" GeV"
46  <<"\n Amplitude = "<<amplitude/millibarn<<" millibarn"
47  <<"ReggeSuppression = "<<100*suppressionfactor<<" %"
48  <<"HadronLifeTime = "<<hadronlifetime<<" s"
49  <<"ReggeModel = "<< reggemodel
50  <<"Mixing = "<< mixing*100 <<" %";
51 
52  checkfraction = 0;
53  n_22 = 0;
54  n_23 = 0;
55 
56  while(getline(process_stream,line)){
57  std::vector<G4String> tokens;
58  //Getting a line
59  ReadAndParse(line,tokens,"#");
60  //Important info
61  G4String incident = tokens[0];
62 
63  G4ParticleDefinition* incidentDef = particleTable->FindParticle(incident);
64  //particleTable->DumpTable();
65  G4int incidentPDG = incidentDef->GetPDGEncoding();
66  known_particles[incidentDef]=true;
67 
68  G4String target = tokens[1];
69  edm::LogInfo("SimG4CoreCustomPhysics")
70  <<"ProcessHelper: Incident "<<incident
71  <<"; Target "<<target;
72 
73  // Making a ReactionProduct
74  ReactionProduct prod;
75  for (size_t i = 2; i != tokens.size();i++){
76  G4String part = tokens[i];
77  if (particleTable->contains(part))
78  {
79  prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
80  } else {
81  G4Exception("G4ProcessHelper", "UnkownParticle", FatalException,
82  "Initialization: The reaction product list contained an unknown particle");
83  }
84  }
85  if (target == "proton")
86  {
87  pReactionMap[incidentPDG].push_back(prod);
88  } else if (target == "neutron") {
89  nReactionMap[incidentPDG].push_back(prod);
90  } else {
91  G4Exception("G4ProcessHelper", "IllegalTarget", FatalException,
92  "Initialization: The reaction product list contained an illegal target particle");
93  }
94  }
95 
96  process_stream.close();
97 
99  CustomParticle* particle = dynamic_cast<CustomParticle*>(part);
100  if(particle) {
101  edm::LogInfo("SimG4CoreCustomPhysics")
102  <<"ProcessHelper: Lifetime of "<<part->GetParticleName()<<" set to "
103  <<particle->GetPDGLifeTime()/s<<" s;"
104  <<" isStable: "<<particle->GetPDGStable();
105  }
106  }
107 }
108 
109 G4bool G4ProcessHelper::ApplicabilityTester(const G4ParticleDefinition& aPart){
110  const G4ParticleDefinition* aP = &aPart;
111  if (known_particles[aP]) return true;
112  return false;
113 }
114 
115 G4double G4ProcessHelper::GetInclusiveCrossSection(const G4DynamicParticle *aParticle,
116  const G4Element *anElement){
117 
118  //We really do need a dedicated class to handle the cross sections. They might not always be constant
119 
120 
121  //Disassemble the PDG-code
122 
123  G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
124  double boost = (aParticle->GetKineticEnergy()+aParticle->GetMass())/aParticle->GetMass();
125  // G4cout<<"thePDGCode: "<<thePDGCode<<G4endl;
126  G4double theXsec = 0;
127  G4String name = aParticle->GetDefinition()->GetParticleName();
128  if(!reggemodel)
129  {
130  //Flat cross section
131  if(CustomPDGParser::s_isRGlueball(thePDGCode)) {
132  theXsec = 24 * millibarn;
133  } else {
134  std::vector<G4int> nq=CustomPDGParser::s_containedQuarks(thePDGCode);
135  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Number of quarks: "<<nq.size()<<G4endl;
136  for (std::vector<G4int>::iterator it = nq.begin();
137  it != nq.end();
138  it++)
139  {
140  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Quarkvector: "<<*it<<G4endl;
141  if (*it == 1 || *it == 2) theXsec += 12 * millibarn;
142  if (*it == 3) theXsec += 6 * millibarn;
143  }
144  }
145  }
146  else
147  { //reggemodel
148  double R = Regge(boost);
149  double P = Pom(boost);
150  if(thePDGCode>0)
151  {
152  if(CustomPDGParser::s_isMesonino(thePDGCode)) theXsec=(P+R)*millibarn;
153  if(CustomPDGParser::s_isSbaryon(thePDGCode)) theXsec=2*P*millibarn;
154  if(CustomPDGParser::s_isRMeson(thePDGCode)||CustomPDGParser::s_isRGlueball(thePDGCode)) theXsec=(R+2*P)*millibarn;
155  if(CustomPDGParser::s_isRBaryon(thePDGCode)) theXsec=3*P*millibarn;
156  }
157  else
158  {
159  if(CustomPDGParser::s_isMesonino(thePDGCode)) theXsec=P*millibarn;
160  if(CustomPDGParser::s_isSbaryon(thePDGCode)) theXsec=(2*(P+R)+30/sqrt(boost))*millibarn;
161  if(CustomPDGParser::s_isRMeson(thePDGCode)||CustomPDGParser::s_isRGlueball(thePDGCode)) theXsec=(R+2*P)*millibarn;
162  if(CustomPDGParser::s_isRBaryon(thePDGCode)) theXsec=3*P*millibarn;
163  }
164  }
165  //Adding resonance
166  if(resonant)
167  {
168  double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass(); //Now total energy
169 
170  e_0 = sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
171  + theProton->GetPDGMass()*theProton->GetPDGMass()
172  + 2.*e_0*theProton->GetPDGMass());
173  // edm::LogInfo("SimG4CoreCustomPhysics")<<e_0/GeV<<G4endl;
174  // edm::LogInfo("SimG4CoreCustomPhysics")<<ek_0/GeV<<" "<<aParticle->GetDefinition()->GetPDGMass()/GeV<<" "<<theProton->GetPDGMass()/GeV<<G4endl;
175  double sqrts=sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
176  + theProton->GetPDGMass()*theProton->GetPDGMass() + 2*aParticle->GetTotalEnergy()*theProton->GetPDGMass());
177 
178  double res_result = amplitude*(gamma*gamma/4.)/((sqrts-e_0)*(sqrts-e_0)+(gamma*gamma/4.));//Non-relativistic Breit Wigner
179 
180  theXsec += res_result;
181  // if(fabs(aParticle->GetKineticEnergy()/GeV-200)<10) std::cout<<sqrts/GeV<<" "<<theXsec/millibarn<<std::endl;
182  }
183 
184  // std::cout<<"Xsec/nucleon: "<<theXsec/millibarn<<"millibarn, Total Xsec: "<<theXsec * anElement->GetN()/millibarn<<" millibarn"<<std::endl;
185  return theXsec * pow(anElement->GetN(),0.7)*1.25;// * 0.523598775598299;
186 }
187 
188 ReactionProduct G4ProcessHelper::GetFinalState(const G4Track& aTrack, G4ParticleDefinition*& aTarget){
189 
190  const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
191 
192  //-----------------------------------------------
193  // Choose n / p as target
194  // and get ReactionProductList pointer
195  //-----------------------------------------------
196 
197  G4Material* aMaterial = aTrack.GetMaterial();
198  const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
199  const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
200 
201  G4double NumberOfProtons=0;
202  G4double NumberOfNucleons=0;
203 
204  for ( size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
205  {
206  //Summing number of protons per unit volume
207  NumberOfProtons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetZ();
208  //Summing nucleons (not neutrons)
209  NumberOfNucleons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetN();
210  }
211 
212  if(G4UniformRand()<NumberOfProtons/NumberOfNucleons)
213  {
214  theReactionMap = &pReactionMap;
215  theTarget = theProton;
216  } else {
217  theReactionMap = &nReactionMap;
218  theTarget = theNeutron;
219  }
220  aTarget = theTarget;
221 
222  G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
223 
224  if(reggemodel
225  &&CustomPDGParser::s_isMesonino(theIncidentPDG)
226  && G4UniformRand()*mixing>0.5
227  &&aDynamicParticle->GetDefinition()->GetPDGCharge()==0.
228  )
229  {
230  // G4cout<<"Oscillating..."<<G4endl;
231  theIncidentPDG *= -1;
232  }
233 
234 
235  bool baryonise=false;
236 
237  if(reggemodel
238  && G4UniformRand()>0.9
239  &&(
240  (CustomPDGParser::s_isMesonino(theIncidentPDG)&&theIncidentPDG>0)
241  ||
242  CustomPDGParser::s_isRMeson(theIncidentPDG)
243  )
244  )
245  baryonise=true;
246 
247  //Making a pointer directly to the ReactionProductList we are looking at. Makes life easier :-)
248  ReactionProductList* aReactionProductList = &((*theReactionMap)[theIncidentPDG]);
249 
250  //-----------------------------------------------
251  // Count processes
252  // kinematic check
253  // compute number of 2 -> 2 and 2 -> 3 processes
254  //-----------------------------------------------
255 
256  G4int N22 = 0; //Number of 2 -> 2 processes
257  G4int N23 = 0; //Number of 2 -> 3 processes. Couldn't think of more informative names
258 
259  //This is the list to be populated
260  ReactionProductList theReactionProductList;
261  std::vector<bool> theChargeChangeList;
262 
263  for (ReactionProductList::iterator prod_it = aReactionProductList->begin();
264  prod_it != aReactionProductList->end();
265  prod_it++){
266  G4int secondaries = prod_it->size();
267  // If the reaction is not possible we will not consider it
268  if(ReactionIsPossible(*prod_it,aDynamicParticle)
269  && (!reggemodel ||
270  (baryonise&&ReactionGivesBaryon(*prod_it))
271  ||
272  (!baryonise&&!ReactionGivesBaryon(*prod_it))
273  ||
274  (CustomPDGParser::s_isSbaryon(theIncidentPDG))
275  ||
276  (CustomPDGParser::s_isRBaryon(theIncidentPDG))
277  )
278  )
279  {
280  // The reaction is possible. Let's store and count it
281  theReactionProductList.push_back(*prod_it);
282  if (secondaries == 2){
283  N22++;
284  } else if (secondaries ==3) {
285  N23++;
286  } else {
287  G4cerr << "ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
288  }
289  } /*else {
290  edm::LogInfo("SimG4CoreCustomPhysics")<<"There was an impossible process"<<G4endl;
291  }*/
292  }
293  // edm::LogInfo("SimG4CoreCustomPhysics")<<"The size of the ReactionProductList is: "<<theReactionProductList.size()<<G4endl;
294 
295  if (theReactionProductList.size()==0) G4Exception("G4ProcessHelper", "NoProcessPossible", FatalException,
296  "GetFinalState: No process could be selected from the given list.");
297 
298  // For the Regge model no phase space considerations. We pick a process at random
299  if(reggemodel)
300  {
301  int n_rps = theReactionProductList.size();
302  int select = (int)( G4UniformRand()*n_rps);
303  // G4cout<<"Possible: "<<n_rps<<", chosen: "<<select<<G4endl;
304  return theReactionProductList[select];
305  }
306 
307  // Fill a probability map. Remember total probability
308  // 2->2 is 0.15*1/n_22 2->3 uses phase space
309  G4double p22 = 0.15;
310  G4double p23 = 1-p22; // :-)
311 
312  std::vector<G4double> Probabilities;
313  std::vector<G4bool> TwotoThreeFlag;
314 
315  G4double CumulatedProbability = 0;
316 
317  // To each ReactionProduct we assign a cumulated probability and a flag
318  // discerning between 2 -> 2 and 2 -> 3
319  for (unsigned int i = 0; i != theReactionProductList.size(); i++){
320  if (theReactionProductList[i].size() == 2) {
321  CumulatedProbability += p22/N22;
322  TwotoThreeFlag.push_back(false);
323  } else {
324  CumulatedProbability += p23/N23;
325  TwotoThreeFlag.push_back(true);
326  }
327  Probabilities.push_back(CumulatedProbability);
328  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Pushing back cumulated probability: "<<CumulatedProbability<<G4endl;
329  }
330 
331  //Renormalising probabilities
332  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Probs: ";
333  for (std::vector<G4double>::iterator it = Probabilities.begin();
334  it != Probabilities.end();
335  it++)
336  {
337  *it /= CumulatedProbability;
338  // edm::LogInfo("SimG4CoreCustomPhysics")<<*it<<" ";
339  }
340  // edm::LogInfo("SimG4CoreCustomPhysics")<<G4endl;
341 
342  // Choosing ReactionProduct
343 
344  G4bool selected = false;
345  G4int tries = 0;
346  // ReactionProductList::iterator prod_it;
347 
348  //Keep looping over the list until we have a choice, or until we have tried 100 times
349  unsigned int i;
350  while(!selected && tries < 100){
351  i=0;
352  G4double dice = G4UniformRand();
353  // edm::LogInfo("SimG4CoreCustomPhysics")<<"What's the dice?"<<dice<<G4endl;
354  while(dice>Probabilities[i] && i<theReactionProductList.size()){
355  // edm::LogInfo("SimG4CoreCustomPhysics")<<"i: "<<i<<G4endl;
356  i++;
357  }
358 
359  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Chosen i: "<<i<<G4endl;
360 
361  if(!TwotoThreeFlag[i]) {
362  // 2 -> 2 processes are chosen immediately
363  selected = true;
364  } else {
365  // 2 -> 3 processes require a phase space lookup
366  if (PhaseSpace(theReactionProductList[i],aDynamicParticle)>G4UniformRand()) selected = true;
367  //selected = true;
368  }
369  // double suppressionfactor=0.5;
370  if(selected&&particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()!=aDynamicParticle->GetDefinition()->GetPDGCharge())
371  {
372  /*
373  edm::LogInfo("SimG4CoreCustomPhysics")<<"Incoming particle "<<aDynamicParticle->GetDefinition()->GetParticleName()
374  <<" has charge "<<aDynamicParticle->GetDefinition()->GetPDGCharge()<<G4endl;
375  edm::LogInfo("SimG4CoreCustomPhysics")<<"Suggested particle "<<particleTable->FindParticle(theReactionProductList[i][0])->GetParticleName()
376  <<" has charge "<<particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()<<G4endl;
377  */
378  if(G4UniformRand()<suppressionfactor) selected = false;
379  }
380  tries++;
381  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Tries: "<<tries<<G4endl;
382  }
383  if(tries>=100) G4cerr<<"Could not select process!!!!"<<G4endl;
384 
385  // edm::LogInfo("SimG4CoreCustomPhysics")<<"So far so good"<<G4endl;
386  // edm::LogInfo("SimG4CoreCustomPhysics")<<"Sec's: "<<theReactionProductList[i].size()<<G4endl;
387 
388  //Updating checkfraction:
389  if (theReactionProductList[i].size()==2) {
390  n_22++;
391  } else {
392  n_23++;
393  }
394 
395  checkfraction = (1.0*n_22)/(n_22+n_23);
396  // edm::LogInfo("SimG4CoreCustomPhysics")<<"n_22: "<<n_22<<" n_23: "<<n_23<<" Checkfraction: "<<checkfraction<<G4endl;
397  // edm::LogInfo("SimG4CoreCustomPhysics") <<"Biig number: "<<n_22+n_23<<G4endl;
398  //Return the chosen ReactionProduct
399  return theReactionProductList[i];
400 }
401 
402 G4double G4ProcessHelper::ReactionProductMass(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
403  // Incident energy:
404  G4double E_incident = aDynamicParticle->GetTotalEnergy();
405  //edm::LogInfo("SimG4CoreCustomPhysics")<<"Total energy: "<<E_incident<<" Kinetic: "<<aDynamicParticle->GetKineticEnergy()<<G4endl;
406  // sqrt(s)= sqrt(m_1^2 + m_2^2 + 2 E_1 m_2)
407  G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
408  G4double m_2 = theTarget->GetPDGMass();
409  //edm::LogInfo("SimG4CoreCustomPhysics")<<"M_R: "<<m_1/GeV<<" GeV, M_np: "<<m_2/GeV<<" GeV"<<G4endl;
410  G4double sqrts = sqrt(m_1*m_1 + m_2*(m_2 + 2 * E_incident));
411  //edm::LogInfo("SimG4CoreCustomPhysics")<<"sqrt(s) = "<<sqrts/GeV<<" GeV"<<G4endl;
412  // Sum of rest masses after reaction:
413  G4double M_after = 0;
414  for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it !=aReaction.end(); r_it++){
415  //edm::LogInfo("SimG4CoreCustomPhysics")<<"Mass contrib: "<<(particleTable->FindParticle(*r_it)->GetPDGMass())/MeV<<" MeV"<<G4endl;
416  M_after += particleTable->FindParticle(*r_it)->GetPDGMass();
417  }
418  //edm::LogInfo("SimG4CoreCustomPhysics")<<"Intending to return this ReactionProductMass: "<<(sqrts - M_after)/MeV<<" MeV"<<G4endl;
419  return sqrts - M_after;
420 }
421 
422 G4bool G4ProcessHelper::ReactionIsPossible(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
423  if (ReactionProductMass(aReaction,aDynamicParticle)>0) return true;
424  return false;
425 }
426 
427 G4bool G4ProcessHelper::ReactionGivesBaryon(const ReactionProduct& aReaction){
428  for (ReactionProduct::const_iterator it = aReaction.begin();it!=aReaction.end();it++)
430  return false;
431 }
432 
433 G4double G4ProcessHelper::PhaseSpace(const ReactionProduct& aReaction,const G4DynamicParticle* aDynamicParticle){
434  G4double qValue = ReactionProductMass(aReaction,aDynamicParticle);
435  G4double phi = sqrt(1+qValue/(2*0.139*GeV))*pow(qValue/(1.1*GeV),3./2.);
436  return (phi/(1+phi));
437 }
438 
439 void G4ProcessHelper::ReadAndParse(const G4String& str,
440  std::vector<G4String>& tokens,
441  const G4String& delimiters)
442 {
443  // Skip delimiters at beginning.
444  G4String::size_type lastPos = str.find_first_not_of(delimiters, 0);
445  // Find first "non-delimiter".
446  G4String::size_type pos = str.find_first_of(delimiters, lastPos);
447 
448  while (G4String::npos != pos || G4String::npos != lastPos)
449  {
450  //Skipping leading / trailing whitespaces
451  G4String temp = str.substr(lastPos, pos - lastPos);
452  while(temp.c_str()[0] == ' ') temp.erase(0,1);
453  while(temp[temp.size()-1] == ' ') temp.erase(temp.size()-1,1);
454  // Found a token, add it to the vector.
455  tokens.push_back(temp);
456  // Skip delimiters. Note the "not_of"
457  lastPos = str.find_first_not_of(delimiters, pos);
458  // Find next "non-delimiter"
459  pos = str.find_first_of(delimiters, lastPos);
460  }
461 }
462 
463 double G4ProcessHelper::Regge(const double boost)
464 {
465  double a=2.165635078566177;
466  double b=0.1467453738547229;
467  double c=-0.9607903711871166;
468  return 1.5*exp(a+b/boost+c*log(boost));
469 }
470 
471 
472 double G4ProcessHelper::Pom(const double boost)
473 {
474  double a=4.138224000651535;
475  double b=1.50377557581421;
476  double c=-0.05449742257808247;
477  double d=0.0008221235048211401;
478  return a + b*sqrt(boost) + c*boost + d*pow(boost,1.5);
479 }
480 
size
Write out results.
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
static bool s_isRGlueball(int pdg)
Definition: CLHEP.h:16
static bool s_isSbaryon(int pdg)
static bool s_isRMeson(int pdg)
uint16_t size_type
static bool s_isMesonino(int pdg)
T sqrt(T t)
Definition: SSEVec.h:18
static bool s_isRBaryon(int pdg)
part
Definition: HCALResponse.h:20
static std::vector< int > s_containedQuarks(int pdg)
double b
Definition: hdecay.h:120
std::pair< OmniClusterRef, TrackingParticleRef > P
static const std::vector< G4ParticleDefinition * > & GetCustomParticles()
double a
Definition: hdecay.h:121
std::string fullPath() const
Definition: FileInPath.cc:184
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40