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