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