1 #include "G4HadReentrentException.hh"
2 #include "G4ProcessManager.hh"
3 #include "G4ParticleTable.hh"
5 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
6 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
8 #include "SimG4Core/CustomPhysics/interface/FullModelReactionDynamics.hh"
12 FullModelHadronicProcess::FullModelHadronicProcess(G4ProcessHelper * aHelper,
const G4String& processName) :
13 G4VDiscreteProcess(processName), theHelper(aHelper)
17 FullModelHadronicProcess::~FullModelHadronicProcess(){}
19 G4bool FullModelHadronicProcess::IsApplicable(
const G4ParticleDefinition& aP)
21 return theHelper->ApplicabilityTester(aP);
24 G4double FullModelHadronicProcess::GetMicroscopicCrossSection(
const G4DynamicParticle *aParticle,
25 const G4Element *anElement,
29 G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
33 G4double HighestEnergyLimit = 10 * TeV ;
34 G4double LowestEnergyLimit = 1 * eV;
36 G4double ParticleEnergy = aParticle->GetKineticEnergy();
39 if (ParticleEnergy > HighestEnergyLimit || ParticleEnergy < LowestEnergyLimit){
48 G4double FullModelHadronicProcess::GetMeanFreePath(
const G4Track& aTrack, G4double, G4ForceCondition*)
50 G4Material *aMaterial = aTrack.GetMaterial();
51 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
54 G4int nElements = aMaterial->GetNumberOfElements();
56 const G4double *theAtomicNumDensityVector =
57 aMaterial->GetAtomicNumDensityVector();
58 G4double aTemp = aMaterial->GetTemperature();
60 for( G4int
i=0;
i<nElements; ++
i )
63 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[
i], aTemp);
64 sigma += theAtomicNumDensityVector[
i] * xSection;
71 G4VParticleChange* FullModelHadronicProcess::PostStepDoIt(
const G4Track& aTrack,
76 const G4TouchableHandle thisTouchable(aTrack.GetTouchableHandle());
79 aParticleChange.Initialize(aTrack);
81 G4DynamicParticle* IncidentRhadron =
const_cast<G4DynamicParticle*
>(aTrack.GetDynamicParticle());
83 const G4ThreeVector aPosition = aTrack.GetPosition();
85 const G4int theIncidentPDG = IncidentRhadron->GetDefinition()->GetPDGEncoding();
86 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
87 std::vector<G4ParticleDefinition*> theParticleDefinitions;
89 G4bool IncidentSurvives =
false;
90 G4bool TargetSurvives =
false;
91 G4Nucleus targetNucleus(aTrack.GetMaterial());
92 G4ParticleDefinition* outgoingRhadron=0;
93 G4ParticleDefinition* outgoingCloud=0;
94 G4ParticleDefinition* outgoingTarget=0;
98 G4ThreeVector p_0 = IncidentRhadron->GetMomentum();
102 G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
105 G4DynamicParticle* cloudParticle =
new G4DynamicParticle();
118 cloudParticle->SetDefinition(CustomIncident->
GetCloud());
120 if(cloudParticle->GetDefinition() == 0)
122 std::cout <<
"FullModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
128 double scale=cloudParticle->GetDefinition()->GetPDGMass()/IncidentRhadron->GetDefinition()->GetPDGMass();
130 G4LorentzVector cloudMomentum;
131 cloudMomentum.setVectM(IncidentRhadron->GetMomentum()*
scale,cloudParticle->GetDefinition()->GetPDGMass());
132 G4LorentzVector gluinoMomentum;
134 gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-
scale),CustomIncident->
GetSpectator()->GetPDGMass());
141 G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
142 G4LorentzVector Cloud4Momentum = cloudMomentum;
144 cloudParticle->Set4Momentum(cloudMomentum);
146 G4DynamicParticle* OrgPart = cloudParticle;
156 double E_0 = IncidentRhadron->GetKineticEnergy();
159 G4double ek = OrgPart->GetKineticEnergy();
160 G4double amas = OrgPart->GetDefinition()->GetPDGMass();
162 G4double tkin = targetNucleus.Cinema( ek );
164 OrgPart->SetKineticEnergy( ek );
165 G4double et = ek + amas;
167 G4double
pp = OrgPart->GetMomentum().mag();
170 G4ThreeVector momentum = OrgPart->GetMomentum();
171 OrgPart->SetMomentum( momentum * (p/pp) );
175 if(ek > 0.) { tkin = targetNucleus.EvaporationEffects( ek ); ek -= tkin; }
176 if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*MeV||ek<=0.) {
178 G4cout<<
"Kinetic energy is sick"<<G4endl;
179 G4cout<<
"Full R-hadron: "<<(ek+gluinoMomentum.e()-gluinoMomentum.m())/MeV<<
" MeV" <<G4endl;
180 G4cout<<
"Quark system: "<<ek/MeV<<
" MeV"<<G4endl;
182 aParticleChange.ProposeTrackStatus( fStopButAlive );
183 return &aParticleChange;
185 OrgPart->SetKineticEnergy( ek );
188 pp = OrgPart->GetMomentum().mag();
192 G4ThreeVector momentum = OrgPart->GetMomentum();
193 OrgPart->SetMomentum( momentum * (p/pp) );
200 G4ParticleDefinition* aTarget;
201 ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
202 G4bool force2to2 =
false;
204 while(rp.size()!=2 && force2to2){
205 rp = theHelper->GetFinalState(aTrack,aTarget);
207 G4int NumberOfSecondaries = rp.size();
211 G4LorentzVector Target4Momentum;
212 Target4Momentum.setVectM(0.,aTarget->GetPDGMass());
214 G4LorentzVector psum_full,psum_cloud;
215 psum_full = FullRhadron4Momentum + Target4Momentum;
216 psum_cloud = Cloud4Momentum + Target4Momentum;
217 G4ThreeVector trafo_full_cms = (-1)*psum_full.boostVector();
218 G4ThreeVector trafo_cloud_cms = (-1)*psum_cloud.boostVector();
228 for(ReactionProduct::iterator it = rp.begin();
232 G4ParticleDefinition* tempDef = theParticleTable->FindParticle(*it);
234 if (tempDef==aTarget) TargetSurvives =
true;
239 outgoingRhadron = tempDef;
248 if(outgoingCloud == 0)
250 std::cout <<
"FullModelHadronicProcess::PostStepDoIt Definition of outgoing particle cloud not available!!" << std::endl;
258 if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
260 if (tempCust==0&&rp.size()==2) outgoingTarget = tempDef;
261 if (tempDef->GetPDGEncoding()==theIncidentPDG) {
262 IncidentSurvives =
true;
264 theParticleDefinitions.push_back(tempDef);
276 if (outgoingTarget==0) outgoingTarget = theParticleTable->FindParticle(rp[1]);
290 if(IncidentSurvives) NumberOfSecondaries--;
291 aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
299 const G4HadProjectile* originalIncident =
new G4HadProjectile(*OrgPart);
303 G4HadProjectile* org2 =
new G4HadProjectile(*OrgPart);
304 G4LorentzRotation trans = org2->GetTrafoToLab();
328 G4DynamicParticle *originalTarget =
new G4DynamicParticle;
329 originalTarget->SetDefinition(aTarget);
331 G4ReactionProduct targetParticle(aTarget);
334 G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
335 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
336 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
346 G4ReactionProduct modifiedOriginal = currentParticle;
348 currentParticle.SetSide( 1 );
349 targetParticle.SetSide( -1 );
350 G4bool incidentHasChanged =
false;
351 if (!IncidentSurvives) incidentHasChanged =
true;
352 G4bool targetHasChanged =
false;
353 if (!TargetSurvives) targetHasChanged =
true;
354 G4bool quasiElastic =
false;
355 if (rp.size()==2) quasiElastic =
true;
356 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;
365 for (G4int
i = 0;
i!=NumberOfSecondaries;
i++){
366 if(theParticleDefinitions[
i]!=aTarget
367 && theParticleDefinitions[
i]!=originalIncident->GetDefinition()
368 && theParticleDefinitions[
i]!=outgoingRhadron
369 && theParticleDefinitions[
i]!=outgoingTarget)
371 G4ReactionProduct* pa =
new G4ReactionProduct;
372 pa->SetDefinition( theParticleDefinitions[
i] );
373 (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
374 vec.SetElement( vecLen++, pa );
380 if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
381 if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );
382 if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
392 CalculateMomenta( vec, vecLen,
393 originalIncident, originalTarget, modifiedOriginal,
394 targetNucleus, currentParticle, targetParticle,
395 incidentHasChanged, targetHasChanged, quasiElastic );
399 G4String cPname = currentParticle.GetDefinition()->GetParticleName();
431 aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
433 G4LorentzVector cloud_p4_new;
445 cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
446 cloud_p4_new *= trans;
448 G4LorentzVector cloud_p4_old_full = Cloud4Momentum;
449 cloud_p4_old_full.boost(trafo_full_cms);
450 G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum;
451 cloud_p4_old_cloud.boost(trafo_cloud_cms);
452 G4LorentzVector cloud_p4_full = cloud_p4_new;
453 cloud_p4_full.boost(trafo_full_cms);
454 G4LorentzVector cloud_p4_cloud = cloud_p4_new;
455 cloud_p4_cloud.boost(trafo_cloud_cms);
457 G4LorentzVector p_g_cms = gluinoMomentum;
458 p_g_cms.boost(trafo_full_cms);
460 G4LorentzVector p4_new;
464 p4_new.setVectM( cloud_p4_new.v() + gluinoMomentum.v(), outgoingRhadron->GetPDGMass() );
468 p_new = p4_new.vect();
470 aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).
m());
472 if( incidentHasChanged )
474 G4DynamicParticle* p0 =
new G4DynamicParticle;
475 p0->SetDefinition( outgoingRhadron );
476 p0->SetMomentum( p_new );
479 G4Track* Track0 =
new G4Track(p0,
480 aTrack.GetGlobalTime(),
482 Track0->SetTouchableHandle(thisTouchable);
483 aParticleChange.AddSecondary(Track0);
489 if(p0->GetKineticEnergy()>e_kin_0) {
490 G4cout<<
"ALAAAAARM!!! (incident changed from "
491 <<IncidentRhadron->GetDefinition()->GetParticleName()
492 <<
" to "<<p0->GetDefinition()->GetParticleName()<<
")"<<G4endl;
493 G4cout<<
"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/GeV<<
" GeV (should be positive)"<<G4endl;
495 if(rp.size()!=3) G4cout<<
"DOUBLE ALAAAAARM!!!"<<G4endl;
499 if(
std::abs(p0->GetKineticEnergy()-e_kin_0)>100*GeV) {
500 G4cout<<
"Diff. too big"<<G4endl;
503 aParticleChange.ProposeTrackStatus( fStopAndKill );
508 G4double p = p_new.mag();
510 aParticleChange.ProposeMomentumDirection( p_new.x()/
p, p_new.y()/
p, p_new.z()/
p );
512 aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
514 G4double aE =
sqrt(p*p+(outgoingRhadron->GetPDGMass()*outgoingRhadron->GetPDGMass()) );
515 e_kin = aE - outgoingRhadron->GetPDGMass();
535 if( targetParticle.GetMass() > 0.0 )
537 G4DynamicParticle *
p1 =
new G4DynamicParticle;
538 p1->SetDefinition( targetParticle.GetDefinition() );
540 G4ThreeVector momentum = targetParticle.GetMomentum();
541 momentum = momentum.rotate(cache,what);
542 p1->SetMomentum( momentum );
543 p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
544 G4Track* Track1 =
new G4Track(p1,
545 aTrack.GetGlobalTime(),
547 Track1->SetTouchableHandle(thisTouchable);
548 aParticleChange.AddSecondary(Track1);
550 G4DynamicParticle *pa;
558 for(
int i=0; i<vecLen; ++
i )
560 pa =
new G4DynamicParticle();
561 pa->SetDefinition( vec[i]->GetDefinition() );
562 pa->SetMomentum( vec[i]->GetMomentum() );
563 pa->Set4Momentum(trans*(pa->Get4Momentum()));
564 G4ThreeVector pvec = pa->GetMomentum();
565 G4Track* Trackn =
new G4Track(pa,
566 aTrack.GetGlobalTime(),
568 Trackn->SetTouchableHandle(thisTouchable);
569 aParticleChange.AddSecondary(Trackn);
584 const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
587 if (theRhadron!=
NULL||IncidentSurvives)
596 E_new = theRhadron->GetKineticEnergy();
613 G4LorentzVector p4_old_full = FullRhadron4Momentum;
614 p4_old_full.boost(trafo_full_cms);
615 G4LorentzVector p4_old_cloud = FullRhadron4Momentum;
616 p4_old_cloud.boost(trafo_cloud_cms);
617 G4LorentzVector p4_full = p4_new;
619 p4_full=p4_full.boost(trafo_full_cms);
621 G4LorentzVector p4_cloud = p4_new;
622 p4_cloud.boost(trafo_cloud_cms);
673 G4double AbsDeltaE = E_0-E_new;
675 if (AbsDeltaE > 10*GeV) {
676 G4cout<<
"Energy loss larger than 10 GeV..."<<G4endl;
677 G4cout<<
"E_0: "<<E_0/GeV<<
" GeV"<<G4endl;
678 G4cout<<
"E_new: "<<E_new/GeV<<
" GeV"<<G4endl;
679 G4cout<<
"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
680 G4cout<<
"x: "<<aPosition.x()/cm<<
" cm"<<G4endl;
683 delete originalIncident;
684 delete originalTarget;
689 ClearNumberOfInteractionLengthLeft();
692 return &aParticleChange;
697 void FullModelHadronicProcess::CalculateMomenta(
698 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
700 const G4HadProjectile *originalIncident,
701 const G4DynamicParticle *originalTarget,
702 G4ReactionProduct &modifiedOriginal,
703 G4Nucleus &targetNucleus,
704 G4ReactionProduct ¤tParticle,
705 G4ReactionProduct &targetParticle,
706 G4bool &incidentHasChanged,
707 G4bool &targetHasChanged,
708 G4bool quasiElastic )
710 FullModelReactionDynamics theReactionDynamics;
713 what = originalIncident->Get4Momentum().vect();
726 theReactionDynamics.TwoBody( vec, vecLen,
727 modifiedOriginal, originalTarget,
728 currentParticle, targetParticle,
729 targetNucleus, targetHasChanged );
735 G4ReactionProduct leadingStrangeParticle;
736 G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
738 leadingStrangeParticle );
745 G4bool finishedGenXPt =
false;
746 G4bool annihilation =
false;
747 if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
748 currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
752 G4double ekcor = 1.0;
753 G4double ek = originalIncident->GetKineticEnergy();
756 const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
757 if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
758 const G4double atomicWeight = targetNucleus.GetN();
759 ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
761 G4double tkin = targetNucleus.Cinema( ek );
764 modifiedOriginal.SetKineticEnergy( ekOrg );
767 const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
768 G4double rand1 = G4UniformRand();
769 G4double rand2 = G4UniformRand();
770 if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/GeV >= 1.0)) &&
771 (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
772 (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
773 (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
774 (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
775 ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
777 theReactionDynamics.GenerateXandPt( vec, vecLen,
778 modifiedOriginal, originalIncident,
779 currentParticle, targetParticle,
780 targetNucleus, incidentHasChanged,
781 targetHasChanged, leadFlag,
782 leadingStrangeParticle );
789 G4bool finishedTwoClu =
false;
790 if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
793 for(G4int i=0; i<vecLen; i++)
delete vec[i];
799 theReactionDynamics.SuppressChargedPions( vec, vecLen,
800 modifiedOriginal, currentParticle,
801 targetParticle, targetNucleus,
802 incidentHasChanged, targetHasChanged );
806 finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
807 modifiedOriginal, originalIncident,
808 currentParticle, targetParticle,
809 targetNucleus, incidentHasChanged,
810 targetHasChanged, leadFlag,
811 leadingStrangeParticle );
813 catch(G4HadReentrentException aC)
816 throw G4HadReentrentException(__FILE__, __LINE__,
"Failing to calculate momenta");
845 theReactionDynamics.TwoBody( vec, vecLen,
846 modifiedOriginal, originalTarget,
847 currentParticle, targetParticle,
848 targetNucleus, targetHasChanged );
852 G4bool FullModelHadronicProcess::MarkLeadingStrangeParticle(
853 const G4ReactionProduct ¤tParticle,
854 const G4ReactionProduct &targetParticle,
855 G4ReactionProduct &leadParticle )
863 if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
864 (currentParticle.GetDefinition() != G4Proton::Proton()) &&
865 (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
868 leadParticle = currentParticle;
870 else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
871 (targetParticle.GetDefinition() != G4Proton::Proton()) &&
872 (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
875 leadParticle = targetParticle;
880 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, G4int &vecLen)
885 for( i=0; i<vecLen; ++
i )
887 G4ThreeVector momentum = vec[
i]->GetMomentum();
888 momentum = momentum.rotate(rotation, what);
889 vec[
i]->SetMomentum(momentum);
893 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
895 G4int nsec = aParticleChange->GetNumberOfSecondaries();
896 if (nsec==0)
return 0;
898 G4bool
found =
false;
899 while (i!=nsec && !found){
902 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found =
true;
906 if(found)
return aParticleChange->GetSecondary(i)->GetDynamicParticle();
G4ParticleDefinition * GetCloud()
static bool s_isRMeson(int pdg)
static bool s_isMesonino(int pdg)
G4ParticleDefinition * GetSpectator()