10 #include"G4ParticleTable.hh" 11 #include "Randomize.hh" 17 using namespace CLHEP;
21 fParticleFactory = ptr;
23 particleTable = G4ParticleTable::GetParticleTable();
25 theProton = particleTable->FindParticle(
"proton");
26 theNeutron = particleTable->FindParticle(
"neutron");
32 std::ifstream process_stream (processDefFilePath.c_str());
38 suppressionfactor = p.
getParameter<
double>(
"reggeSuppression");
39 hadronlifetime = p.
getParameter<
double>(
"hadronLifeTime");
44 <<
"ProcessHelper: Read in physics parameters:" 45 <<
"\n Resonant = "<< resonant
46 <<
"\n ResonanceEnergy = "<<ek_0/
GeV<<
" GeV" 47 <<
"\n Gamma = "<<gamma/
GeV<<
" GeV" 48 <<
"\n Amplitude = "<<amplitude/millibarn<<
" millibarn" 49 <<
"ReggeSuppression = "<<100*suppressionfactor<<
" %" 50 <<
"HadronLifeTime = "<<hadronlifetime<<
" s" 51 <<
"ReggeModel = "<< reggemodel
52 <<
"Mixing = "<< mixing*100 <<
" %";
58 while(getline(process_stream,line)){
59 std::vector<G4String> tokens;
61 ReadAndParse(line,tokens,
"#");
63 G4String incident = tokens[0];
65 G4ParticleDefinition* incidentDef = particleTable->FindParticle(incident);
67 G4int incidentPDG = incidentDef->GetPDGEncoding();
68 known_particles[incidentDef]=
true;
70 G4String
target = tokens[1];
72 <<
"ProcessHelper: Incident "<<incident
77 for (
size_t i = 2;
i != tokens.size();
i++){
78 G4String
part = tokens[
i];
79 if (particleTable->contains(part))
81 prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
83 G4Exception(
"G4ProcessHelper",
"UnkownParticle", FatalException,
84 "Initialization: The reaction product list contained an unknown particle");
87 if (target ==
"proton")
89 pReactionMap[incidentPDG].push_back(prod);
90 }
else if (target ==
"neutron") {
91 nReactionMap[incidentPDG].push_back(prod);
93 G4Exception(
"G4ProcessHelper",
"IllegalTarget", FatalException,
94 "Initialization: The reaction product list contained an illegal target particle");
98 process_stream.close();
100 for(
auto part : fParticleFactory->GetCustomParticles()) {
104 <<
"ProcessHelper: Lifetime of "<<
part->GetParticleName()<<
" set to " 105 <<particle->GetPDGLifeTime()/
s<<
" s;" 106 <<
" isStable: "<<particle->GetPDGStable();
115 const G4ParticleDefinition* aP = &aPart;
116 if (known_particles[aP])
return true;
121 const G4Element *anElement){
128 G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
129 double boost = (aParticle->GetKineticEnergy()+aParticle->GetMass())/aParticle->GetMass();
131 G4double theXsec = 0;
132 G4String
name = aParticle->GetDefinition()->GetParticleName();
137 theXsec = 24 * millibarn;
141 for (std::vector<G4int>::iterator it = nq.begin();
146 if (*it == 1 || *it == 2) theXsec += 12 * millibarn;
147 if (*it == 3) theXsec += 6 * millibarn;
153 double R = Regge(boost);
154 double P = Pom(boost);
173 double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass();
175 e_0 =
sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
176 + theProton->GetPDGMass()*theProton->GetPDGMass()
177 + 2.*e_0*theProton->GetPDGMass());
180 double sqrts=
sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
181 + theProton->GetPDGMass()*theProton->GetPDGMass() + 2*aParticle->GetTotalEnergy()*theProton->GetPDGMass());
185 theXsec += res_result;
190 return theXsec *
pow(anElement->GetN(),0.7)*1.25;
195 const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
202 G4Material* aMaterial = aTrack.GetMaterial();
203 const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
204 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
206 G4double NumberOfProtons=0;
207 G4double NumberOfNucleons=0;
209 for (
size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
212 NumberOfProtons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetZ();
214 NumberOfNucleons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetN();
217 if(G4UniformRand()<NumberOfProtons/NumberOfNucleons)
219 theReactionMap = &pReactionMap;
220 theTarget = theProton;
222 theReactionMap = &nReactionMap;
223 theTarget = theNeutron;
227 G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
231 && G4UniformRand()*
mixing>0.5
232 &&aDynamicParticle->GetDefinition()->GetPDGCharge()==0.
236 theIncidentPDG *= -1;
240 bool baryonise=
false;
243 && G4UniformRand()>0.9
266 std::vector<bool> theChargeChangeList;
268 for (ReactionProductList::iterator prod_it = aReactionProductList->begin();
269 prod_it != aReactionProductList->end();
271 G4int secondaries = prod_it->size();
273 if(ReactionIsPossible(*prod_it,aDynamicParticle)
275 (baryonise&&ReactionGivesBaryon(*prod_it))
277 (!baryonise&&!ReactionGivesBaryon(*prod_it))
286 theReactionProductList.push_back(*prod_it);
287 if (secondaries == 2){
289 }
else if (secondaries ==3) {
292 G4cerr <<
"ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
300 if (theReactionProductList.empty()) G4Exception(
"G4ProcessHelper",
"NoProcessPossible", FatalException,
301 "GetFinalState: No process could be selected from the given list.");
306 int n_rps = theReactionProductList.size();
307 int select = (
int)( G4UniformRand()*n_rps);
309 return theReactionProductList[
select];
315 G4double p23 = 1-p22;
317 std::vector<G4double> Probabilities;
318 std::vector<G4bool> TwotoThreeFlag;
320 G4double CumulatedProbability = 0;
324 for (
unsigned int i = 0;
i != theReactionProductList.size();
i++){
325 if (theReactionProductList[
i].
size() == 2) {
326 CumulatedProbability += p22/N22;
327 TwotoThreeFlag.push_back(
false);
329 CumulatedProbability += p23/N23;
330 TwotoThreeFlag.push_back(
true);
332 Probabilities.push_back(CumulatedProbability);
338 for (std::vector<G4double>::iterator it = Probabilities.begin();
339 it != Probabilities.end();
342 *it /= CumulatedProbability;
349 G4bool selected =
false;
355 while(!selected && tries < 100){
357 G4double dice = G4UniformRand();
359 while(dice>Probabilities[i] && i<theReactionProductList.size()){
366 if(!TwotoThreeFlag[i]) {
371 if (PhaseSpace(theReactionProductList[i],aDynamicParticle)>G4UniformRand()) selected =
true;
375 if(selected&&particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()!=aDynamicParticle->GetDefinition()->GetPDGCharge())
383 if(G4UniformRand()<suppressionfactor) selected =
false;
388 if(tries>=100)
G4cerr<<
"Could not select process!!!!"<<G4endl;
394 if (theReactionProductList[i].
size()==2) {
400 checkfraction = (1.0*n_22)/(n_22+n_23);
404 return theReactionProductList[
i];
409 G4double E_incident = aDynamicParticle->GetTotalEnergy();
412 G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
413 G4double m_2 = theTarget->GetPDGMass();
415 G4double sqrts =
sqrt(m_1*m_1 + m_2*(m_2 + 2 * E_incident));
418 G4double M_after = 0;
419 for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it !=aReaction.end(); r_it++){
421 M_after += particleTable->FindParticle(*r_it)->GetPDGMass();
424 return sqrts - M_after;
428 if (ReactionProductMass(aReaction,aDynamicParticle)>0)
return true;
433 for (ReactionProduct::const_iterator it = aReaction.begin();it!=aReaction.end();it++)
439 G4double qValue = ReactionProductMass(aReaction,aDynamicParticle);
440 G4double phi =
sqrt(1+qValue/(2*0.139*
GeV))*
pow(qValue/(1.1*
GeV),3./2.);
441 return (phi/(1+phi));
445 std::vector<G4String>& tokens,
446 const G4String& delimiters)
453 while (G4String::npos != pos || G4String::npos != lastPos)
456 G4String
temp = str.substr(lastPos, pos - lastPos);
457 while(temp.c_str()[0] ==
' ') temp.erase(0,1);
458 while(temp[temp.size()-1] ==
' ') temp.erase(temp.size()-1,1);
460 tokens.push_back(temp);
462 lastPos = str.find_first_not_of(delimiters, pos);
464 pos = str.find_first_of(delimiters, lastPos);
470 double a=2.165635078566177;
471 double b=0.1467453738547229;
472 double c=-0.9607903711871166;
473 return 1.5*
exp(a+b/boost+c*
log(boost));
479 double a=4.138224000651535;
480 double b=1.50377557581421;
481 double c=-0.05449742257808247;
482 double d=0.0008221235048211401;
483 return a + b*
sqrt(boost) + c*boost + d*
pow(boost,1.5);
G4double ReactionProductMass(const ReactionProduct &aReaction, const G4DynamicParticle *aDynamicParticle)
T getParameter(std::string const &) const
static bool s_isRGlueball(int pdg)
std::vector< ReactionProduct > ReactionProductList
G4bool ApplicabilityTester(const G4ParticleDefinition &aPart)
static bool s_isSbaryon(int pdg)
static bool s_isRMeson(int pdg)
G4double PhaseSpace(const ReactionProduct &aReaction, const G4DynamicParticle *aDynamicParticle)
std::vector< G4int > ReactionProduct
static bool s_isMesonino(int pdg)
void ReadAndParse(const G4String &str, std::vector< G4String > &tokens, const G4String &delimiters=" ")
ReactionProduct GetFinalState(const G4Track &aTrack, G4ParticleDefinition *&aTarget)
G4double GetInclusiveCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement)
static bool s_isRBaryon(int pdg)
static std::vector< int > s_containedQuarks(int pdg)
std::pair< OmniClusterRef, TrackingParticleRef > P
G4double Regge(const double boost)
G4ProcessHelper(const edm::ParameterSet &p, CustomParticleFactory *ptr)
std::string fullPath() const
G4bool ReactionIsPossible(const ReactionProduct &aReaction, const G4DynamicParticle *aDynamicParticle)
Power< A, B >::type pow(const A &a, const B &b)
G4bool ReactionGivesBaryon(const ReactionProduct &aReaction)
G4double Pom(const double boost)