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"
11 using namespace CLHEP;
13 FullModelHadronicProcess::FullModelHadronicProcess(G4ProcessHelper * aHelper,
15 G4VDiscreteProcess(processName), theHelper(aHelper)
18 FullModelHadronicProcess::~FullModelHadronicProcess()
21 G4bool FullModelHadronicProcess::IsApplicable(
const G4ParticleDefinition& aP)
23 return theHelper->ApplicabilityTester(aP);
26 G4double FullModelHadronicProcess::GetMicroscopicCrossSection(
const G4DynamicParticle *aParticle,
27 const G4Element *anElement,
31 G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
36 G4double FullModelHadronicProcess::GetMeanFreePath(
const G4Track& aTrack, G4double,
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 )
52 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[
i], aTemp);
53 sigma += theAtomicNumDensityVector[
i] * xSection;
55 G4double res = DBL_MAX;
56 if(sigma > 0.0) { res = 1./sigma; }
60 G4VParticleChange* FullModelHadronicProcess::PostStepDoIt(
const G4Track& aTrack,
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=0;
79 G4ParticleDefinition* outgoingCloud=0;
80 G4ParticleDefinition* outgoingTarget=0;
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() == 0)
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 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;
155 ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
156 G4bool force2to2 =
false;
158 while(rp.size()!=2 && force2to2){
159 rp = theHelper->GetFinalState(aTrack,aTarget);
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;
182 outgoingRhadron = tempDef;
185 if(outgoingCloud == 0) {
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==0&&rp.size()==2) outgoingTarget = tempDef;
196 if (tempDef->GetPDGEncoding()==theIncidentPDG) {
197 IncidentSurvives =
true;
199 theParticleDefinitions.push_back(tempDef);
203 if (outgoingTarget==0) 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();
252 G4DynamicParticle *originalTarget =
new G4DynamicParticle;
253 originalTarget->SetDefinition(aTarget);
255 G4ReactionProduct targetParticle(aTarget);
257 G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
258 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
259 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
268 G4ReactionProduct modifiedOriginal = currentParticle;
270 currentParticle.SetSide( 1 );
271 targetParticle.SetSide( -1 );
272 G4bool incidentHasChanged =
false;
273 if (!IncidentSurvives) incidentHasChanged =
true;
274 G4bool targetHasChanged =
false;
275 if (!TargetSurvives) targetHasChanged =
true;
276 G4bool quasiElastic =
false;
277 if (rp.size()==2) quasiElastic =
true;
278 G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> vec;
285 for (G4int
i = 0;
i!=NumberOfSecondaries;
i++){
286 if(theParticleDefinitions[
i]!=aTarget
287 && theParticleDefinitions[
i]!=originalIncident->GetDefinition()
288 && theParticleDefinitions[
i]!=outgoingRhadron
289 && theParticleDefinitions[
i]!=outgoingTarget)
291 G4ReactionProduct*
pa =
new G4ReactionProduct;
292 pa->SetDefinition( theParticleDefinitions[
i] );
293 (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
294 vec.SetElement( vecLen++, pa );
298 if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
299 if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );
300 if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
309 CalculateMomenta( vec, vecLen,
310 originalIncident, originalTarget, modifiedOriginal,
311 targetNucleus, currentParticle, targetParticle,
312 incidentHasChanged, targetHasChanged, quasiElastic );
316 G4String cPname = currentParticle.GetDefinition()->GetParticleName();
344 aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
346 G4LorentzVector cloud_p4_new;
358 cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
359 cloud_p4_new *= trans;
361 G4LorentzVector cloud_p4_old_full = Cloud4Momentum;
362 cloud_p4_old_full.boost(trafo_full_cms);
363 G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum;
364 cloud_p4_old_cloud.boost(trafo_cloud_cms);
365 G4LorentzVector cloud_p4_full = cloud_p4_new;
366 cloud_p4_full.boost(trafo_full_cms);
367 G4LorentzVector cloud_p4_cloud = cloud_p4_new;
368 cloud_p4_cloud.boost(trafo_cloud_cms);
370 G4LorentzVector p_g_cms = gluinoMomentum;
371 p_g_cms.boost(trafo_full_cms);
373 double e = cloud_p4_new.e() + gluinoMomentum.e();
374 if(outgoingRhadron) e += outgoingRhadron->GetPDGMass();
375 G4LorentzVector p4_new( cloud_p4_new.v() + gluinoMomentum.v(),
e );
379 G4ThreeVector p_new = p4_new.vect();
381 aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).
m());
383 if( incidentHasChanged )
385 G4DynamicParticle* p0 =
new G4DynamicParticle;
386 p0->SetDefinition( outgoingRhadron );
387 p0->SetMomentum( p_new );
390 G4Track* Track0 =
new G4Track(p0,
391 aTrack.GetGlobalTime(),
393 Track0->SetTouchableHandle(thisTouchable);
394 aParticleChange.AddSecondary(Track0);
400 if(p0->GetKineticEnergy()>e_kin_0) {
401 G4cout<<
"ALAAAAARM!!! (incident changed from "
402 <<IncidentRhadron->GetDefinition()->GetParticleName()
403 <<
" to "<<p0->GetDefinition()->GetParticleName()<<
")"<<G4endl;
404 G4cout<<
"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/
GeV
405 <<
" GeV (should be positive)"<<G4endl;
407 if(rp.size()!=3) G4cout<<
"DOUBLE ALAAAAARM!!!"<<G4endl;
411 if(
std::abs(p0->GetKineticEnergy()-e_kin_0)>100*
GeV) {
412 G4cout<<
"Diff. too big"<<G4endl;
414 aParticleChange.ProposeTrackStatus( fStopAndKill );
418 G4double p = p_new.mag();
420 aParticleChange.ProposeMomentumDirection( p_new.x()/
p, p_new.y()/
p, p_new.z()/
p );
422 aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
424 G4double aE =
sqrt(p*p+(outgoingRhadron->GetPDGMass()*outgoingRhadron->GetPDGMass()) );
425 e_kin = aE - outgoingRhadron->GetPDGMass();
429 if( targetParticle.GetMass() > 0.0 )
431 G4DynamicParticle *
p1 =
new G4DynamicParticle;
432 p1->SetDefinition( targetParticle.GetDefinition() );
434 G4ThreeVector momentum = targetParticle.GetMomentum();
435 momentum = momentum.rotate(cache,what);
436 p1->SetMomentum( momentum );
437 p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
438 G4Track* Track1 =
new G4Track(p1,
439 aTrack.GetGlobalTime(),
441 Track1->SetTouchableHandle(thisTouchable);
442 aParticleChange.AddSecondary(Track1);
444 G4DynamicParticle *
pa;
450 for(
int i=0; i<vecLen; ++
i )
452 pa =
new G4DynamicParticle();
453 pa->SetDefinition( vec[i]->GetDefinition() );
454 pa->SetMomentum( vec[i]->GetMomentum() );
455 pa->Set4Momentum(trans*(pa->Get4Momentum()));
456 G4ThreeVector pvec = pa->GetMomentum();
457 G4Track* Trackn =
new G4Track(pa,
458 aTrack.GetGlobalTime(),
460 Trackn->SetTouchableHandle(thisTouchable);
461 aParticleChange.AddSecondary(Trackn);
467 const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
469 if (theRhadron!=
NULL||IncidentSurvives) {
471 if(IncidentSurvives) {
474 E_new = theRhadron->GetKineticEnergy();
490 G4LorentzVector p4_old_full = FullRhadron4Momentum;
491 p4_old_full.boost(trafo_full_cms);
492 G4LorentzVector p4_old_cloud = FullRhadron4Momentum;
493 p4_old_cloud.boost(trafo_cloud_cms);
494 G4LorentzVector p4_full = p4_new;
496 p4_full=p4_full.boost(trafo_full_cms);
498 G4LorentzVector p4_cloud = p4_new;
499 p4_cloud.boost(trafo_cloud_cms);
501 G4double AbsDeltaE = E_0-E_new;
503 if (AbsDeltaE > 10*
GeV) {
504 G4cout<<
"Energy loss larger than 10 GeV..."<<G4endl;
505 G4cout<<
"E_0: "<<E_0/
GeV<<
" GeV"<<G4endl;
506 G4cout<<
"E_new: "<<E_new/
GeV<<
" GeV"<<G4endl;
507 G4cout<<
"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
508 G4cout<<
"x: "<<aPosition.x()/cm<<
" cm"<<G4endl;
511 delete originalIncident;
512 delete originalTarget;
517 ClearNumberOfInteractionLengthLeft();
519 return &aParticleChange;
523 void FullModelHadronicProcess::CalculateMomenta(
524 G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
526 const G4HadProjectile *originalIncident,
527 const G4DynamicParticle *originalTarget,
528 G4ReactionProduct &modifiedOriginal,
529 G4Nucleus &targetNucleus,
530 G4ReactionProduct ¤tParticle,
531 G4ReactionProduct &targetParticle,
532 G4bool &incidentHasChanged,
533 G4bool &targetHasChanged,
534 G4bool quasiElastic )
536 FullModelReactionDynamics theReactionDynamics;
539 what = originalIncident->Get4Momentum().vect();
544 theReactionDynamics.TwoBody( vec, vecLen,
545 modifiedOriginal, originalTarget,
546 currentParticle, targetParticle,
547 targetNucleus, targetHasChanged );
553 G4ReactionProduct leadingStrangeParticle;
554 G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
556 leadingStrangeParticle );
561 G4bool finishedGenXPt =
false;
562 G4bool annihilation =
false;
563 if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
564 currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
568 G4double ekcor = 1.0;
569 G4double ek = originalIncident->GetKineticEnergy();
572 const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
573 if( ek > 1.0*
GeV )ekcor = 1./(ek/
GeV);
574 const G4double atomicWeight = G4double(targetNucleus.GetN_asInt());
575 ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
577 G4double tkin = targetNucleus.Cinema( ek );
580 modifiedOriginal.SetKineticEnergy( ekOrg );
583 const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
584 G4double rand1 = G4UniformRand();
585 G4double rand2 = G4UniformRand();
586 if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/
GeV >= 1.0)) &&
587 (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
588 (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
589 (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
590 (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
591 ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
593 theReactionDynamics.GenerateXandPt( vec, vecLen,
594 modifiedOriginal, originalIncident,
595 currentParticle, targetParticle,
596 targetNucleus, incidentHasChanged,
597 targetHasChanged, leadFlag,
598 leadingStrangeParticle );
605 G4bool finishedTwoClu =
false;
606 if( modifiedOriginal.GetTotalMomentum()/
MeV < 1.0 )
609 for(G4int i=0; i<vecLen; i++)
delete vec[i];
615 theReactionDynamics.SuppressChargedPions( vec, vecLen,
616 modifiedOriginal, currentParticle,
617 targetParticle, targetNucleus,
618 incidentHasChanged, targetHasChanged );
622 finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
623 modifiedOriginal, originalIncident,
624 currentParticle, targetParticle,
625 targetNucleus, incidentHasChanged,
626 targetHasChanged, leadFlag,
627 leadingStrangeParticle );
629 catch(G4HadReentrentException aC)
632 throw G4HadReentrentException(__FILE__, __LINE__,
"Failing to calculate momenta");
652 theReactionDynamics.TwoBody( vec, vecLen,
653 modifiedOriginal, originalTarget,
654 currentParticle, targetParticle,
655 targetNucleus, targetHasChanged );
658 G4bool FullModelHadronicProcess::MarkLeadingStrangeParticle(
659 const G4ReactionProduct ¤tParticle,
660 const G4ReactionProduct &targetParticle,
661 G4ReactionProduct &leadParticle )
669 if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
670 (currentParticle.GetDefinition() != G4Proton::Proton()) &&
671 (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
674 leadParticle = currentParticle;
676 else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
677 (targetParticle.GetDefinition() != G4Proton::Proton()) &&
678 (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
681 leadParticle = targetParticle;
686 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
692 for( i=0; i<vecLen; ++
i )
694 G4ThreeVector momentum = vec[
i]->GetMomentum();
695 momentum = momentum.rotate(rotation, what);
696 vec[
i]->SetMomentum(momentum);
700 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
702 G4int nsec = aParticleChange->GetNumberOfSecondaries();
703 if (nsec==0)
return 0;
705 G4bool
found =
false;
706 while (i!=nsec && !found){
709 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found =
true;
713 if(found)
return aParticleChange->GetSecondary(i)->GetDynamicParticle();
G4ParticleDefinition * GetCloud()
static bool s_isRMeson(int pdg)
static bool s_isMesonino(int pdg)
Abs< T >::type abs(const T &t)
G4ParticleDefinition * GetSpectator()