10 #include "G4ParticleTable.hh" 11 #include "Randomize.hh" 17 using namespace CLHEP;
20 fParticleFactory = ptr;
22 particleTable = G4ParticleTable::GetParticleTable();
24 theProton = particleTable->FindParticle(
"proton");
25 theNeutron = particleTable->FindParticle(
"neutron");
31 std::ifstream process_stream(processDefFilePath.c_str());
37 suppressionfactor = p.
getParameter<
double>(
"reggeSuppression");
38 hadronlifetime = p.
getParameter<
double>(
"hadronLifeTime");
42 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"ProcessHelper: Read in physics parameters:" 43 <<
"\n Resonant = " << resonant <<
"\n ResonanceEnergy = " << ek_0 /
GeV 45 <<
"\n Gamma = " << gamma /
GeV <<
" GeV" 46 <<
"\n Amplitude = " << amplitude / millibarn <<
" millibarn" 47 <<
"ReggeSuppression = " << 100 * suppressionfactor <<
" %" 48 <<
"HadronLifeTime = " << hadronlifetime <<
" s" 49 <<
"ReggeModel = " << reggemodel <<
"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];
68 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"ProcessHelper: Incident " << incident <<
"; Target " <<
target;
72 for (
size_t i = 2;
i != tokens.size();
i++) {
73 G4String
part = tokens[
i];
74 if (particleTable->contains(part)) {
75 prod.push_back(particleTable->FindParticle(part)->GetPDGEncoding());
77 G4Exception(
"G4ProcessHelper",
80 "Initialization: The reaction product list contained an unknown particle");
83 if (target ==
"proton") {
84 pReactionMap[incidentPDG].push_back(prod);
85 }
else if (target ==
"neutron") {
86 nReactionMap[incidentPDG].push_back(prod);
88 G4Exception(
"G4ProcessHelper",
91 "Initialization: The reaction product list contained an illegal target particle");
95 process_stream.close();
97 for (
auto part : fParticleFactory->GetCustomParticles()) {
100 edm::LogInfo(
"SimG4CoreCustomPhysics") <<
"ProcessHelper: Lifetime of " <<
part->GetParticleName() <<
" set to " 101 << particle->GetPDGLifeTime() /
s <<
" s;" 102 <<
" isStable: " << particle->GetPDGStable();
110 const G4ParticleDefinition* aP = &aPart;
111 if (known_particles[aP])
121 G4int thePDGCode = aParticle->GetDefinition()->GetPDGEncoding();
122 double boost = (aParticle->GetKineticEnergy() + aParticle->GetMass()) / aParticle->GetMass();
124 G4double theXsec = 0;
125 G4String
name = aParticle->GetDefinition()->GetParticleName();
129 theXsec = 24 * millibarn;
133 for (std::vector<G4int>::iterator it = nq.begin(); it != nq.end(); it++) {
135 if (*it == 1 || *it == 2)
136 theXsec += 12 * millibarn;
138 theXsec += 6 * millibarn;
142 double R = Regge(boost);
143 double P = Pom(boost);
144 if (thePDGCode > 0) {
146 theXsec = (P +
R) * millibarn;
148 theXsec = 2 * P * millibarn;
150 theXsec = (R + 2 * P) * millibarn;
152 theXsec = 3 * P * millibarn;
155 theXsec = P * millibarn;
157 theXsec = (2 * (P + R) + 30 /
sqrt(boost)) * millibarn;
159 theXsec = (R + 2 * P) * millibarn;
161 theXsec = 3 * P * millibarn;
166 double e_0 = ek_0 + aParticle->GetDefinition()->GetPDGMass();
168 e_0 =
sqrt(aParticle->GetDefinition()->GetPDGMass() * aParticle->GetDefinition()->GetPDGMass() +
169 theProton->GetPDGMass() * theProton->GetPDGMass() + 2. * e_0 * theProton->GetPDGMass());
172 double sqrts =
sqrt(aParticle->GetDefinition()->GetPDGMass() * aParticle->GetDefinition()->GetPDGMass() +
173 theProton->GetPDGMass() * theProton->GetPDGMass() +
174 2 * aParticle->GetTotalEnergy() * theProton->GetPDGMass());
177 ((sqrts - e_0) * (sqrts - e_0) + (
gamma *
gamma / 4.));
179 theXsec += res_result;
184 return theXsec *
pow(anElement->GetN(), 0.7) * 1.25;
188 const G4DynamicParticle* aDynamicParticle = aTrack.GetDynamicParticle();
195 G4Material* aMaterial = aTrack.GetMaterial();
196 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
197 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
199 G4double NumberOfProtons = 0;
200 G4double NumberOfNucleons = 0;
202 for (
size_t elm = 0; elm < aMaterial->GetNumberOfElements(); elm++) {
204 NumberOfProtons += NbOfAtomsPerVolume[elm] * (*theElementVector)[elm]->GetZ();
206 NumberOfNucleons += NbOfAtomsPerVolume[elm] * (*theElementVector)[elm]->GetN();
209 if (G4UniformRand() < NumberOfProtons / NumberOfNucleons) {
210 theReactionMap = &pReactionMap;
211 theTarget = theProton;
213 theReactionMap = &nReactionMap;
214 theTarget = theNeutron;
218 G4int theIncidentPDG = aDynamicParticle->GetDefinition()->GetPDGEncoding();
221 aDynamicParticle->GetDefinition()->GetPDGCharge() == 0.) {
223 theIncidentPDG *= -1;
226 bool baryonise =
false;
228 if (reggemodel && G4UniformRand() > 0.9 &&
247 std::vector<bool> theChargeChangeList;
249 for (ReactionProductList::iterator prod_it = aReactionProductList->begin(); prod_it != aReactionProductList->end();
251 G4int secondaries = prod_it->size();
253 if (ReactionIsPossible(*prod_it, aDynamicParticle) &&
254 (!reggemodel || (baryonise && ReactionGivesBaryon(*prod_it)) ||
258 theReactionProductList.push_back(*prod_it);
259 if (secondaries == 2) {
261 }
else if (secondaries == 3) {
264 G4cerr <<
"ReactionProduct has unsupported number of secondaries: " << secondaries << G4endl;
272 if (theReactionProductList.empty())
273 G4Exception(
"G4ProcessHelper",
276 "GetFinalState: No process could be selected from the given list.");
280 int n_rps = theReactionProductList.size();
281 int select = (
int)(G4UniformRand() * n_rps);
283 return theReactionProductList[
select];
289 G4double p23 = 1 - p22;
291 std::vector<G4double> Probabilities;
292 std::vector<G4bool> TwotoThreeFlag;
294 G4double CumulatedProbability = 0;
298 for (
unsigned int i = 0;
i != theReactionProductList.size();
i++) {
299 if (theReactionProductList[
i].
size() == 2) {
300 CumulatedProbability += p22 / N22;
301 TwotoThreeFlag.push_back(
false);
303 CumulatedProbability += p23 / N23;
304 TwotoThreeFlag.push_back(
true);
306 Probabilities.push_back(CumulatedProbability);
312 for (std::vector<G4double>::iterator it = Probabilities.begin(); it != Probabilities.end(); it++) {
313 *it /= CumulatedProbability;
320 G4bool selected =
false;
326 while (!selected && tries < 100) {
328 G4double dice = G4UniformRand();
330 while (dice > Probabilities[i] && i < theReactionProductList.size()) {
337 if (!TwotoThreeFlag[i]) {
342 if (PhaseSpace(theReactionProductList[i], aDynamicParticle) > G4UniformRand())
347 if (selected && particleTable->FindParticle(theReactionProductList[i][0])->GetPDGCharge() !=
348 aDynamicParticle->GetDefinition()->GetPDGCharge()) {
355 if (G4UniformRand() < suppressionfactor)
362 G4cerr <<
"Could not select process!!!!" << G4endl;
368 if (theReactionProductList[i].
size() == 2) {
374 checkfraction = (1.0 * n_22) / (n_22 + n_23);
378 return theReactionProductList[
i];
382 const G4DynamicParticle* aDynamicParticle) {
384 G4double E_incident = aDynamicParticle->GetTotalEnergy();
387 G4double m_1 = aDynamicParticle->GetDefinition()->GetPDGMass();
388 G4double m_2 = theTarget->GetPDGMass();
390 G4double sqrts =
sqrt(m_1 * m_1 + m_2 * (m_2 + 2 * E_incident));
393 G4double M_after = 0;
394 for (ReactionProduct::const_iterator r_it = aReaction.begin(); r_it != aReaction.end(); r_it++) {
396 M_after += particleTable->FindParticle(*r_it)->GetPDGMass();
399 return sqrts - M_after;
403 const G4DynamicParticle* aDynamicParticle) {
404 if (ReactionProductMass(aReaction, aDynamicParticle) > 0)
410 for (ReactionProduct::const_iterator it = aReaction.begin(); it != aReaction.end(); it++)
417 G4double qValue = ReactionProductMass(aReaction, aDynamicParticle);
418 G4double phi =
sqrt(1 + qValue / (2 * 0.139 *
GeV)) *
pow(qValue / (1.1 *
GeV), 3. / 2.);
419 return (phi / (1 + phi));
428 while (G4String::npos != pos || G4String::npos != lastPos) {
430 G4String
temp = str.substr(lastPos, pos - lastPos);
431 while (temp.c_str()[0] ==
' ')
433 while (temp[temp.size() - 1] ==
' ')
434 temp.erase(temp.size() - 1, 1);
436 tokens.push_back(temp);
438 lastPos = str.find_first_not_of(delimiters, pos);
440 pos = str.find_first_of(delimiters, lastPos);
445 double a = 2.165635078566177;
446 double b = 0.1467453738547229;
447 double c = -0.9607903711871166;
448 return 1.5 *
exp(a + b / boost + c *
log(boost));
452 double a = 4.138224000651535;
453 double b = 1.50377557581421;
454 double c = -0.05449742257808247;
455 double d = 0.0008221235048211401;
456 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)
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)
std::vector< ReactionProduct > ReactionProductList
Power< A, B >::type pow(const A &a, const B &b)
G4bool ReactionGivesBaryon(const ReactionProduct &aReaction)
G4double Pom(const double boost)