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 using namespace CLHEP;
16 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper,
const G4String& processName) :
17 G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(
true)
27 G4bool ToyModelHadronicProcess::IsApplicable(
const G4ParticleDefinition& aP)
29 return m_helper->applicabilityTester(aP);
32 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(
const G4DynamicParticle *particle,
33 const G4Element *element,
37 G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
40 G4double highestEnergyLimit = 10 * TeV ;
41 G4double lowestEnergyLimit = 1 * eV;
42 G4double particleEnergy = particle->GetKineticEnergy();
44 if (particleEnergy > highestEnergyLimit || particleEnergy < lowestEnergyLimit){
45 if(m_verboseLevel >= 1)
std::cout <<
"ToyModelHadronicProcess: Energy out of bounds [" <<
46 lowestEnergyLimit / MeV <<
"MeV , " <<
47 highestEnergyLimit / MeV <<
"MeV ] while it is " << particleEnergy/MeV <<
51 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection << std::endl;
52 return inclusiveCrossSection;
58 G4double ToyModelHadronicProcess::GetMeanFreePath(
const G4Track& aTrack, G4double, G4ForceCondition*)
60 G4Material *aMaterial = aTrack.GetMaterial();
61 const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
64 G4int nElements = aMaterial->GetNumberOfElements();
66 const G4double *theAtomicNumDensityVector =
67 aMaterial->GetAtomicNumDensityVector();
68 G4double aTemp = aMaterial->GetTemperature();
70 for( G4int
i=0;
i<nElements; ++
i )
73 GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[
i], aTemp);
74 sigma += theAtomicNumDensityVector[
i] * xSection;
82 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(
const G4Track& track,
86 const G4TouchableHandle thisTouchable(track.GetTouchableHandle());
89 G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
90 m_particleChange.Initialize(track);
92 const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle();
95 const G4ThreeVector aPosition = track.GetPosition();
100 G4DynamicParticle* cloudParticle =
new G4DynamicParticle();
105 if(CustomIncident==0)
107 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
109 cloudParticle->SetDefinition(CustomIncident->
GetCloud());
113 double scale=cloudParticle->GetDefinition()->GetPDGMass()/incidentRHadron->GetDefinition()->GetPDGMass();
114 G4LorentzVector cloudMomentum;
115 cloudMomentum.setVectM(incidentRHadron->GetMomentum()*
scale,cloudParticle->GetDefinition()->GetPDGMass());
116 G4LorentzVector gluinoMomentum;
117 gluinoMomentum.setVectM(incidentRHadron->GetMomentum()*(1.-
scale),CustomIncident->
GetSpectator()->GetPDGMass());
119 cloudParticle->Set4Momentum(cloudMomentum);
121 const G4DynamicParticle *incidentParticle;
122 if(! m_detachCloud) incidentParticle = incidentRHadron;
123 else incidentParticle= cloudParticle;
125 if(m_verboseLevel >= 3)
127 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt After scaling " << std::endl;
128 std::cout <<
" RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()<<
" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() / GeV
129 <<
" P= " << incidentRHadron->Get4Momentum() / GeV<<
" mass= "<< incidentRHadron->Get4Momentum().m() / GeV <<std::endl;
130 std::cout <<
" Cloud "<<cloudParticle->GetDefinition()->GetParticleName()<<
" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() / GeV
131 <<
" P= " << cloudParticle->Get4Momentum() / GeV<<
" mass= "<< cloudParticle->Get4Momentum().m() / GeV <<std::endl;
132 std::cout <<
" Incident pdgmass= " << incidentParticle->GetDefinition()->GetPDGMass() / GeV
133 <<
" P= " << incidentParticle->Get4Momentum() / GeV<<
" mass= "<< incidentParticle->Get4Momentum().m() / GeV <<std::endl;
136 const G4ThreeVector
position = track.GetPosition();
137 const G4int incidentParticlePDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
138 std::vector<G4ParticleDefinition*> newParticles;
139 std::vector<G4ParticleDefinition*> newParticlesRHadron;
141 G4bool incidentSurvives =
false;
144 G4ParticleDefinition* targetParticle;
145 HadronicProcessHelper::ReactionProduct reactionProduct = m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
150 for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
151 it != reactionProduct.end() ;
154 G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
156 if (productDefinition->GetPDGEncoding()==incidentParticlePDG)
157 incidentSurvives =
true;
159 newParticlesRHadron.push_back(productDefinition);
164 if( cProd!=0 && m_detachCloud)
165 productDefinition=cProd->
GetCloud();
167 newParticles.push_back(productDefinition);
170 int numberOfSecondaries = reactionProduct.size();
172 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt N secondaries: "
173 << numberOfSecondaries << std::endl;
181 const G4LorentzVector incident4Momentum = incidentParticle->Get4Momentum();
182 const G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
183 const G4LorentzVector sum4Momentum = incident4Momentum + target4Momentum;
184 const G4ThreeVector cmBoost = sum4Momentum.boostVector();
185 const G4LorentzVector cm4Momentum = sum4Momentum.rest4Vector();
187 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << std::endl
188 <<
" Boost = " << cmBoost / GeV<< std::endl
189 <<
" 4P CM = " << cm4Momentum / GeV << std::endl
190 <<
" 4P Inc = " << incident4Momentum / GeV << std::endl
191 <<
" 4P Targ = " << target4Momentum / GeV<< std::endl;
194 const G4double phi_p = 2*
pi*G4UniformRand()-
pi ;
195 const G4double theta_p = pi*G4UniformRand() ;
196 const G4ThreeVector randomDirection(
sin(theta_p)*
cos(phi_p),
201 std::vector<G4double>
m;
202 std::vector<G4LorentzVector> fourMomenta;
205 for(
int ip=0;ip<numberOfSecondaries;ip++)
207 m.push_back(newParticles[ip]->GetPDGMass());
211 if (numberOfSecondaries==2){
215 G4double
energy = cm4Momentum.e();
222 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]
223 - 2*(m[0]*m[0] + m[1]*m[1])*energy*energy -2*m[0]*m[0]*m[1]*m[1]);
224 p[0] = cmMomentum * randomDirection;
227 if(m_verboseLevel >= 2)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt 2->2: " << std::endl
228 <<
" Pcm(GeV) = " << cmMomentum / GeV << std::endl;
230 for(
int ip=0;ip<2;ip++)
233 G4double
e =
sqrt(p[ip].
mag2() + m[ip]*m[ip]);
235 fourMomenta.push_back(G4LorentzVector(p[ip],e));
237 fourMomenta[ip].boost(cmBoost);
239 if(m_verboseLevel >= 2)
std::cout <<
" particle " << ip <<
" Plab(GeV) = " << fourMomenta[ip] /GeV << std::endl;
242 }
else if (numberOfSecondaries==3) {
246 for (std::vector<G4double>::iterator it=m.begin();it!=m.end();it++)
247 fourMomenta.push_back(G4LorentzVector(0,0,0,*it));
250 KinCalc.
doDecay(cm4Momentum, fourMomenta[0], fourMomenta[1], fourMomenta[2] );
253 CLHEP::HepRotation
rotation(randomDirection,G4UniformRand()*2*pi);
254 for (std::vector<G4LorentzVector>::iterator it = fourMomenta.begin();
255 it!=fourMomenta.end();
259 (*it).boost(cmBoost);
261 if(m_verboseLevel >= 3) G4cout<<
"Momentum-check: "<<incident4Momentum /GeV<<
" GeV vs "
262 << (fourMomenta[0]+fourMomenta[1]+fourMomenta[2])/GeV<<G4endl;
268 if(incidentSurvives){
270 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries-1);
271 if(m_verboseLevel >= 3)
std::cout <<
"Incident survives: set num secondaries to " << numberOfSecondaries-1 << std::endl;
275 m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
276 m_particleChange.ProposeTrackStatus(fStopAndKill);
277 if(m_verboseLevel >= 3)
std::cout <<
"Incident does not survive: stopAndKill + set num secondaries to " << numberOfSecondaries << std::endl;
281 for (
int ip=0; ip <numberOfSecondaries;ip++)
283 if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition())
288 if(m_verboseLevel >= 3)
289 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
290 fourMomenta[ip]/GeV <<
"(m="<< fourMomenta[ip].m()/ GeV<<
") + " <<gluinoMomentum/GeV<<std::endl; ;
291 G4LorentzVector p4_new = gluinoMomentum+fourMomenta[ip];
292 G4ThreeVector momentum = p4_new.vect();
293 double rhMass=newParticlesRHadron[ip]->GetPDGMass() ;
296 fourMomenta[ip].setVectM(momentum,rhMass);
300 if(m_verboseLevel >= 3)
301 std::cout <<
" = " << fourMomenta[ip]/GeV <<
"(m="<< fourMomenta[ip].m() / GeV<<
") vs. "<<rhMass/GeV
304 m_particleChange.ProposeMomentumDirection(fourMomenta[ip].vect()/fourMomenta[ip].vect().
mag());
305 m_particleChange.ProposeEnergy(fourMomenta[ip].
e()-fourMomenta[ip].
mag());
306 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Propose momentum " << fourMomenta[ip]/GeV << std::endl;
309 G4DynamicParticle* productDynParticle =
new G4DynamicParticle();
311 productDynParticle->SetDefinition(newParticlesRHadron[ip]);
313 if(newParticlesRHadron[ip]!=newParticles[ip] && m_detachCloud )
315 G4LorentzVector p4_new;
316 p4_new.setVectM(gluinoMomentum.vect()+fourMomenta[ip].vect(),productDynParticle->GetDefinition()->GetPDGMass());
318 productDynParticle->Set4Momentum(p4_new);
322 if(m_verboseLevel >= 3)
323 std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
325 <<
"(m="<< fourMomenta[ip].m()/ GeV<<
") + "
327 <<
" = " << productDynParticle->Get4Momentum()/GeV
328 <<
" => "<<productDynParticle->Get4Momentum().m()/GeV
329 <<
" vs. "<<productDynParticle->GetDefinition()->GetPDGMass()/GeV
334 productDynParticle->Set4Momentum(fourMomenta[ip]);
337 G4Track* productTrack =
new G4Track(productDynParticle,
338 track.GetGlobalTime(),
340 productTrack->SetTouchableHandle(thisTouchable);
342 if(m_verboseLevel >= 3)
std::cout <<
"ToyModelHadronicProcess::PostStepDoIt Add secondary with 4-Momentum " << productDynParticle->Get4Momentum()/GeV << std::endl;
343 m_particleChange.AddSecondary(productTrack);
356 ClearNumberOfInteractionLengthLeft();
358 return &m_particleChange;
362 const G4DynamicParticle* ToyModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
364 G4int nsec = aParticleChange->GetNumberOfSecondaries();
365 if (nsec==0)
return 0;
367 G4bool
found =
false;
368 while (i!=nsec && !found){
369 if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found =
true;
373 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()