1 #include "Randomize.hh" 2 #include "G4ParticleTable.hh" 9 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh" 14 using namespace CLHEP;
18 m_particleTable = G4ParticleTable::GetParticleTable();
19 m_proton = m_particleTable->FindParticle(
"proton");
20 m_neutron = m_particleTable->FindParticle(
"neutron");
24 std::ifstream processListStream (fileName.c_str());
26 while(getline(processListStream,line)){
27 std::vector<G4String> tokens;
30 m_readAndParse(line,tokens,
"#");
33 G4String incident = tokens[0];
34 G4ParticleDefinition* incidentDef = m_particleTable->FindParticle(incident);
35 G4int incidentPDG = incidentDef->GetPDGEncoding();
36 m_knownParticles[incidentDef]=
true;
38 G4String
target = tokens[1];
43 for (
size_t i = 2;
i != tokens.size();
i++){
44 G4String
part = tokens[
i];
45 if (m_particleTable->contains(part))
47 prod.push_back(m_particleTable->FindParticle(part)->GetPDGEncoding());
49 G4Exception(
"HadronicProcessHelper",
"UnkownParticle", FatalException,
50 "Initialization: The reaction product list contained an unknown particle");
53 if (target ==
"proton")
55 m_protonReactionMap[incidentPDG].push_back(prod);
56 }
else if (target ==
"neutron") {
57 m_neutronReactionMap[incidentPDG].push_back(prod);
59 G4Exception(
"HadronicProcessHelper",
"IllegalTarget", FatalException,
60 "Initialization: The reaction product list contained an illegal target particle");
64 processListStream.close();
71 G4bool HadronicProcessHelper::applicabilityTester(
const G4ParticleDefinition& aPart){
72 const G4ParticleDefinition* aP = &aPart;
73 if (m_knownParticles[aP])
return true;
77 G4double HadronicProcessHelper::inclusiveCrossSection(
const G4DynamicParticle *particle,
78 const G4Element *element){
82 G4int pdgCode = particle->GetDefinition()->GetPDGEncoding();
90 G4double totalNucleonCrossSection = 0;
92 for (std::vector<G4int>::iterator it = quarks.begin();
97 if (*it == 1 || *it == 2) totalNucleonCrossSection += 12 * millibarn;
99 if (*it == 3) totalNucleonCrossSection += 6 * millibarn;
104 return totalNucleonCrossSection *
pow(element->GetN(),0.7)*1.25;
107 HadronicProcessHelper::ReactionProduct HadronicProcessHelper::finalState(
const G4DynamicParticle* incidentDynamicParticle,
108 const G4Material *material, G4ParticleDefinition*& target){
116 ReactionMap* m_reactionMap;
118 const G4ElementVector* elementVector = material->GetElementVector() ;
119 const G4double* numberOfAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
121 G4double numberOfProtons=0;
122 G4double numberOfNucleons=0;
125 for (
size_t elm=0 ; elm < material->GetNumberOfElements() ; elm++ )
128 numberOfProtons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetZ();
130 numberOfNucleons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetN();
134 if(G4UniformRand()<numberOfProtons/numberOfNucleons)
136 m_reactionMap = &m_protonReactionMap;
139 m_reactionMap = &m_neutronReactionMap;
143 const G4int incidentPDG = incidentDynamicParticle->GetDefinition()->GetPDGEncoding();
146 ReactionProductList* reactionProductList = &((*m_reactionMap)[incidentPDG]);
158 ReactionProductList goodReactionProductList;
160 for (ReactionProductList::iterator prod_it = reactionProductList->begin();
161 prod_it != reactionProductList->end();
163 G4int secondaries = prod_it->size();
165 if(m_reactionIsPossible(*prod_it,incidentDynamicParticle,target)){
167 goodReactionProductList.push_back(*prod_it);
168 if (secondaries == 2){
170 }
else if (secondaries ==3) {
173 G4cerr <<
"ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
181 if (goodReactionProductList.empty()) G4Exception(
"HadronicProcessHelper",
"NoProcessPossible", FatalException,
182 "GetFinalState: No process could be selected from the given list.");
185 G4double prob22 = 0.15;
186 G4double prob23 = 1-prob22;
188 std::vector<G4double> probabilities;
189 std::vector<G4bool> twoToThreeFlag;
191 G4double cumulatedProbability = 0;
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);
201 cumulatedProbability += prob23/good23;
202 twoToThreeFlag.push_back(
true);
204 probabilities.push_back(cumulatedProbability);
208 for (std::vector<G4double>::iterator it = probabilities.begin();
209 it != probabilities.end();
212 *it /= cumulatedProbability;
216 G4bool selected =
false;
222 while(!selected && tries < 100){
224 G4double dice = G4UniformRand();
226 while(dice>probabilities[i] && i<numberOfReactions) i++;
228 if(twoToThreeFlag[i]) {
230 if (m_phaseSpace(goodReactionProductList[i],incidentDynamicParticle,target)>G4UniformRand()) selected =
true;
237 if(tries>=100)
G4cerr<<
"Could not select process!!!!"<<G4endl;
242 if (goodReactionProductList[i].
size()==2) {
247 m_checkFraction = (1.0*m_n22)/(m_n22+m_n23);
250 return goodReactionProductList[
i];
253 G4double HadronicProcessHelper::m_reactionProductMass(
const ReactionProduct& reactionProd,
254 const G4DynamicParticle* incidentDynamicParticle,G4ParticleDefinition* target){
256 G4double incidentEnergy = incidentDynamicParticle->GetTotalEnergy();
257 G4double m_1 = incidentDynamicParticle->GetDefinition()->GetPDGMass();
258 G4double m_2 = target->GetPDGMass();
260 G4double sqrtS =
sqrt(m_1*m_1 + m_2*(m_2 + 2 * incidentEnergy));
262 G4double productsMass = 0;
264 for (ReactionProduct::const_iterator r_it = reactionProd.begin(); r_it !=reactionProd.end(); r_it++){
266 productsMass += m_particleTable->FindParticle(*r_it)->GetPDGMass();
269 return sqrtS - productsMass;
272 G4bool HadronicProcessHelper::m_reactionIsPossible(
const ReactionProduct& aReaction,
273 const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
274 if (m_reactionProductMass(aReaction,aDynamicParticle,target)>0)
return true;
278 G4double HadronicProcessHelper::m_phaseSpace(
const ReactionProduct& aReaction,
279 const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
280 G4double qValue = m_reactionProductMass(aReaction,aDynamicParticle,target);
282 return (phi/(1+phi));
285 void HadronicProcessHelper::m_readAndParse(
const G4String&
str,
286 std::vector<G4String>& tokens,
287 const G4String& delimiters)
294 while (G4String::npos != pos || G4String::npos != lastPos)
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);
301 tokens.push_back(temp);
303 lastPos = str.find_first_not_of(delimiters, pos);
305 pos = str.find_first_of(delimiters, lastPos);
static bool s_isRGlueball(int pdg)
static std::vector< int > s_containedQuarks(int pdg)
Power< A, B >::type pow(const A &a, const B &b)