CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HadronicProcessHelper.cc
Go to the documentation of this file.
1 #include "Randomize.hh"
2 #include "G4ParticleTable.hh"
3 
4 #include <iostream>
5 #include <fstream>
6 #include <string>
7 
9 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
10 
13 
14 using namespace CLHEP;
15 
16 HadronicProcessHelper::HadronicProcessHelper(const std::string & fileName){
17 
18  m_particleTable = G4ParticleTable::GetParticleTable();
19  m_proton = m_particleTable->FindParticle("proton");
20  m_neutron = m_particleTable->FindParticle("neutron");
21 
22  G4String line;
23 
24  std::ifstream processListStream (fileName.c_str());
25 
26  while(getline(processListStream,line)){
27  std::vector<G4String> tokens;
28 
29  //Getting a line
30  m_readAndParse(line,tokens,"#");
31 
32  //Important info
33  G4String incident = tokens[0];
34  G4ParticleDefinition* incidentDef = m_particleTable->FindParticle(incident);
35  G4int incidentPDG = incidentDef->GetPDGEncoding();
36  m_knownParticles[incidentDef]=true;
37  // G4cout<<"Incident: "<<incident<<G4endl;
38  G4String target = tokens[1];
39  // G4cout<<"Target: "<<target<<G4endl;
40 
41  // Making a ReactionProduct
42  ReactionProduct prod;
43  for (size_t i = 2; i != tokens.size();i++){
44  G4String part = tokens[i];
45  if (m_particleTable->contains(part))
46  {
47  prod.push_back(m_particleTable->FindParticle(part)->GetPDGEncoding());
48  } else {
49  G4cout<<"Particle: "<<part<<" is unknown."<<G4endl;
50  G4Exception("HadronicProcessHelper", "UnkownParticle", FatalException,
51  "Initialization: The reaction product list contained an unknown particle");
52  }
53  }
54  if (target == "proton")
55  {
56  m_protonReactionMap[incidentPDG].push_back(prod);
57  } else if (target == "neutron") {
58  m_neutronReactionMap[incidentPDG].push_back(prod);
59  } else {
60  G4Exception("HadronicProcessHelper", "IllegalTarget", FatalException,
61  "Initialization: The reaction product list contained an illegal target particle");
62  }
63  }
64 
65  processListStream.close();
66 
67 
68  m_checkFraction = 0;
69  m_n22 = 0;
70  m_n23 = 0;
71 
72 
73 }
74 
75 
76 G4bool HadronicProcessHelper::applicabilityTester(const G4ParticleDefinition& aPart){
77  const G4ParticleDefinition* aP = &aPart;
78  if (m_knownParticles[aP]) return true;
79  return false;
80 }
81 
82 G4double HadronicProcessHelper::inclusiveCrossSection(const G4DynamicParticle *particle,
83  const G4Element *element){
84 
85  //We really do need a dedicated class to handle the cross sections. They might not always be constant
86 
87  G4int pdgCode = particle->GetDefinition()->GetPDGEncoding();
88 
89  //24mb for gluino-balls
90  if(CustomPDGParser::s_isRGlueball(pdgCode)) return 24 * millibarn * element->GetN();
91 
92  //get quark vector
93  std::vector<G4int> quarks=CustomPDGParser::s_containedQuarks(pdgCode);
94 
95  G4double totalNucleonCrossSection = 0;
96 
97  for (std::vector<G4int>::iterator it = quarks.begin();
98  it != quarks.end();
99  it++)
100  {
101  // 12mb for each 'up' or 'down'
102  if (*it == 1 || *it == 2) totalNucleonCrossSection += 12 * millibarn;
103  // 6mb for each 'strange'
104  if (*it == 3) totalNucleonCrossSection += 6 * millibarn;
105  }
106 
107  //Convert to xsec per nucleus
108  // return totalNucleonCrossSection * element->GetN();
109  return totalNucleonCrossSection * pow(element->GetN(),0.7)*1.25;// * 0.523598775598299;
110 }
111 
112 HadronicProcessHelper::ReactionProduct HadronicProcessHelper::finalState(const G4DynamicParticle* incidentDynamicParticle,
113  const G4Material *material, G4ParticleDefinition*& target){
114 
115 // const G4DynamicParticle* incidentDynamicParticle = track.GetDynamicParticle();
116 
117  //-----------------------------------------------
118  // Choose n / p as target
119  // and get ReactionProductList pointer
120  //-----------------------------------------------
121  ReactionMap* m_reactionMap;
122  //G4Material* material = track.GetMaterial();
123  const G4ElementVector* elementVector = material->GetElementVector() ;
124  const G4double* numberOfAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
125 
126  G4double numberOfProtons=0;
127  G4double numberOfNucleons=0;
128 
129  //Loop on elements
130  for ( size_t elm=0 ; elm < material->GetNumberOfElements() ; elm++ )
131  {
132  //Summing number of protons per unit volume
133  numberOfProtons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetZ();
134  //Summing nucleons (not neutrons)
135  numberOfNucleons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetN();
136  }
137 
138  //random decision of the target
139  if(G4UniformRand()<numberOfProtons/numberOfNucleons)
140  { //target is a proton
141  m_reactionMap = &m_protonReactionMap;
142  target = m_proton;
143  } else { //target is a neutron
144  m_reactionMap = &m_neutronReactionMap;
145  target = m_neutron;
146  }
147 
148  const G4int incidentPDG = incidentDynamicParticle->GetDefinition()->GetPDGEncoding();
149 
150  //Making a pointer directly to the ReactionProductList we are looking at. Makes life easier :-)
151  ReactionProductList* reactionProductList = &((*m_reactionMap)[incidentPDG]);
152 
153  //-----------------------------------------------
154  // Count processes
155  // kinematic check
156  // compute number of 2 -> 2 and 2 -> 3 processes
157  //-----------------------------------------------
158 
159  G4int good22 = 0; //Number of 2 -> 2 processes that are possible
160  G4int good23 = 0; //Number of 2 -> 3 processes that are possible
161 
162  //This is the list to be populated
163  ReactionProductList goodReactionProductList;
164 
165  for (ReactionProductList::iterator prod_it = reactionProductList->begin();
166  prod_it != reactionProductList->end();
167  prod_it++){
168  G4int secondaries = prod_it->size();
169  // If the reaction is not possible we will not consider it
170  if(m_reactionIsPossible(*prod_it,incidentDynamicParticle,target)){
171  // The reaction is possible. Let's store and count it
172  goodReactionProductList.push_back(*prod_it);
173  if (secondaries == 2){
174  good22++;
175  } else if (secondaries ==3) {
176  good23++;
177  } else {
178  G4cerr << "ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
179  }
180  } /*else {
181  G4cout<<"There was an impossible process"<<G4endl;
182  }*/
183  }
184  // G4cout<<"The size of the ReactionProductList is: "<<theReactionProductList.size()<<G4endl;
185 
186  if (goodReactionProductList.size()==0) G4Exception("HadronicProcessHelper", "NoProcessPossible", FatalException,
187  "GetFinalState: No process could be selected from the given list.");
188  // Fill a probability map. Remember total probability
189  // 2->2 is 0.15*1/n_22 2->3 uses phase space
190  G4double prob22 = 0.15;
191  G4double prob23 = 1-prob22; // :-)
192 
193  std::vector<G4double> probabilities;
194  std::vector<G4bool> twoToThreeFlag;
195 
196  G4double cumulatedProbability = 0;
197 
198  // To each ReactionProduct we assign a cumulated probability and a flag
199  // discerning between 2 -> 2 and 2 -> 3
200  size_t numberOfReactions = goodReactionProductList.size();
201  for (size_t i = 0; i != numberOfReactions; i++){
202  if (goodReactionProductList[i].size() == 2) {
203  cumulatedProbability += prob22/good22;
204  twoToThreeFlag.push_back(false);
205  } else {
206  cumulatedProbability += prob23/good23;
207  twoToThreeFlag.push_back(true);
208  }
209  probabilities.push_back(cumulatedProbability);
210  }
211 
212  //Normalising probabilities to 1
213  for (std::vector<G4double>::iterator it = probabilities.begin();
214  it != probabilities.end();
215  it++)
216  {
217  *it /= cumulatedProbability;
218  }
219 
220  // Choosing ReactionProduct
221  G4bool selected = false;
222  G4int tries = 0;
223  // ReactionProductList::iterator prod_it;
224 
225  //Keep looping over the list until we have a choice, or until we have tried 100 times
226  size_t i;
227  while(!selected && tries < 100){
228  i=0;
229  G4double dice = G4UniformRand();
230  //select the process using the dice
231  while(dice>probabilities[i] && i<numberOfReactions) i++;
232 
233  if(twoToThreeFlag[i]) {
234  // 2 -> 3 processes require a phase space lookup
235  if (m_phaseSpace(goodReactionProductList[i],incidentDynamicParticle,target)>G4UniformRand()) selected = true;
236  } else {
237  // 2 -> 2 processes are chosen immediately
238  selected = true;
239  }
240  tries++;
241  }
242  if(tries>=100) G4cerr<<"Could not select process!!!!"<<G4endl;
243 
244 
245  //Debugging stuff
246  //Updating checkfraction:
247  if (goodReactionProductList[i].size()==2) {
248  m_n22++;
249  } else {
250  m_n23++;
251  }
252  m_checkFraction = (1.0*m_n22)/(m_n22+m_n23);
253 
254  //return the selected productList
255  return goodReactionProductList[i];
256 }
257 
258 G4double HadronicProcessHelper::m_reactionProductMass(const ReactionProduct& reactionProd,
259  const G4DynamicParticle* incidentDynamicParticle,G4ParticleDefinition* target){
260  // Incident energy:
261  G4double incidentEnergy = incidentDynamicParticle->GetTotalEnergy();
262  G4double m_1 = incidentDynamicParticle->GetDefinition()->GetPDGMass();
263  G4double m_2 = target->GetPDGMass();
264  //square root of "s"
265  G4double sqrtS = sqrt(m_1*m_1 + m_2*(m_2 + 2 * incidentEnergy));
266  // Sum of rest masses after reaction:
267  G4double productsMass = 0;
268  //Loop on reaction producs
269  for (ReactionProduct::const_iterator r_it = reactionProd.begin(); r_it !=reactionProd.end(); r_it++){
270  //Sum the masses of the products
271  productsMass += m_particleTable->FindParticle(*r_it)->GetPDGMass();
272  }
273  //the result is square root of "s" minus the masses of the products
274  return sqrtS - productsMass;
275 }
276 
277 G4bool HadronicProcessHelper::m_reactionIsPossible(const ReactionProduct& aReaction,
278  const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
279  if (m_reactionProductMass(aReaction,aDynamicParticle,target)>0) return true;
280  return false;
281 }
282 
283 G4double HadronicProcessHelper::m_phaseSpace(const ReactionProduct& aReaction,
284  const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
285  G4double qValue = m_reactionProductMass(aReaction,aDynamicParticle,target);
286  G4double phi = sqrt(1+qValue/(2*0.139*GeV))*pow(qValue/(1.1*GeV),3./2.);
287  return (phi/(1+phi));
288 }
289 
290 void HadronicProcessHelper::m_readAndParse(const G4String& str,
291  std::vector<G4String>& tokens,
292  const G4String& delimiters)
293 {
294  // Skip delimiters at beginning.
295  G4String::size_type lastPos = str.find_first_not_of(delimiters, 0);
296  // Find first "non-delimiter".
297  G4String::size_type pos = str.find_first_of(delimiters, lastPos);
298 
299  while (G4String::npos != pos || G4String::npos != lastPos)
300  {
301  //Skipping leading / trailing whitespaces
302  G4String temp = str.substr(lastPos, pos - lastPos);
303  while(temp.c_str()[0] == ' ') temp.erase(0,1);
304  while(temp[temp.size()-1] == ' ') temp.erase(temp.size()-1,1);
305  // Found a token, add it to the vector.
306  tokens.push_back(temp);
307  // Skip delimiters. Note the "not_of"
308  lastPos = str.find_first_not_of(delimiters, pos);
309  // Find next "non-delimiter"
310  pos = str.find_first_of(delimiters, lastPos);
311  }
312 }
int i
Definition: DBlmapReader.cc:9
const double GeV
Definition: MathUtil.h:16
static bool s_isRGlueball(int pdg)
uint16_t size_type
T sqrt(T t)
Definition: SSEVec.h:48
part
Definition: HCALResponse.h:20
static std::vector< int > s_containedQuarks(int pdg)
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: DDAxes.h:10