1 #include "G4HadReentrentException.hh" 2 #include "G4ProcessManager.hh" 3 #include "G4ParticleTable.hh" 11 using namespace CLHEP;
15 G4VDiscreteProcess(processName), theHelper(aHelper)
27 const G4Element *anElement,
39 G4Material *aMaterial = aTrack.GetMaterial();
40 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
43 G4int nElements = aMaterial->GetNumberOfElements();
45 const G4double *theAtomicNumDensityVector =
46 aMaterial->GetAtomicNumDensityVector();
47 G4double aTemp = aMaterial->GetTemperature();
49 for( G4int
i=0;
i<nElements; ++
i )
53 sigma += theAtomicNumDensityVector[
i] * xSection;
55 G4double
res = DBL_MAX;
56 if(sigma > 0.0) { res = 1./sigma; }
64 const G4TouchableHandle& thisTouchable(aTrack.GetTouchableHandle());
67 aParticleChange.Initialize(aTrack);
68 const G4DynamicParticle* IncidentRhadron = aTrack.GetDynamicParticle();
70 const G4ThreeVector& aPosition = aTrack.GetPosition();
72 const G4int theIncidentPDG = IncidentRhadron->GetDefinition()->GetPDGEncoding();
73 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
74 std::vector<G4ParticleDefinition*> theParticleDefinitions;
75 G4bool IncidentSurvives =
false;
76 G4bool TargetSurvives =
false;
77 G4Nucleus targetNucleus(aTrack.GetMaterial());
78 G4ParticleDefinition* outgoingRhadron=
nullptr;
79 G4ParticleDefinition* outgoingCloud=
nullptr;
80 G4ParticleDefinition* outgoingTarget=
nullptr;
82 G4ThreeVector p_0 = IncidentRhadron->GetMomentum();
83 G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
86 G4DynamicParticle* cloudParticle =
new G4DynamicParticle();
93 cloudParticle->SetDefinition(CustomIncident->
GetCloud());
95 if(cloudParticle->GetDefinition() ==
nullptr)
97 G4cout <<
"FullModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << G4endl;
104 double scale=cloudParticle->GetDefinition()->GetPDGMass()
105 /IncidentRhadron->GetDefinition()->GetPDGMass();
107 G4LorentzVector cloudMomentum(IncidentRhadron->GetMomentum()*
scale,
108 cloudParticle->GetDefinition()->GetPDGMass());
109 G4LorentzVector gluinoMomentum(IncidentRhadron->GetMomentum()*(1.-
scale),
113 G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
114 const G4LorentzVector& Cloud4Momentum = cloudMomentum;
116 cloudParticle->Set4Momentum(cloudMomentum);
118 G4DynamicParticle* OrgPart = cloudParticle;
130 double E_0 = IncidentRhadron->GetKineticEnergy();
131 G4double ek = OrgPart->GetKineticEnergy();
132 G4double amas = OrgPart->GetDefinition()->GetPDGMass();
133 G4ThreeVector
dir = (OrgPart->GetMomentum()).
unit();
134 G4double tkin = targetNucleus.Cinema( ek );
138 tkin = targetNucleus.EvaporationEffects( ek );
141 if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*
MeV||ek<=0.) {
143 G4cout<<
"Kinetic energy is sick"<<G4endl;
144 G4cout<<
"Full R-hadron: "<<(ek+gluinoMomentum.e()-gluinoMomentum.m())/
MeV<<
" MeV" <<G4endl;
145 G4cout<<
"Quark system: "<<ek/
MeV<<
" MeV"<<G4endl;
146 aParticleChange.ProposeTrackStatus( fStopButAlive );
147 return &aParticleChange;
149 OrgPart->SetKineticEnergy( ek );
151 OrgPart->SetMomentum(dir * p);
154 G4ParticleDefinition* aTarget;
156 G4bool force2to2 =
false;
158 while(rp.size()!=2 && force2to2){
161 G4int NumberOfSecondaries = rp.size();
165 G4LorentzVector Target4Momentum(0.,0.,0.,aTarget->GetPDGMass());
167 G4LorentzVector psum_full = FullRhadron4Momentum + Target4Momentum;
168 G4LorentzVector psum_cloud = Cloud4Momentum + Target4Momentum;
169 G4ThreeVector trafo_full_cms = (-1)*psum_full.boostVector();
170 G4ThreeVector trafo_cloud_cms = (-1)*psum_cloud.boostVector();
174 for(ReactionProduct::iterator it = rp.begin();
175 it != rp.end(); ++it) {
176 G4ParticleDefinition* tempDef = theParticleTable->FindParticle(*it);
178 if (tempDef==aTarget) TargetSurvives =
true;
181 if (tempCust!=
nullptr) {
182 outgoingRhadron = tempDef;
185 if(outgoingCloud ==
nullptr) {
186 G4cout <<
"FullModelHadronicProcess::PostStepDoIt Definition of outgoing particle cloud not available!" << G4endl;
194 if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
195 if (tempCust==
nullptr&&rp.size()==2) outgoingTarget = tempDef;
196 if (tempDef->GetPDGEncoding()==theIncidentPDG) {
197 IncidentSurvives =
true;
199 theParticleDefinitions.push_back(tempDef);
203 if (outgoingTarget==
nullptr) outgoingTarget = theParticleTable->FindParticle(rp[1]);
217 if(IncidentSurvives) NumberOfSecondaries--;
218 aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
225 const G4HadProjectile* originalIncident =
new G4HadProjectile(*OrgPart);
228 G4HadProjectile* org2 =
new G4HadProjectile(*OrgPart);
229 G4LorentzRotation trans = org2->GetTrafoToLab();
234 G4DynamicParticle *originalTarget =
new G4DynamicParticle;
235 originalTarget->SetDefinition(aTarget);
237 G4ReactionProduct targetParticle(aTarget);
239 G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
240 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
241 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
250 G4ReactionProduct modifiedOriginal = currentParticle;
252 currentParticle.SetSide( 1 );
253 targetParticle.SetSide( -1 );
254 G4bool incidentHasChanged =
false;
255 if (!IncidentSurvives) incidentHasChanged =
true;
256 G4bool targetHasChanged =
false;
257 if (!TargetSurvives) targetHasChanged =
true;
258 G4bool quasiElastic =
false;
259 if (rp.size()==2) quasiElastic =
true;
260 G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> vec;
267 for (G4int
i = 0;
i!=NumberOfSecondaries;
i++){
268 if(theParticleDefinitions[
i]!=aTarget
269 && theParticleDefinitions[
i]!=originalIncident->GetDefinition()
270 && theParticleDefinitions[
i]!=outgoingRhadron
271 && theParticleDefinitions[
i]!=outgoingTarget)
273 G4ReactionProduct*
pa =
new G4ReactionProduct;
274 pa->SetDefinition( theParticleDefinitions[
i] );
275 (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
276 vec.SetElement( vecLen++, pa );
280 if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
281 if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );
282 if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
292 originalIncident, originalTarget, modifiedOriginal,
293 targetNucleus, currentParticle, targetParticle,
294 incidentHasChanged, targetHasChanged, quasiElastic );
298 G4String cPname = currentParticle.GetDefinition()->GetParticleName();
325 aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
327 G4LorentzVector cloud_p4_new;
339 cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
340 cloud_p4_new *= trans;
342 G4LorentzVector cloud_p4_old_full = Cloud4Momentum;
343 cloud_p4_old_full.boost(trafo_full_cms);
344 G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum;
345 cloud_p4_old_cloud.boost(trafo_cloud_cms);
346 G4LorentzVector cloud_p4_full = cloud_p4_new;
347 cloud_p4_full.boost(trafo_full_cms);
348 G4LorentzVector cloud_p4_cloud = cloud_p4_new;
349 cloud_p4_cloud.boost(trafo_cloud_cms);
351 G4LorentzVector p_g_cms = gluinoMomentum;
352 p_g_cms.boost(trafo_full_cms);
354 double e = cloud_p4_new.e() + gluinoMomentum.e();
355 if(outgoingRhadron) e += outgoingRhadron->GetPDGMass();
356 G4LorentzVector p4_new( cloud_p4_new.v() + gluinoMomentum.v(),
e );
360 G4ThreeVector p_new = p4_new.vect();
362 aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).
m());
364 if( incidentHasChanged )
366 G4DynamicParticle* p0 =
new G4DynamicParticle;
367 p0->SetDefinition( outgoingRhadron );
368 p0->SetMomentum( p_new );
371 G4Track* Track0 =
new G4Track(p0,
372 aTrack.GetGlobalTime(),
374 Track0->SetTouchableHandle(thisTouchable);
375 aParticleChange.AddSecondary(Track0);
381 if(p0->GetKineticEnergy()>e_kin_0) {
382 G4cout<<
"ALAAAAARM!!! (incident changed from " 383 <<IncidentRhadron->GetDefinition()->GetParticleName()
384 <<
" to "<<p0->GetDefinition()->GetParticleName()<<
")"<<G4endl;
385 G4cout<<
"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/
GeV 386 <<
" GeV (should be positive)"<<G4endl;
388 if(rp.size()!=3)
G4cout<<
"DOUBLE ALAAAAARM!!!"<<G4endl;
392 if(
std::abs(p0->GetKineticEnergy()-e_kin_0)>100*
GeV) {
393 G4cout<<
"Diff. too big"<<G4endl;
395 aParticleChange.ProposeTrackStatus( fStopAndKill );
399 G4double p = p_new.mag();
401 aParticleChange.ProposeMomentumDirection( p_new.x()/
p, p_new.y()/
p, p_new.z()/
p );
403 aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
407 if( targetParticle.GetMass() > 0.0 )
409 G4DynamicParticle *
p1 =
new G4DynamicParticle;
410 p1->SetDefinition( targetParticle.GetDefinition() );
412 G4ThreeVector momentum = targetParticle.GetMomentum();
414 p1->SetMomentum( momentum );
415 p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
416 G4Track* Track1 =
new G4Track(p1,
417 aTrack.GetGlobalTime(),
419 Track1->SetTouchableHandle(thisTouchable);
420 aParticleChange.AddSecondary(Track1);
422 G4DynamicParticle *
pa;
428 for(
int i=0;
i<vecLen; ++
i )
430 pa =
new G4DynamicParticle();
431 pa->SetDefinition( vec[
i]->GetDefinition() );
432 pa->SetMomentum( vec[
i]->GetMomentum() );
433 pa->Set4Momentum(trans*(pa->Get4Momentum()));
434 G4ThreeVector pvec = pa->GetMomentum();
435 G4Track* Trackn =
new G4Track(pa,
436 aTrack.GetGlobalTime(),
438 Trackn->SetTouchableHandle(thisTouchable);
439 aParticleChange.AddSecondary(Trackn);
445 const G4DynamicParticle* theRhadron =
FindRhadron(&aParticleChange);
447 if (theRhadron!=
nullptr||IncidentSurvives) {
449 if(IncidentSurvives) {
452 E_new = theRhadron->GetKineticEnergy();
468 G4LorentzVector p4_old_full = FullRhadron4Momentum;
469 p4_old_full.boost(trafo_full_cms);
470 G4LorentzVector p4_old_cloud = FullRhadron4Momentum;
471 p4_old_cloud.boost(trafo_cloud_cms);
472 G4LorentzVector p4_full = p4_new;
474 p4_full=p4_full.boost(trafo_full_cms);
476 G4LorentzVector p4_cloud = p4_new;
477 p4_cloud.boost(trafo_cloud_cms);
479 G4double AbsDeltaE = E_0-E_new;
481 if (AbsDeltaE > 10*
GeV) {
482 G4cout<<
"Energy loss larger than 10 GeV..."<<G4endl;
483 G4cout<<
"E_0: "<<E_0/
GeV<<
" GeV"<<G4endl;
484 G4cout<<
"E_new: "<<E_new/
GeV<<
" GeV"<<G4endl;
485 G4cout<<
"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
486 G4cout<<
"x: "<<aPosition.x()/cm<<
" cm"<<G4endl;
489 delete originalIncident;
490 delete originalTarget;
495 ClearNumberOfInteractionLengthLeft();
497 return &aParticleChange;
502 G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
504 const G4HadProjectile *originalIncident,
505 const G4DynamicParticle *originalTarget,
506 G4ReactionProduct &modifiedOriginal,
507 G4Nucleus &targetNucleus,
508 G4ReactionProduct ¤tParticle,
509 G4ReactionProduct &targetParticle,
510 G4bool &incidentHasChanged,
511 G4bool &targetHasChanged,
512 G4bool quasiElastic )
517 what = originalIncident->Get4Momentum().vect();
522 theReactionDynamics.
TwoBody( vec, vecLen,
523 modifiedOriginal, originalTarget,
524 currentParticle, targetParticle,
525 targetNucleus, targetHasChanged );
531 G4ReactionProduct leadingStrangeParticle;
534 leadingStrangeParticle );
539 G4bool finishedGenXPt =
false;
540 G4bool annihilation =
false;
541 if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
542 currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
546 G4double ekcor = 1.0;
547 G4double ek = originalIncident->GetKineticEnergy();
550 const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
551 if( ek > 1.0*
GeV )ekcor = 1./(ek/
GeV);
552 const G4double atomicWeight = G4double(targetNucleus.GetN_asInt());
553 ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
555 G4double tkin = targetNucleus.Cinema( ek );
558 modifiedOriginal.SetKineticEnergy( ekOrg );
561 const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
562 G4double rand1 = G4UniformRand();
563 G4double rand2 = G4UniformRand();
564 if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/
GeV >= 1.0)) &&
565 (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
566 (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
567 (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
568 (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
569 ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
572 modifiedOriginal, originalIncident,
573 currentParticle, targetParticle,
574 targetNucleus, incidentHasChanged,
575 targetHasChanged, leadFlag,
576 leadingStrangeParticle );
583 G4bool finishedTwoClu =
false;
584 if( modifiedOriginal.GetTotalMomentum()/
MeV < 1.0 )
587 for(G4int
i=0;
i<vecLen;
i++)
delete vec[
i];
594 modifiedOriginal, currentParticle,
595 targetParticle, targetNucleus,
596 incidentHasChanged, targetHasChanged );
600 finishedTwoClu = theReactionDynamics.
TwoCluster( vec, vecLen,
601 modifiedOriginal, originalIncident,
602 currentParticle, targetParticle,
603 targetNucleus, incidentHasChanged,
604 targetHasChanged, leadFlag,
605 leadingStrangeParticle );
607 catch(G4HadReentrentException &aC)
610 throw G4HadReentrentException(__FILE__, __LINE__,
"Failing to calculate momenta");
630 theReactionDynamics.
TwoBody( vec, vecLen,
631 modifiedOriginal, originalTarget,
632 currentParticle, targetParticle,
633 targetNucleus, targetHasChanged );
637 const G4ReactionProduct ¤tParticle,
638 const G4ReactionProduct &targetParticle,
639 G4ReactionProduct &leadParticle )
647 if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
648 (currentParticle.GetDefinition() != G4Proton::Proton()) &&
649 (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
652 leadParticle = currentParticle;
654 else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
655 (targetParticle.GetDefinition() != G4Proton::Proton()) &&
656 (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
659 leadParticle = targetParticle;
670 for( i=0; i<vecLen; ++
i )
672 G4ThreeVector momentum = vec[
i]->GetMomentum();
673 momentum = momentum.rotate(rotation,
what);
674 vec[
i]->SetMomentum(momentum);
680 G4int nsec = aParticleChange->GetNumberOfSecondaries();
681 if (nsec==0)
return nullptr;
683 G4bool
found =
false;
684 while (i!=nsec && !found){
687 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=
nullptr) found =
true;
691 if(found)
return aParticleChange->GetSecondary(i)->GetDynamicParticle();
G4bool MarkLeadingStrangeParticle(const G4ReactionProduct ¤tParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
G4bool ApplicabilityTester(const G4ParticleDefinition &aPart)
void TwoBody(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4DynamicParticle *originalTarget, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &targetHasChanged)
G4ParticleDefinition * GetCloud()
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
static bool s_isRMeson(int pdg)
virtual G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, G4double aTemp)
std::vector< G4int > ReactionProduct
static bool s_isMesonino(int pdg)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
~FullModelHadronicProcess() override
Abs< T >::type abs(const T &t)
ReactionProduct GetFinalState(const G4Track &aTrack, G4ParticleDefinition *&aTarget)
G4bool TwoCluster(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
G4double GetInclusiveCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement)
G4ParticleDefinition * GetSpectator()
G4bool GenerateXandPt(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
void CalculateMomenta(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
FullModelHadronicProcess(G4ProcessHelper *aHelper, const G4String &processName="FullModelHadronicProcess")
G4ProcessHelper * theHelper
void SuppressChargedPions(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, const G4ReactionProduct &modifiedOriginal, G4ReactionProduct ¤tParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged)
void Rotate(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen)
const G4DynamicParticle * FindRhadron(G4ParticleChange *)
G4bool IsApplicable(const G4ParticleDefinition &aP) override