CMS 3D CMS Logo

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