2 #include "G4ProcessManager.hh"
3 #include "G4ParticleTable.hh"
5 #include "G4InelasticInteraction.hh"
6 #include "Randomize.hh"
8 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
11 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
14 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper,
const G4String& processName) :
15 G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(
true)
25 G4bool ToyModelHadronicProcess::IsApplicable(
const G4ParticleDefinition& aP)
27 return m_helper->applicabilityTester(aP);
30 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(
const G4DynamicParticle *particle,
31 const G4Element *element,
35 G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
38 G4double highestEnergyLimit = 10 * TeV ;
39 G4double lowestEnergyLimit = 1 * eV;
40 G4double particleEnergy = particle->GetKineticEnergy();
42 if (particleEnergy > highestEnergyLimit || particleEnergy < lowestEnergyLimit){
43 if(m_verboseLevel >= 1)
std::cout <<
"ToyModelHadronicProcess: Energy out of bounds [" <<
44 lowestEnergyLimit / MeV <<
"MeV , " <<
45 highestEnergyLimit / MeV <<
"MeV ] while it is " << particleEnergy/MeV <<
49 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection << std::endl;
50 return inclusiveCrossSection;
56 G4double ToyModelHadronicProcess::GetMeanFreePath(
const G4Track& aTrack, G4double, G4ForceCondition*)
58 G4Material *aMaterial = aTrack.GetMaterial();
59 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
62 G4int nElements = aMaterial->GetNumberOfElements();
64 const G4double *theAtomicNumDensityVector =
65 aMaterial->GetAtomicNumDensityVector();
66 G4double aTemp = aMaterial->GetTemperature();
68 for( G4int
i=0;
i<nElements; ++
i )
71 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[
i], aTemp);
72 sigma += theAtomicNumDensityVector[
i] * xSection;
80 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(
const G4Track& track,
84 const G4TouchableHandle thisTouchable(track.GetTouchableHandle());
87 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
88 m_particleChange.Initialize(track);
90 const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle();
93 const G4ThreeVector aPosition = track.GetPosition();
98 G4DynamicParticle* cloudParticle =
new G4DynamicParticle();
103 if(CustomIncident==0)
105 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
107 cloudParticle->SetDefinition(CustomIncident->
GetCloud());
111 double scale=cloudParticle->GetDefinition()->GetPDGMass()/incidentRHadron->GetDefinition()->GetPDGMass();
112 G4LorentzVector cloudMomentum;
113 cloudMomentum.setVectM(incidentRHadron->GetMomentum()*
scale,cloudParticle->GetDefinition()->GetPDGMass());
114 G4LorentzVector gluinoMomentum;
115 gluinoMomentum.setVectM(incidentRHadron->GetMomentum()*(1.-
scale),CustomIncident->
GetSpectator()->GetPDGMass());
117 cloudParticle->Set4Momentum(cloudMomentum);
119 const G4DynamicParticle *incidentParticle;
120 if(! m_detachCloud) incidentParticle = incidentRHadron;
121 else incidentParticle= cloudParticle;
123 if(m_verboseLevel >= 3)
125 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt After scaling " << std::endl;
126 std::cout <<
" RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()<<
" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() / GeV
127 <<
" P= " << incidentRHadron->Get4Momentum() / GeV<<
" mass= "<< incidentRHadron->Get4Momentum().m() / GeV <<std::endl;
128 std::cout <<
" Cloud "<<cloudParticle->GetDefinition()->GetParticleName()<<
" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() / GeV
129 <<
" P= " << cloudParticle->Get4Momentum() / GeV<<
" mass= "<< cloudParticle->Get4Momentum().m() / GeV <<std::endl;
130 std::cout <<
" Incident pdgmass= " << incidentParticle->GetDefinition()->GetPDGMass() / GeV
131 <<
" P= " << incidentParticle->Get4Momentum() / GeV<<
" mass= "<< incidentParticle->Get4Momentum().m() / GeV <<std::endl;
134 const G4ThreeVector
position = track.GetPosition();
135 const G4int incidentParticlePDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
136 std::vector<G4ParticleDefinition*> newParticles;
137 std::vector<G4ParticleDefinition*> newParticlesRHadron;
139 G4bool incidentSurvives =
false;
142 G4ParticleDefinition* targetParticle;
143 HadronicProcessHelper::ReactionProduct reactionProduct = m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
148 for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
149 it != reactionProduct.end() ;
152 G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
154 if (productDefinition->GetPDGEncoding()==incidentParticlePDG)
155 incidentSurvives =
true;
157 newParticlesRHadron.push_back(productDefinition);
162 if( cProd!=0 && m_detachCloud)
163 productDefinition=cProd->
GetCloud();
165 newParticles.push_back(productDefinition);
168 int numberOfSecondaries = reactionProduct.size();
170 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt N secondaries: "
171 << numberOfSecondaries << std::endl;
179 const G4LorentzVector incident4Momentum = incidentParticle->Get4Momentum();
180 const G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
181 const G4LorentzVector sum4Momentum = incident4Momentum + target4Momentum;
182 const G4ThreeVector cmBoost = sum4Momentum.boostVector();
183 const G4LorentzVector cm4Momentum = sum4Momentum.rest4Vector();
185 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << std::endl
186 <<
" Boost = " << cmBoost / GeV<< std::endl
187 <<
" 4P CM = " << cm4Momentum / GeV << std::endl
188 <<
" 4P Inc = " << incident4Momentum / GeV << std::endl
189 <<
" 4P Targ = " << target4Momentum / GeV<< std::endl;
192 const G4double phi_p = 2*
pi*G4UniformRand()-
pi ;
193 const G4double theta_p = pi*G4UniformRand() ;
194 const G4ThreeVector randomDirection(
sin(theta_p)*
cos(phi_p),
199 std::vector<G4double>
m;
200 std::vector<G4LorentzVector> fourMomenta;
203 for(
int ip=0;ip<numberOfSecondaries;ip++)
205 m.push_back(newParticles[ip]->GetPDGMass());
209 if (numberOfSecondaries==2){
213 G4double
energy = cm4Momentum.e();
220 G4double cmMomentum = 1/(2*
energy)*
sqrt(energy*energy*energy*energy + m[0]*m[0]*m[0]*m[0] + m[1]*m[1]*m[1]*m[1]
221 - 2*(m[0]*m[0] + m[1]*m[1])*energy*energy -2*m[0]*m[0]*m[1]*m[1]);
222 p[0] = cmMomentum * randomDirection;
225 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt 2->2: " << std::endl
226 <<
" Pcm(GeV) = " << cmMomentum / GeV << std::endl;
228 for(
int ip=0;ip<2;ip++)
231 G4double
e =
sqrt(p[ip].
mag2() + m[ip]*m[ip]);
233 fourMomenta.push_back(G4LorentzVector(p[ip],e));
235 fourMomenta[ip].boost(cmBoost);
237 if(m_verboseLevel >= 2)
std::cout <<
" particle " << ip <<
" Plab(GeV) = " << fourMomenta[ip] /GeV << std::endl;
240 }
else if (numberOfSecondaries==3) {
244 for (std::vector<G4double>::iterator it=m.begin();it!=m.end();it++)
245 fourMomenta.push_back(G4LorentzVector(0,0,0,*it));
248 KinCalc.
doDecay(cm4Momentum, fourMomenta[0], fourMomenta[1], fourMomenta[2] );
251 CLHEP::HepRotation
rotation(randomDirection,G4UniformRand()*2*pi);
252 for (std::vector<G4LorentzVector>::iterator it = fourMomenta.begin();
253 it!=fourMomenta.end();
257 (*it).boost(cmBoost);
259 if(m_verboseLevel >= 3) G4cout<<
"Momentum-check: "<<incident4Momentum /GeV<<
" GeV vs "
260 << (fourMomenta[0]+fourMomenta[1]+fourMomenta[2])/GeV<<G4endl;
266 if(incidentSurvives){
268 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries-1);
269 if(m_verboseLevel >= 3)
std::cout <<
"Incident survives: set num secondaries to " << numberOfSecondaries-1 << std::endl;
273 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
274 m_particleChange.ProposeTrackStatus(fStopAndKill);
275 if(m_verboseLevel >= 3)
std::cout <<
"Incident does not survive: stopAndKill + set num secondaries to " << numberOfSecondaries << std::endl;
279 for (
int ip=0; ip <numberOfSecondaries;ip++)
281 if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition())
286 if(m_verboseLevel >= 3)
287 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
288 fourMomenta[ip]/GeV <<
"(m="<< fourMomenta[ip].m()/ GeV<<
") + " <<gluinoMomentum/GeV<<std::endl; ;
289 G4LorentzVector p4_new = gluinoMomentum+fourMomenta[ip];
290 G4ThreeVector momentum = p4_new.vect();
291 double rhMass=newParticlesRHadron[ip]->GetPDGMass() ;
294 fourMomenta[ip].setVectM(momentum,rhMass);
298 if(m_verboseLevel >= 3)
299 std::cout <<
" = " << fourMomenta[ip]/GeV <<
"(m="<< fourMomenta[ip].m() / GeV<<
") vs. "<<rhMass/GeV
302 m_particleChange.ProposeMomentumDirection(fourMomenta[ip].vect()/fourMomenta[ip].vect().
mag());
303 m_particleChange.ProposeEnergy(fourMomenta[ip].
e()-fourMomenta[ip].
mag());
304 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Propose momentum " << fourMomenta[ip]/GeV << std::endl;
307 G4DynamicParticle* productDynParticle =
new G4DynamicParticle();
309 productDynParticle->SetDefinition(newParticlesRHadron[ip]);
311 if(newParticlesRHadron[ip]!=newParticles[ip] && m_detachCloud )
313 G4LorentzVector p4_new;
314 p4_new.setVectM(gluinoMomentum.vect()+fourMomenta[ip].vect(),productDynParticle->GetDefinition()->GetPDGMass());
316 productDynParticle->Set4Momentum(p4_new);
320 if(m_verboseLevel >= 3)
321 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
323 <<
"(m="<< fourMomenta[ip].m()/ GeV<<
") + "
325 <<
" = " << productDynParticle->Get4Momentum()/GeV
326 <<
" => "<<productDynParticle->Get4Momentum().m()/GeV
327 <<
" vs. "<<productDynParticle->GetDefinition()->GetPDGMass()/GeV
332 productDynParticle->Set4Momentum(fourMomenta[ip]);
335 G4Track* productTrack =
new G4Track(productDynParticle,
336 track.GetGlobalTime(),
338 productTrack->SetTouchableHandle(thisTouchable);
340 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add secondary with 4-Momentum " << productDynParticle->Get4Momentum()/GeV << std::endl;
341 m_particleChange.AddSecondary(productTrack);
354 ClearNumberOfInteractionLengthLeft();
356 return &m_particleChange;
360 const G4DynamicParticle* ToyModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
362 G4int nsec = aParticleChange->GetNumberOfSecondaries();
363 if (nsec==0)
return 0;
365 G4bool
found =
false;
366 while (i!=nsec && !found){
367 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found =
true;
371 if(found)
return aParticleChange->GetSecondary(i)->GetDynamicParticle();
void doDecay(const G4LorentzVector &mother, G4LorentzVector &daughter1, G4LorentzVector &daughter2, G4LorentzVector &daughter3)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Sin< T >::type sin(const T &t)
G4ParticleDefinition * GetCloud()
static int position[TOTALCHAMBERS][3]
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
Cos< T >::type cos(const T &t)
G4ParticleDefinition * GetSpectator()