1 #include "Randomize.hh"
2 #include "G4ParticleTable.hh"
9 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
14 HadronicProcessHelper::HadronicProcessHelper(
const std::string &
fileName){
16 m_particleTable = G4ParticleTable::GetParticleTable();
17 m_proton = m_particleTable->FindParticle(
"proton");
18 m_neutron = m_particleTable->FindParticle(
"neutron");
22 std::ifstream processListStream (fileName.c_str());
24 while(getline(processListStream,line)){
25 std::vector<G4String> tokens;
28 m_readAndParse(line,tokens,
"#");
31 G4String incident = tokens[0];
32 G4ParticleDefinition* incidentDef = m_particleTable->FindParticle(incident);
33 G4int incidentPDG = incidentDef->GetPDGEncoding();
34 m_knownParticles[incidentDef]=
true;
36 G4String
target = tokens[1];
41 for (
size_t i = 2;
i != tokens.size();
i++){
42 G4String
part = tokens[
i];
43 if (m_particleTable->contains(part))
45 prod.push_back(m_particleTable->FindParticle(part)->GetPDGEncoding());
47 G4cout<<
"Particle: "<<part<<
" is unknown."<<G4endl;
48 G4Exception(
"HadronicProcessHelper",
"UnkownParticle", FatalException,
49 "Initialization: The reaction product list contained an unknown particle");
52 if (target ==
"proton")
54 m_protonReactionMap[incidentPDG].push_back(prod);
55 }
else if (target ==
"neutron") {
56 m_neutronReactionMap[incidentPDG].push_back(prod);
58 G4Exception(
"HadronicProcessHelper",
"IllegalTarget", FatalException,
59 "Initialization: The reaction product list contained an illegal target particle");
63 processListStream.close();
74 G4bool HadronicProcessHelper::applicabilityTester(
const G4ParticleDefinition& aPart){
75 const G4ParticleDefinition* aP = &aPart;
76 if (m_knownParticles[aP])
return true;
80 G4double HadronicProcessHelper::inclusiveCrossSection(
const G4DynamicParticle *particle,
81 const G4Element *element){
85 G4int pdgCode = particle->GetDefinition()->GetPDGEncoding();
93 G4double totalNucleonCrossSection = 0;
95 for (std::vector<G4int>::iterator it = quarks.begin();
100 if (*it == 1 || *it == 2) totalNucleonCrossSection += 12 * millibarn;
102 if (*it == 3) totalNucleonCrossSection += 6 * millibarn;
107 return totalNucleonCrossSection *
pow(element->GetN(),0.7)*1.25;
110 HadronicProcessHelper::ReactionProduct HadronicProcessHelper::finalState(
const G4DynamicParticle* incidentDynamicParticle,
111 const G4Material *material, G4ParticleDefinition*& target){
119 ReactionMap* m_reactionMap;
121 const G4ElementVector* elementVector = material->GetElementVector() ;
122 const G4double* numberOfAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
124 G4double numberOfProtons=0;
125 G4double numberOfNucleons=0;
128 for (
size_t elm=0 ;
elm < material->GetNumberOfElements() ;
elm++ )
131 numberOfProtons += numberOfAtomsPerVolume[
elm]*(*elementVector)[
elm]->GetZ();
133 numberOfNucleons += numberOfAtomsPerVolume[
elm]*(*elementVector)[
elm]->GetN();
137 if(G4UniformRand()<numberOfProtons/numberOfNucleons)
139 m_reactionMap = &m_protonReactionMap;
142 m_reactionMap = &m_neutronReactionMap;
146 const G4int incidentPDG = incidentDynamicParticle->GetDefinition()->GetPDGEncoding();
149 ReactionProductList* reactionProductList = &((*m_reactionMap)[incidentPDG]);
161 ReactionProductList goodReactionProductList;
163 for (ReactionProductList::iterator prod_it = reactionProductList->begin();
164 prod_it != reactionProductList->end();
166 G4int secondaries = prod_it->size();
168 if(m_reactionIsPossible(*prod_it,incidentDynamicParticle,target)){
170 goodReactionProductList.push_back(*prod_it);
171 if (secondaries == 2){
173 }
else if (secondaries ==3) {
176 G4cerr <<
"ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
184 if (goodReactionProductList.size()==0) G4Exception(
"HadronicProcessHelper",
"NoProcessPossible", FatalException,
185 "GetFinalState: No process could be selected from the given list.");
188 G4double prob22 = 0.15;
189 G4double prob23 = 1-prob22;
191 std::vector<G4double> probabilities;
192 std::vector<G4bool> twoToThreeFlag;
194 G4double cumulatedProbability = 0;
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);
204 cumulatedProbability += prob23/good23;
205 twoToThreeFlag.push_back(
true);
207 probabilities.push_back(cumulatedProbability);
211 for (std::vector<G4double>::iterator it = probabilities.begin();
212 it != probabilities.end();
215 *it /= cumulatedProbability;
219 G4bool selected =
false;
225 while(!selected && tries < 100){
227 G4double dice = G4UniformRand();
229 while(dice>probabilities[i] && i<numberOfReactions) i++;
231 if(twoToThreeFlag[i]) {
233 if (m_phaseSpace(goodReactionProductList[i],incidentDynamicParticle,target)>G4UniformRand()) selected =
true;
240 if(tries>=100) G4cerr<<
"Could not select process!!!!"<<G4endl;
245 if (goodReactionProductList[i].
size()==2) {
250 m_checkFraction = (1.0*m_n22)/(m_n22+m_n23);
253 return goodReactionProductList[
i];
256 G4double HadronicProcessHelper::m_reactionProductMass(
const ReactionProduct& reactionProd,
257 const G4DynamicParticle* incidentDynamicParticle,G4ParticleDefinition* target){
259 G4double incidentEnergy = incidentDynamicParticle->GetTotalEnergy();
260 G4double m_1 = incidentDynamicParticle->GetDefinition()->GetPDGMass();
261 G4double m_2 = target->GetPDGMass();
263 G4double sqrtS =
sqrt(m_1*m_1 + m_2*(m_2 + 2 * incidentEnergy));
265 G4double productsMass = 0;
267 for (ReactionProduct::const_iterator r_it = reactionProd.begin(); r_it !=reactionProd.end(); r_it++){
269 productsMass += m_particleTable->FindParticle(*r_it)->GetPDGMass();
272 return sqrtS - productsMass;
275 G4bool HadronicProcessHelper::m_reactionIsPossible(
const ReactionProduct& aReaction,
276 const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
277 if (m_reactionProductMass(aReaction,aDynamicParticle,target)>0)
return true;
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));
288 void HadronicProcessHelper::m_readAndParse(
const G4String& str,
289 std::vector<G4String>& tokens,
290 const G4String& delimiters)
297 while (G4String::npos != pos || G4String::npos != lastPos)
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);
304 tokens.push_back(temp);
306 lastPos = str.find_first_not_of(delimiters, pos);
308 pos = str.find_first_of(delimiters, lastPos);
static bool s_isRGlueball(int pdg)
static std::vector< int > s_containedQuarks(int pdg)
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)