CMS 3D CMS Logo

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