2 #include "G4ProcessManager.hh" 3 #include "G4ParticleTable.hh" 5 #include "G4FermiPhaseSpaceDecay.hh" 7 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh" 10 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh" 15 using namespace CLHEP;
17 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper,
19 G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(
true)
21 m_decay =
new G4FermiPhaseSpaceDecay();
24 G4bool ToyModelHadronicProcess::IsApplicable(
const G4ParticleDefinition& aP)
26 return m_helper->applicabilityTester(aP);
29 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(
const G4DynamicParticle *particle,
30 const G4Element *element,
34 G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
36 if(m_verboseLevel >= 3) {
37 G4cout <<
"ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection
40 return inclusiveCrossSection;
43 G4double ToyModelHadronicProcess::GetMeanFreePath(
const G4Track& aTrack, G4double,
46 G4Material *aMaterial = aTrack.GetMaterial();
47 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
50 G4int nElements = aMaterial->GetNumberOfElements();
52 const G4double *theAtomicNumDensityVector =
53 aMaterial->GetAtomicNumDensityVector();
54 G4double aTemp = aMaterial->GetTemperature();
56 for( G4int
i=0;
i<nElements; ++
i )
60 sigma += theAtomicNumDensityVector[
i] * xSection;
63 if(sigma > 0.0) { x = 1./sigma; }
67 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(
const G4Track&
track,
70 const G4TouchableHandle thisTouchable(track.GetTouchableHandle());
73 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
74 m_particleChange.Initialize(track);
77 const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle();
78 const G4ThreeVector aPosition = track.GetPosition();
83 G4DynamicParticle* cloudParticle =
new G4DynamicParticle();
86 if(CustomIncident==0) {
87 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!" 90 cloudParticle->SetDefinition(CustomIncident->
GetCloud());
94 double scale = cloudParticle->GetDefinition()->GetPDGMass()
95 /incidentRHadron->GetDefinition()->GetPDGMass();
96 G4LorentzVector cloudMomentum(incidentRHadron->GetMomentum()*
scale,
97 cloudParticle->GetDefinition()->GetPDGMass());
98 G4LorentzVector gluinoMomentum(incidentRHadron->GetMomentum()*(1.-
scale),
101 cloudParticle->Set4Momentum(cloudMomentum);
103 const G4DynamicParticle *incidentParticle;
104 if(! m_detachCloud) { incidentParticle = incidentRHadron; }
105 else { incidentParticle = cloudParticle; }
107 if(m_verboseLevel >= 3)
109 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt After scaling " << G4endl;
110 G4cout <<
" RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()
111 <<
" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() /
GeV 112 <<
" P= " << incidentRHadron->Get4Momentum() /
GeV<<
" mass= " 113 << incidentRHadron->Get4Momentum().m() /
GeV <<G4endl;
114 G4cout <<
" Cloud "<<cloudParticle->GetDefinition()->GetParticleName()
115 <<
" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() /
GeV 116 <<
" P= " << cloudParticle->Get4Momentum() /
GeV<<
" mass= " 117 << cloudParticle->Get4Momentum().m() /
GeV <<G4endl;
118 G4cout <<
" Incident pdgmass= " 119 << incidentParticle->GetDefinition()->GetPDGMass() /
GeV 120 <<
" P= " << incidentParticle->Get4Momentum() /
GeV<<
" mass= " 121 << incidentParticle->Get4Momentum().m() /
GeV <<G4endl;
124 const G4ThreeVector
position = track.GetPosition();
125 std::vector<G4ParticleDefinition*> newParticles;
126 std::vector<G4ParticleDefinition*> newParticlesRHadron;
129 G4ParticleDefinition* targetParticle;
130 HadronicProcessHelper::ReactionProduct reactionProduct =
131 m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
134 for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
135 it != reactionProduct.end() ;
138 G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
139 newParticlesRHadron.push_back(productDefinition);
142 if( cProd!=0 && m_detachCloud)
143 productDefinition=cProd->
GetCloud();
145 newParticles.push_back(productDefinition);
148 int numberOfSecondaries = reactionProduct.size();
150 if(m_verboseLevel >= 2) {
151 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt N secondaries: " 152 << numberOfSecondaries << G4endl;
156 G4LorentzVector sum4Momentum = incidentParticle->Get4Momentum();
157 G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
158 sum4Momentum += target4Momentum;
159 G4ThreeVector cmBoost = sum4Momentum.boostVector();
160 G4double M = sum4Momentum.m();
162 if(m_verboseLevel >= 2) {
163 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << G4endl
164 <<
" Boost = " << cmBoost /
GeV<< G4endl
165 <<
" 4P Lab = " << sum4Momentum /
GeV << G4endl
166 <<
" Ecm = " << M /
GeV << G4endl;
168 std::vector<G4double>
m;
171 for(
int ip=0; ip<numberOfSecondaries; ++ip)
173 m.push_back(newParticles[ip]->GetPDGMass());
176 std::vector<G4LorentzVector*>* fourMomenta = m_decay->Decay(M, m);
181 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
182 m_particleChange.ProposeTrackStatus(fStopAndKill);
183 if(m_verboseLevel >= 3) {
184 G4cout <<
"Incident does not survive: stopAndKill + set num secondaries to " 185 << numberOfSecondaries << G4endl;
188 for (
int ip=0; ip <numberOfSecondaries; ++ip) {
190 (*fourMomenta)[ip]->boost(cmBoost);
193 if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition()) {
195 if(m_verboseLevel >= 3) {
196 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " 197 << *((*fourMomenta)[ip])/
GeV <<
"(m="<< (*fourMomenta)[ip]->m()/
GeV 198 <<
") + " <<gluinoMomentum/
GeV<<G4endl;
200 G4LorentzVector p4_new = gluinoMomentum + *((*fourMomenta)[ip]);
201 G4ThreeVector momentum = p4_new.vect();
202 double rhMass = newParticlesRHadron[ip]->GetPDGMass();
203 (*fourMomenta)[ip]->setVectM(momentum,rhMass);
204 if(m_verboseLevel >= 3) {
205 G4cout <<
" = " << *((*fourMomenta)[ip])/
GeV <<
"(m="<< (*fourMomenta)[ip]->m() /
GeV 206 <<
") vs. "<<rhMass/
GeV 212 G4DynamicParticle* productDynParticle =
new G4DynamicParticle();
214 productDynParticle->SetDefinition(newParticlesRHadron[ip]);
216 productDynParticle->Set4Momentum(*((*fourMomenta)[ip]));
218 if(m_verboseLevel >= 3) {
219 G4cout <<
"ToyModelHadronicProcess::PostStepDoIt: " << *((*fourMomenta)[ip])/
GeV 220 <<
"(m="<< (*fourMomenta)[ip]->m()/
GeV<<
") " 221 << newParticlesRHadron[ip]->GetParticleName()
225 G4Track* productTrack =
new G4Track(productDynParticle,
226 track.GetGlobalTime(),
228 productTrack->SetTouchableHandle(thisTouchable);
229 m_particleChange.AddSecondary(productTrack);
230 delete (*fourMomenta)[ip];
235 ClearNumberOfInteractionLengthLeft();
236 return &m_particleChange;
G4ParticleDefinition * GetCloud()
G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, G4double aTemp)
G4ParticleDefinition * GetSpectator()
static int position[264][3]