1 #include"G4ParticleTable.hh"
2 #include "Randomize.hh"
12 #include"SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
16 using namespace CLHEP;
20 particleTable = G4ParticleTable::GetParticleTable();
22 theProton = particleTable->FindParticle(
"proton");
23 theNeutron = particleTable->FindParticle(
"neutron");
29 std::ifstream process_stream (processDefFilePath.c_str());
34 amplitude = p.
getParameter<
double>(
"amplitude")*millibarn;
35 suppressionfactor = p.
getParameter<
double>(
"reggeSuppression");
36 hadronlifetime = p.
getParameter<
double>(
"hadronLifeTime");
41 <<
"ProcessHelper: Read in physics parameters:"
42 <<
"\n Resonant = "<< resonant
43 <<
"\n ResonanceEnergy = "<<ek_0/
GeV<<
" GeV"
44 <<
"\n Gamma = "<<gamma/
GeV<<
" GeV"
45 <<
"\n Amplitude = "<<amplitude/millibarn<<
" millibarn"
46 <<
"ReggeSuppression = "<<100*suppressionfactor<<
" %"
47 <<
"HadronLifeTime = "<<hadronlifetime<<
" s"
48 <<
"ReggeModel = "<< reggemodel
49 <<
"Mixing = "<< mixing*100 <<
" %";
55 while(getline(process_stream,line)){
56 std::vector<G4String> tokens;
58 ReadAndParse(line,tokens,
"#");
60 G4String incident = tokens[0];
62 G4ParticleDefinition* incidentDef = particleTable->FindParticle(incident);
64 G4int incidentPDG = incidentDef->GetPDGEncoding();
65 known_particles[incidentDef]=
true;
67 G4String
target = tokens[1];
69 <<
"ProcessHelper: Incident "<<incident
74 for (
size_t i = 2;
i != tokens.size();
i++){
75 G4String
part = tokens[
i];
76 if (particleTable->contains(part))
78 prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
80 G4Exception(
"G4ProcessHelper",
"UnkownParticle", FatalException,
81 "Initialization: The reaction product list contained an unknown particle");
84 if (target ==
"proton")
86 pReactionMap[incidentPDG].push_back(prod);
87 }
else if (target ==
"neutron") {
88 nReactionMap[incidentPDG].push_back(prod);
90 G4Exception(
"G4ProcessHelper",
"IllegalTarget", FatalException,
91 "Initialization: The reaction product list contained an illegal target particle");
95 process_stream.close();
97 G4ParticleTable::G4PTblDicIterator* theParticleIterator;
98 theParticleIterator = particleTable->GetIterator();
100 theParticleIterator->reset();
101 while( (*theParticleIterator)() ){
104 G4DecayTable*
table = theParticleIterator->value()->GetDecayTable();
105 if(particle!=0&&table!=0&&name.find(
"cloud")>name.size()&&hadronlifetime > 0)
107 particle->SetPDGLifeTime(hadronlifetime*
s);
108 particle->SetPDGStable(
false);
110 <<
"ProcessHelper: Lifetime of "<<name<<
" set to "
111 <<particle->GetPDGLifeTime()/s<<
" s;"
112 <<
" isStable: "<<particle->GetPDGStable();
115 theParticleIterator->reset();
118 G4bool G4ProcessHelper::ApplicabilityTester(
const G4ParticleDefinition& aPart){
119 const G4ParticleDefinition* aP = &aPart;
120 if (known_particles[aP])
return true;
124 G4double G4ProcessHelper::GetInclusiveCrossSection(
const G4DynamicParticle *aParticle,
125 const G4Element *anElement){
132 G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
133 double boost = (aParticle->GetKineticEnergy()+aParticle->GetMass())/aParticle->GetMass();
135 G4double theXsec = 0;
136 G4String name = aParticle->GetDefinition()->GetParticleName();
141 theXsec = 24 * millibarn;
145 for (std::vector<G4int>::iterator it = nq.begin();
150 if (*it == 1 || *it == 2) theXsec += 12 * millibarn;
151 if (*it == 3) theXsec += 6 * millibarn;
157 double R = Regge(boost);
158 double P = Pom(boost);
177 double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass();
179 e_0 =
sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
180 + theProton->GetPDGMass()*theProton->GetPDGMass()
181 + 2.*e_0*theProton->GetPDGMass());
184 double sqrts=
sqrt(aParticle->GetDefinition()->GetPDGMass()*aParticle->GetDefinition()->GetPDGMass()
185 + theProton->GetPDGMass()*theProton->GetPDGMass() + 2*aParticle->GetTotalEnergy()*theProton->GetPDGMass());
187 double res_result = amplitude*(gamma*gamma/4.)/((sqrts-e_0)*(sqrts-e_0)+(gamma*gamma/4.));
189 theXsec += res_result;
194 return theXsec *
pow(anElement->GetN(),0.7)*1.25;
197 ReactionProduct G4ProcessHelper::GetFinalState(
const G4Track& aTrack, G4ParticleDefinition*& aTarget){
199 const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
206 G4Material* aMaterial = aTrack.GetMaterial();
207 const G4ElementVector* theElementVector = aMaterial->GetElementVector() ;
208 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
210 G4double NumberOfProtons=0;
211 G4double NumberOfNucleons=0;
213 for (
size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ )
216 NumberOfProtons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetZ();
218 NumberOfNucleons += NbOfAtomsPerVolume[elm]*(*theElementVector)[elm]->GetN();
221 if(G4UniformRand()<NumberOfProtons/NumberOfNucleons)
223 theReactionMap = &pReactionMap;
224 theTarget = theProton;
226 theReactionMap = &nReactionMap;
227 theTarget = theNeutron;
231 G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
235 && G4UniformRand()*mixing>0.5
236 &&aDynamicParticle->GetDefinition()->GetPDGCharge()==0.
240 theIncidentPDG *= -1;
244 bool baryonise=
false;
247 && G4UniformRand()>0.9
257 ReactionProductList* aReactionProductList = &((*theReactionMap)[theIncidentPDG]);
269 ReactionProductList theReactionProductList;
270 std::vector<bool> theChargeChangeList;
272 for (ReactionProductList::iterator prod_it = aReactionProductList->begin();
273 prod_it != aReactionProductList->end();
275 G4int secondaries = prod_it->size();
277 if(ReactionIsPossible(*prod_it,aDynamicParticle)
279 (baryonise&&ReactionGivesBaryon(*prod_it))
281 (!baryonise&&!ReactionGivesBaryon(*prod_it))
290 theReactionProductList.push_back(*prod_it);
291 if (secondaries == 2){
293 }
else if (secondaries ==3) {
296 G4cerr <<
"ReactionProduct has unsupported number of secondaries: "<<secondaries<<G4endl;
304 if (theReactionProductList.size()==0) G4Exception(
"G4ProcessHelper",
"NoProcessPossible", FatalException,
305 "GetFinalState: No process could be selected from the given list.");
310 int n_rps = theReactionProductList.size();
311 int select = (int)( G4UniformRand()*n_rps);
313 return theReactionProductList[
select];
319 G4double p23 = 1-p22;
321 std::vector<G4double> Probabilities;
322 std::vector<G4bool> TwotoThreeFlag;
324 G4double CumulatedProbability = 0;
328 for (
unsigned int i = 0;
i != theReactionProductList.size();
i++){
329 if (theReactionProductList[
i].
size() == 2) {
330 CumulatedProbability += p22/N22;
331 TwotoThreeFlag.push_back(
false);
333 CumulatedProbability += p23/N23;
334 TwotoThreeFlag.push_back(
true);
336 Probabilities.push_back(CumulatedProbability);
342 for (std::vector<G4double>::iterator it = Probabilities.begin();
343 it != Probabilities.end();
346 *it /= CumulatedProbability;
353 G4bool selected =
false;
359 while(!selected && tries < 100){
361 G4double dice = G4UniformRand();
363 while(dice>Probabilities[i] && i<theReactionProductList.size()){
370 if(!TwotoThreeFlag[i]) {
375 if (PhaseSpace(theReactionProductList[i],aDynamicParticle)>G4UniformRand()) selected =
true;
379 if(selected&&particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge()!=aDynamicParticle->GetDefinition()->GetPDGCharge())
387 if(G4UniformRand()<suppressionfactor) selected =
false;
392 if(tries>=100) G4cerr<<
"Could not select process!!!!"<<G4endl;
398 if (theReactionProductList[i].
size()==2) {
404 checkfraction = (1.0*n_22)/(n_22+n_23);
408 return theReactionProductList[
i];
411 G4double G4ProcessHelper::ReactionProductMass(
const ReactionProduct& aReaction,
const G4DynamicParticle* aDynamicParticle){
413 G4double E_incident = aDynamicParticle->GetTotalEnergy();
416 G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
417 G4double m_2 = theTarget->GetPDGMass();
419 G4double sqrts =
sqrt(m_1*m_1 + m_2*(m_2 + 2 * E_incident));
422 G4double M_after = 0;
423 for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it !=aReaction.end(); r_it++){
425 M_after += particleTable->FindParticle(*r_it)->GetPDGMass();
428 return sqrts - M_after;
431 G4bool G4ProcessHelper::ReactionIsPossible(
const ReactionProduct& aReaction,
const G4DynamicParticle* aDynamicParticle){
432 if (ReactionProductMass(aReaction,aDynamicParticle)>0)
return true;
436 G4bool G4ProcessHelper::ReactionGivesBaryon(
const ReactionProduct& aReaction){
437 for (ReactionProduct::const_iterator it = aReaction.begin();it!=aReaction.end();it++)
442 G4double G4ProcessHelper::PhaseSpace(
const ReactionProduct& aReaction,
const G4DynamicParticle* aDynamicParticle){
443 G4double qValue = ReactionProductMass(aReaction,aDynamicParticle);
445 return (phi/(1+phi));
448 void G4ProcessHelper::ReadAndParse(
const G4String& str,
449 std::vector<G4String>& tokens,
450 const G4String& delimiters)
457 while (G4String::npos != pos || G4String::npos != lastPos)
460 G4String
temp = str.substr(lastPos, pos - lastPos);
461 while(temp.c_str()[0] ==
' ') temp.erase(0,1);
462 while(temp[temp.size()-1] ==
' ') temp.erase(temp.size()-1,1);
464 tokens.push_back(temp);
466 lastPos = str.find_first_not_of(delimiters, pos);
468 pos = str.find_first_of(delimiters, lastPos);
472 double G4ProcessHelper::Regge(
const double boost)
474 double a=2.165635078566177;
475 double b=0.1467453738547229;
476 double c=-0.9607903711871166;
477 return 1.5*
exp(a+b/boost+c*
log(boost));
481 double G4ProcessHelper::Pom(
const double boost)
483 double a=4.138224000651535;
484 double b=1.50377557581421;
485 double c=-0.05449742257808247;
486 double d=0.0008221235048211401;
487 return a + b*
sqrt(boost) + c*boost + d*
pow(boost,1.5);
T getParameter(std::string const &) const
static bool s_isRGlueball(int pdg)
static bool s_isSbaryon(int pdg)
static bool s_isRMeson(int pdg)
static bool s_isMesonino(int pdg)
static bool s_isRBaryon(int pdg)
static std::vector< int > s_containedQuarks(int pdg)
Geom::Phi< T > phi() const
std::string fullPath() const
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)