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 G4cout<<
"Particle: "<<part<<
" is unknown."<<G4endl;
50 G4Exception(
"HadronicProcessHelper",
"UnkownParticle", FatalException,
51 "Initialization: The reaction product list contained an unknown particle");
54 if (target ==
"proton")
56 m_protonReactionMap[incidentPDG].push_back(prod);
57 }
else if (target ==
"neutron") {
58 m_neutronReactionMap[incidentPDG].push_back(prod);
60 G4Exception(
"HadronicProcessHelper",
"IllegalTarget", FatalException,
61 "Initialization: The reaction product list contained an illegal target particle");
65 processListStream.close();
76 G4bool HadronicProcessHelper::applicabilityTester(
const G4ParticleDefinition& aPart){
77 const G4ParticleDefinition* aP = &aPart;
78 if (m_knownParticles[aP])
return true;
82 G4double HadronicProcessHelper::inclusiveCrossSection(
const G4DynamicParticle *particle,
83 const G4Element *element){
87 G4int pdgCode = particle->GetDefinition()->GetPDGEncoding();
95 G4double totalNucleonCrossSection = 0;
97 for (std::vector<G4int>::iterator it = quarks.begin();
102 if (*it == 1 || *it == 2) totalNucleonCrossSection += 12 * millibarn;
104 if (*it == 3) totalNucleonCrossSection += 6 * millibarn;
109 return totalNucleonCrossSection *
pow(element->GetN(),0.7)*1.25;
112 HadronicProcessHelper::ReactionProduct HadronicProcessHelper::finalState(
const G4DynamicParticle* incidentDynamicParticle,
113 const G4Material *material, G4ParticleDefinition*& target){
121 ReactionMap* m_reactionMap;
123 const G4ElementVector* elementVector = material->GetElementVector() ;
124 const G4double* numberOfAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
126 G4double numberOfProtons=0;
127 G4double numberOfNucleons=0;
130 for (
size_t elm=0 ; elm < material->GetNumberOfElements() ; elm++ )
133 numberOfProtons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetZ();
135 numberOfNucleons += numberOfAtomsPerVolume[elm]*(*elementVector)[elm]->GetN();
139 if(G4UniformRand()<numberOfProtons/numberOfNucleons)
141 m_reactionMap = &m_protonReactionMap;
144 m_reactionMap = &m_neutronReactionMap;
148 const G4int incidentPDG = incidentDynamicParticle->GetDefinition()->GetPDGEncoding();
151 ReactionProductList* reactionProductList = &((*m_reactionMap)[incidentPDG]);
163 ReactionProductList goodReactionProductList;
165 for (ReactionProductList::iterator prod_it = reactionProductList->begin();
166 prod_it != reactionProductList->end();
168 G4int secondaries = prod_it->size();
170 if(m_reactionIsPossible(*prod_it,incidentDynamicParticle,target)){
172 goodReactionProductList.push_back(*prod_it);
173 if (secondaries == 2){
175 }
else if (secondaries ==3) {
178 G4cerr <<
"ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
186 if (goodReactionProductList.size()==0) G4Exception(
"HadronicProcessHelper",
"NoProcessPossible", FatalException,
187 "GetFinalState: No process could be selected from the given list.");
190 G4double prob22 = 0.15;
191 G4double prob23 = 1-prob22;
193 std::vector<G4double> probabilities;
194 std::vector<G4bool> twoToThreeFlag;
196 G4double cumulatedProbability = 0;
200 size_t numberOfReactions = goodReactionProductList.size();
201 for (
size_t i = 0;
i != numberOfReactions;
i++){
202 if (goodReactionProductList[
i].
size() == 2) {
203 cumulatedProbability += prob22/good22;
204 twoToThreeFlag.push_back(
false);
206 cumulatedProbability += prob23/good23;
207 twoToThreeFlag.push_back(
true);
209 probabilities.push_back(cumulatedProbability);
213 for (std::vector<G4double>::iterator it = probabilities.begin();
214 it != probabilities.end();
217 *it /= cumulatedProbability;
221 G4bool selected =
false;
227 while(!selected && tries < 100){
229 G4double dice = G4UniformRand();
231 while(dice>probabilities[i] && i<numberOfReactions) i++;
233 if(twoToThreeFlag[i]) {
235 if (m_phaseSpace(goodReactionProductList[i],incidentDynamicParticle,target)>G4UniformRand()) selected =
true;
242 if(tries>=100) G4cerr<<
"Could not select process!!!!"<<G4endl;
247 if (goodReactionProductList[i].
size()==2) {
252 m_checkFraction = (1.0*m_n22)/(m_n22+m_n23);
255 return goodReactionProductList[
i];
258 G4double HadronicProcessHelper::m_reactionProductMass(
const ReactionProduct& reactionProd,
259 const G4DynamicParticle* incidentDynamicParticle,G4ParticleDefinition* target){
261 G4double incidentEnergy = incidentDynamicParticle->GetTotalEnergy();
262 G4double m_1 = incidentDynamicParticle->GetDefinition()->GetPDGMass();
263 G4double m_2 = target->GetPDGMass();
265 G4double sqrtS =
sqrt(m_1*m_1 + m_2*(m_2 + 2 * incidentEnergy));
267 G4double productsMass = 0;
269 for (ReactionProduct::const_iterator r_it = reactionProd.begin(); r_it !=reactionProd.end(); r_it++){
271 productsMass += m_particleTable->FindParticle(*r_it)->GetPDGMass();
274 return sqrtS - productsMass;
277 G4bool HadronicProcessHelper::m_reactionIsPossible(
const ReactionProduct& aReaction,
278 const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
279 if (m_reactionProductMass(aReaction,aDynamicParticle,target)>0)
return true;
283 G4double HadronicProcessHelper::m_phaseSpace(
const ReactionProduct& aReaction,
284 const G4DynamicParticle* aDynamicParticle,G4ParticleDefinition* target){
285 G4double qValue = m_reactionProductMass(aReaction,aDynamicParticle,target);
287 return (phi/(1+phi));
290 void HadronicProcessHelper::m_readAndParse(
const G4String& str,
291 std::vector<G4String>& tokens,
292 const G4String& delimiters)
299 while (G4String::npos != pos || G4String::npos != lastPos)
302 G4String
temp = str.substr(lastPos, pos - lastPos);
303 while(temp.c_str()[0] ==
' ') temp.erase(0,1);
304 while(temp[temp.size()-1] ==
' ') temp.erase(temp.size()-1,1);
306 tokens.push_back(temp);
308 lastPos = str.find_first_not_of(delimiters, pos);
310 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)