CMS 3D CMS Logo

FullModelHadronicProcess.cc
Go to the documentation of this file.
1 #include "G4HadReentrentException.hh"
2 #include "G4ProcessManager.hh"
3 #include "G4ParticleTable.hh"
4 
10 
11 using namespace CLHEP;
12 
14  const G4String& processName) :
15  G4VDiscreteProcess(processName), theHelper(aHelper)
16 {}
17 
19 {}
20 
21 G4bool FullModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
22 {
23  return theHelper->ApplicabilityTester(aP);
24 }
25 
26 G4double FullModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
27  const G4Element *anElement,
28  G4double aTemp)
29 {
30  //Get the cross section for this particle/element combination from the ProcessHelper
31  G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
32  // G4cout<<"Returned cross section from helper was: "<<InclXsec/millibarn<<" millibarn"<<G4endl;
33  return InclXsec;
34 }
35 
36 G4double FullModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double,
37  G4ForceCondition*)
38 {
39  G4Material *aMaterial = aTrack.GetMaterial();
40  const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
41  G4double sigma = 0.0;
42 
43  G4int nElements = aMaterial->GetNumberOfElements();
44 
45  const G4double *theAtomicNumDensityVector =
46  aMaterial->GetAtomicNumDensityVector();
47  G4double aTemp = aMaterial->GetTemperature();
48 
49  for( G4int i=0; i<nElements; ++i )
50  {
51  G4double xSection =
52  GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
53  sigma += theAtomicNumDensityVector[i] * xSection;
54  }
55  G4double res = DBL_MAX;
56  if(sigma > 0.0) { res = 1./sigma; }
57  return res;
58 }
59 
60 G4VParticleChange* FullModelHadronicProcess::PostStepDoIt(const G4Track& aTrack,
61  const G4Step& aStep)
62 {
63  // G4cout<<"**** Entering FullModelHadronicProcess::PostStepDoIt ******"<<G4endl;
64  const G4TouchableHandle& thisTouchable(aTrack.GetTouchableHandle());
65 
66  // A little setting up
67  aParticleChange.Initialize(aTrack);
68  const G4DynamicParticle* IncidentRhadron = aTrack.GetDynamicParticle();
69  CustomParticle* CustomIncident = static_cast<CustomParticle*>(IncidentRhadron->GetDefinition());
70  const G4ThreeVector& aPosition = aTrack.GetPosition();
71  // G4cout<<"G: "<<aStep.GetStepLength()/cm<<G4endl;
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=nullptr;
79  G4ParticleDefinition* outgoingCloud=nullptr;
80  G4ParticleDefinition* outgoingTarget=nullptr;
81 
82  G4ThreeVector p_0 = IncidentRhadron->GetMomentum();
83  G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
84  // G4cout<<e_kin_0/GeV<<G4endl;
85 
86  G4DynamicParticle* cloudParticle = new G4DynamicParticle();
87  /*
88  if(CustomPDGParser::s_isRMeson(theIncidentPDG))
89  G4cout<<"Rmeson"<<G4endl;
90  if(CustomPDGParser::s_isRBaryon(theIncidentPDG))
91  G4cout<<"Rbaryon"<<G4endl;
92  */
93  cloudParticle->SetDefinition(CustomIncident->GetCloud());
94 
95  if(cloudParticle->GetDefinition() == nullptr)
96  {
97  G4cout << "FullModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << G4endl;
98  }
99  /*
100  G4cout<<"Incoming particle was "<<IncidentRhadron->GetDefinition()->GetParticleName()
101  <<". Corresponding cloud is "<<cloudParticle->GetDefinition()->GetParticleName()<<G4endl;
102  G4cout<<"Kinetic energy was: "<<IncidentRhadron->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
103  */
104  double scale=cloudParticle->GetDefinition()->GetPDGMass()
105  /IncidentRhadron->GetDefinition()->GetPDGMass();
106  // G4cout<<"Mass ratio: "<<scale<<G4endl;
107  G4LorentzVector cloudMomentum(IncidentRhadron->GetMomentum()*scale,
108  cloudParticle->GetDefinition()->GetPDGMass());
109  G4LorentzVector gluinoMomentum(IncidentRhadron->GetMomentum()*(1.-scale),
110  CustomIncident->GetSpectator()->GetPDGMass());
111 
112  //These two for getting CMS transforms later (histogramming purposes...)
113  G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
114  const G4LorentzVector& Cloud4Momentum = cloudMomentum;
115 
116  cloudParticle->Set4Momentum(cloudMomentum);
117 
118  G4DynamicParticle* OrgPart = cloudParticle;
119 
120  /*
121  G4cout<<"Original momentum: "<<IncidentRhadron->Get4Momentum().v().mag()/GeV
122  <<" GeV, corresponding to gamma: "
123  <<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
124 
125  G4cout<<"Cloud momentum: "<<cloudParticle->Get4Momentum().v().mag()/GeV
126  <<" GeV, corresponding to gamma: "
127  <<cloudParticle->GetTotalEnergy()/cloudParticle->GetDefinition()->GetPDGMass()<<G4endl;
128  */
129 
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 );
135  ek += tkin;
136 
137  // calculate black track energies
138  tkin = targetNucleus.EvaporationEffects( ek );
139  ek -= tkin;
140 
141  if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*MeV||ek<=0.) {
142  //Very rare event...
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 ); // AR_NEWCODE_IMPORT
147  return &aParticleChange;
148  }
149  OrgPart->SetKineticEnergy( ek );
150  G4double p = std::sqrt(ek*(ek + 2*amas));
151  OrgPart->SetMomentum(dir * p);
152 
153  //Get the final state particles
154  G4ParticleDefinition* aTarget;
155  ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
156  G4bool force2to2 = false;
157  // G4cout<<"Trying to get final state..."<<G4endl;
158  while(rp.size()!=2 && force2to2){
159  rp = theHelper->GetFinalState(aTrack,aTarget);
160  }
161  G4int NumberOfSecondaries = rp.size();
162  // G4cout<<"Multiplicity of selected final state: "<<rp.size()<<G4endl;
163 
164  //Getting CMS transforms. Boosting is done at histogram filling
165  G4LorentzVector Target4Momentum(0.,0.,0.,aTarget->GetPDGMass());
166 
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();
171 
172  // OK Let's make some particles :-)
173  // We're not using them for anything yet, I know, but let's make sure the machinery is there
174  for(ReactionProduct::iterator it = rp.begin();
175  it != rp.end(); ++it) {
176  G4ParticleDefinition* tempDef = theParticleTable->FindParticle(*it);
177  CustomParticle* tempCust = dynamic_cast<CustomParticle*>(tempDef);
178  if (tempDef==aTarget) TargetSurvives = true;
179 
180  // if (tempDef->GetParticleType()=="rhadron")
181  if (tempCust!=nullptr) {
182  outgoingRhadron = tempDef;
183  //Setting outgoing cloud definition
184  outgoingCloud=tempCust->GetCloud();
185  if(outgoingCloud == nullptr) {
186  G4cout << "FullModelHadronicProcess::PostStepDoIt Definition of outgoing particle cloud not available!" << G4endl;
187  }
188  /*
189  G4cout<<"Outgoing Rhadron is: "<<outgoingRhadron->GetParticleName()<<G4endl;
190  G4cout<<"Outgoing cloud is: "<<outgoingCloud->GetParticleName()<<G4endl;
191  */
192  }
193 
194  if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
195  if (tempCust==nullptr&&rp.size()==2) outgoingTarget = tempDef;
196  if (tempDef->GetPDGEncoding()==theIncidentPDG) {
197  IncidentSurvives = true;
198  } else {
199  theParticleDefinitions.push_back(tempDef);
200  }
201  }
202 
203  if (outgoingTarget==nullptr) outgoingTarget = theParticleTable->FindParticle(rp[1]);
204 
205  // A little debug information
206  /*
207  G4cout<<"The particles coming out of this reaction will be: ";
208  for (std::vector<G4DynamicParticle*>::iterator it = theDynamicParticles.begin();
209  it != theDynamicParticles.end();
210  it++){
211  G4cout<< (*it)->GetDefinition()->GetParticleName()<<" ";
212  }
213  G4cout<<G4endl;
214  */
215  // If the incident particle survives it is not a "secondary", although
216  // it would probably be easier to fStopAndKill the track every time.
217  if(IncidentSurvives) NumberOfSecondaries--;
218  aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
219 
220  // ADAPTED FROM G4LEPionMinusInelastic::ApplyYourself
221  // It is bloody ugly, but such is the way of cut 'n' paste
222 
223  // Set up the incident
224  // This is where rotation to z-axis is done
225  const G4HadProjectile* originalIncident = new G4HadProjectile(*OrgPart);
226 
227  //Maybe I need to calculate trafo from R-hadron... Bollocks! Labframe does not move. CMS does.
228  G4HadProjectile* org2 = new G4HadProjectile(*OrgPart);
229  G4LorentzRotation trans = org2->GetTrafoToLab();
230  delete org2;
231 
232  // create the target particle
233 
234  G4DynamicParticle *originalTarget = new G4DynamicParticle;
235  originalTarget->SetDefinition(aTarget);
236 
237  G4ReactionProduct targetParticle(aTarget);
238 
239  G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
240  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
241  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
242 
243  /*
244  G4cout<<"After creation:"<<G4endl;
245  G4cout<<"currentParticle: "<<currentParticle.GetMomentum()/GeV<<" GeV vs. "<<OrgPart->Get4Momentum()/GeV<<" GeV"<<G4endl;
246  G4cout<<"targetParticle: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
247  G4cout<<"Fourmomentum from originalIncident: "<<originalIncident->Get4Momentum()<<" vs "<<OrgPart->Get4Momentum()<<G4endl;
248  */
249 
250  G4ReactionProduct modifiedOriginal = currentParticle;
251 
252  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
253  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
254  G4bool incidentHasChanged = false;
255  if (!IncidentSurvives) incidentHasChanged = true; //I wonder if I am supposed to do this...
256  G4bool targetHasChanged = false;
257  if (!TargetSurvives) targetHasChanged = true; //Ditto here
258  G4bool quasiElastic = false;
259  if (rp.size()==2) quasiElastic = true; //Oh well...
260  G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> vec; // vec will contain the secondary particles
261  G4int vecLen = 0;
262  vec.Initialize( 0 );
263 
264  // I hope my understanding of "secondary" is correct here
265  // I think that it entails only what is not a surviving incident of target
266 
267  for (G4int i = 0; i!=NumberOfSecondaries;i++){
268  if(theParticleDefinitions[i]!=aTarget
269  && theParticleDefinitions[i]!=originalIncident->GetDefinition()
270  && theParticleDefinitions[i]!=outgoingRhadron
271  && theParticleDefinitions[i]!=outgoingTarget)
272  {
273  G4ReactionProduct* pa = new G4ReactionProduct;
274  pa->SetDefinition( theParticleDefinitions[i] );
275  (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
276  vec.SetElement( vecLen++, pa );
277  }
278  }
279 
280  if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
281  if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );//Is this correct? It solves the "free energy" checking in ReactionDynamics
282  if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
283 
284  // G4cout<<"Calling CalculateMomenta... "<<G4endl;
285  /*
286  G4cout<<"Current particle starts as: "<<currentParticle.GetDefinition()->GetParticleName()<<G4endl;
287  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
288  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
289  */
290 
291  CalculateMomenta( vec, vecLen,
292  originalIncident, originalTarget, modifiedOriginal,
293  targetNucleus, currentParticle, targetParticle,
294  incidentHasChanged, targetHasChanged, quasiElastic );
295 
296  // G4cout <<"Cloud loss: "<<(e_temp-currentParticle.GetKineticEnergy())/GeV<<" GeV"<<G4endl;
297 
298  G4String cPname = currentParticle.GetDefinition()->GetParticleName();
299 
300  // if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud")
301  // {
302  /*
303  G4cout<<"Current particle is now: "<<cPname <<G4endl;
304  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
305  G4cout<<"and kinetic energy: "<<currentParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
306  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
307  G4cout<<"Modified original is: "<<modifiedOriginal.GetDefinition()->GetParticleName()<<G4endl;
308  G4cout<<"with momentum: "<<modifiedOriginal.GetMomentum()/GeV<<" GeV"<<G4endl;
309  G4cout<<"and kinetic energy: "<<modifiedOriginal.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
310  G4cout<<"May be killed?: "<<modifiedOriginal.GetMayBeKilled()<<G4endl;
311  G4cout<<"Target particle is: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
312  G4cout<<"with momentum: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
313  G4cout<<"and kinetic energy: "<<targetParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
314  G4cout<<"May be killed?: "<<targetParticle.GetMayBeKilled()<<G4endl;
315  G4cout<<"incidentHasChanged: "<<incidentHasChanged<<G4endl;
316  G4cout<<"targetHasChanged: "<<targetHasChanged<<G4endl;
317  G4cout<<"Particles in vec:"<<G4endl;
318  for(int i=0; i<vecLen; ++i )
319  {
320  G4cout<< vec[i]->GetDefinition()->GetParticleName()<<G4endl;
321  }
322  */
323  // G4cout<<"Done!"<<G4endl;
324 
325  aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
326  G4double e_kin=0;
327  G4LorentzVector cloud_p4_new; //Cloud 4-momentum in lab after collision
328  // n++;
329  // G4cout << n << G4endl;
330  /*
331  if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud") {
332  G4cout<<"Cloud deleted!!! AAARRRRGGGHHH!!!"<<G4endl;
333  G4cout<<"Cloud name: "<<cPname<<G4endl;
334  G4cout<<"E_kin_0: "<<e_kin_0/GeV<<" GeV"<<G4endl;
335  // G4cout<<"n: "<<n<<G4endl;
336  // n=0;
337  }
338  */
339  cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
340  cloud_p4_new *= trans;
341 
342  G4LorentzVector cloud_p4_old_full = Cloud4Momentum; //quark system in CMS BEFORE collision
343  cloud_p4_old_full.boost(trafo_full_cms);
344  G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum; //quark system in cloud CMS BEFORE collision
345  cloud_p4_old_cloud.boost(trafo_cloud_cms);
346  G4LorentzVector cloud_p4_full = cloud_p4_new; //quark system in CMS AFTER collision
347  cloud_p4_full.boost(trafo_full_cms);
348  G4LorentzVector cloud_p4_cloud = cloud_p4_new; //quark system in cloud CMS AFTER collision
349  cloud_p4_cloud.boost(trafo_cloud_cms);
350 
351  G4LorentzVector p_g_cms = gluinoMomentum; //gluino in CMS BEFORE collision
352  p_g_cms.boost(trafo_full_cms);
353 
354  double e = cloud_p4_new.e() + gluinoMomentum.e();
355  if(outgoingRhadron) e += outgoingRhadron->GetPDGMass();
356  G4LorentzVector p4_new( cloud_p4_new.v() + gluinoMomentum.v(), e );
357  // G4cout<<"P4-diff: "<<(p4_new-cloud_p4_new-gluinoMomentum)/GeV<<", magnitude: "
358  // <<(p4_new-cloud_p4_new-gluinoMomentum).m()/MeV<<" MeV" <<G4endl;
359 
360  G4ThreeVector p_new = p4_new.vect();
361 
362  aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).m());
363 
364  if( incidentHasChanged )
365  {
366  G4DynamicParticle* p0 = new G4DynamicParticle;
367  p0->SetDefinition( outgoingRhadron );
368  p0->SetMomentum( p_new );
369 
370  // May need to run SetDefinitionAndUpdateE here...
371  G4Track* Track0 = new G4Track(p0,
372  aTrack.GetGlobalTime(),
373  aPosition);
374  Track0->SetTouchableHandle(thisTouchable);
375  aParticleChange.AddSecondary(Track0);
376  /*
377  G4cout<<"Adding a particle "<<p0->GetDefinition()->GetParticleName()<<G4endl;
378  G4cout<<"with momentum: "<<p0->GetMomentum()/GeV<<" GeV"<<G4endl;
379  G4cout<<"and kinetic energy: "<<p0->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
380  */
381  if(p0->GetKineticEnergy()>e_kin_0) {
382  G4cout<<"ALAAAAARM!!! (incident changed from "
383  <<IncidentRhadron->GetDefinition()->GetParticleName()
384  <<" to "<<p0->GetDefinition()->GetParticleName()<<")"<<G4endl;
385  G4cout<<"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/GeV
386  <<" GeV (should be positive)"<<G4endl;
387  //Turns out problem is only in 2 -> 3 (Won't fix 2 -> 2 energy deposition)
388  if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
389  } /*else {
390  G4cout<<"NO ALAAAAARM!!!"<<G4endl;
391  }*/
392  if(std::abs(p0->GetKineticEnergy()-e_kin_0)>100*GeV) {
393  G4cout<<"Diff. too big"<<G4endl;
394  }
395  aParticleChange.ProposeTrackStatus( fStopAndKill );
396  }
397  else
398  {
399  G4double p = p_new.mag();
400  if( p > DBL_MIN )
401  aParticleChange.ProposeMomentumDirection( p_new.x()/p, p_new.y()/p, p_new.z()/p );
402  else
403  aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
404  }
405 
406  // return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
407  if( targetParticle.GetMass() > 0.0 ) // targetParticle can be eliminated in TwoBody
408  {
409  G4DynamicParticle *p1 = new G4DynamicParticle;
410  p1->SetDefinition( targetParticle.GetDefinition() );
411  //G4cout<<"Target secondary: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
412  G4ThreeVector momentum = targetParticle.GetMomentum();
413  momentum = momentum.rotate(cache,what);
414  p1->SetMomentum( momentum );
415  p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
416  G4Track* Track1 = new G4Track(p1,
417  aTrack.GetGlobalTime(),
418  aPosition);
419  Track1->SetTouchableHandle(thisTouchable);
420  aParticleChange.AddSecondary(Track1);
421  }
422  G4DynamicParticle *pa;
423  /*
424  G4cout<<"vecLen: "<<vecLen<<G4endl;
425  G4cout<<"#sec's: "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
426  */
427 
428  for(int i=0; i<vecLen; ++i )
429  {
430  pa = new G4DynamicParticle();
431  pa->SetDefinition( vec[i]->GetDefinition() );
432  pa->SetMomentum( vec[i]->GetMomentum() );
433  pa->Set4Momentum(trans*(pa->Get4Momentum()));
434  G4ThreeVector pvec = pa->GetMomentum();
435  G4Track* Trackn = new G4Track(pa,
436  aTrack.GetGlobalTime(),
437  aPosition);
438  Trackn->SetTouchableHandle(thisTouchable);
439  aParticleChange.AddSecondary(Trackn);
440 
441  delete vec[i];
442  }
443 
444  // Histogram filling
445  const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
446 
447  if (theRhadron!=nullptr||IncidentSurvives) {
448  double E_new;
449  if(IncidentSurvives) {
450  E_new = e_kin;
451  } else {
452  E_new = theRhadron->GetKineticEnergy();
453  if(CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
454  !=CustomPDGParser::s_isRMeson(theIncidentPDG)
455  ||
456  CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
457  !=CustomPDGParser::s_isMesonino(theIncidentPDG)
458  ) {
459  G4cout<<"Rm: "<<CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
460  <<" vs: "<<CustomPDGParser::s_isRMeson(theIncidentPDG)<<G4endl;
461  G4cout<<"Sm: "<<CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
462  <<" vs: "<<CustomPDGParser::s_isMesonino(theIncidentPDG)<<G4endl;
463 
464  }
465  }
466 
467  //Calculating relevant scattering angles.
468  G4LorentzVector p4_old_full = FullRhadron4Momentum; //R-hadron in CMS BEFORE collision
469  p4_old_full.boost(trafo_full_cms);
470  G4LorentzVector p4_old_cloud = FullRhadron4Momentum; //R-hadron in cloud CMS BEFORE collision
471  p4_old_cloud.boost(trafo_cloud_cms);
472  G4LorentzVector p4_full = p4_new; //R-hadron in CMS AFTER collision
473  // G4cout<<p4_full.v()/GeV<<G4endl;
474  p4_full=p4_full.boost(trafo_full_cms);
475  // G4cout<<p4_full.m()<<" / "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).m()<<G4endl;
476  G4LorentzVector p4_cloud = p4_new; //R-hadron in cloud CMS AFTER collision
477  p4_cloud.boost(trafo_cloud_cms);
478 
479  G4double AbsDeltaE = E_0-E_new;
480  // G4cout <<"Energy loss: "<<AbsDeltaE/GeV<<G4endl;
481  if (AbsDeltaE > 10*GeV) {
482  G4cout<<"Energy loss larger than 10 GeV..."<<G4endl;
483  G4cout<<"E_0: "<<E_0/GeV<<" GeV"<<G4endl;
484  G4cout<<"E_new: "<<E_new/GeV<<" GeV"<<G4endl;
485  G4cout<<"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
486  G4cout<<"x: "<<aPosition.x()/cm<<" cm"<<G4endl;
487  }
488  }
489  delete originalIncident;
490  delete originalTarget;
491  // aParticleChange.DumpInfo();
492  // G4cout << "Exiting FullModelHadronicProcess::PostStepDoIt"<<G4endl;
493 
494  //clear interaction length
495  ClearNumberOfInteractionLengthLeft();
496 
497  return &aParticleChange;
498 }
499 
500 
502  G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
503  G4int &vecLen,
504  const G4HadProjectile *originalIncident, // the original incident particle
505  const G4DynamicParticle *originalTarget,
506  G4ReactionProduct &modifiedOriginal, // Fermi motion and evap. effects included
507  G4Nucleus &targetNucleus,
508  G4ReactionProduct &currentParticle,
509  G4ReactionProduct &targetParticle,
510  G4bool &incidentHasChanged,
511  G4bool &targetHasChanged,
512  G4bool quasiElastic )
513 {
514  FullModelReactionDynamics theReactionDynamics;
515 
516  cache = 0;
517  what = originalIncident->Get4Momentum().vect();
518 
519  if( quasiElastic )
520  {
521  // G4cout<<"We are calling TwoBody..."<<G4endl;
522  theReactionDynamics.TwoBody( vec, vecLen,
523  modifiedOriginal, originalTarget,
524  currentParticle, targetParticle,
525  targetNucleus, targetHasChanged );
526 
527  return;
528  }
529 
530  //If ProduceStrangeParticlePairs is commented out, let's cut this one as well
531  G4ReactionProduct leadingStrangeParticle;
532  G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
533  targetParticle,
534  leadingStrangeParticle );
535 
536  //
537  // Note: the number of secondaries can be reduced in GenerateXandPt and TwoCluster
538  //
539  G4bool finishedGenXPt = false;
540  G4bool annihilation = false;
541  if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
542  currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
543  {
544  // original was an anti-particle and annihilation has taken place
545  annihilation = true;
546  G4double ekcor = 1.0;
547  G4double ek = originalIncident->GetKineticEnergy();
548  G4double ekOrg = ek;
549 
550  const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
551  if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
552  const G4double atomicWeight = G4double(targetNucleus.GetN_asInt());
553  ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
554 
555  G4double tkin = targetNucleus.Cinema( ek );
556  //ek += tkin;
557  ekOrg += tkin;
558  modifiedOriginal.SetKineticEnergy( ekOrg );
559  }
560 
561  const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
562  G4double rand1 = G4UniformRand();
563  G4double rand2 = G4UniformRand();
564  if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/GeV >= 1.0)) &&
565  (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
566  (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
567  (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
568  (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
569  ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
570  finishedGenXPt =
571  theReactionDynamics.GenerateXandPt( vec, vecLen,
572  modifiedOriginal, originalIncident,
573  currentParticle, targetParticle,
574  targetNucleus, incidentHasChanged,
575  targetHasChanged, leadFlag,
576  leadingStrangeParticle );
577  if( finishedGenXPt )
578  {
579  Rotate(vec, vecLen);
580  return;
581  }
582 
583  G4bool finishedTwoClu = false;
584  if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
585  {
586 
587  for(G4int i=0; i<vecLen; i++) delete vec[i];
588  vecLen = 0;
589  }
590  else
591  {
592 
593  theReactionDynamics.SuppressChargedPions( vec, vecLen,
594  modifiedOriginal, currentParticle,
595  targetParticle, targetNucleus,
596  incidentHasChanged, targetHasChanged );
597 
598  try
599  {
600  finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
601  modifiedOriginal, originalIncident,
602  currentParticle, targetParticle,
603  targetNucleus, incidentHasChanged,
604  targetHasChanged, leadFlag,
605  leadingStrangeParticle );
606  }
607  catch(G4HadReentrentException &aC)
608  {
609  aC.Report(G4cout);
610  throw G4HadReentrentException(__FILE__, __LINE__, "Failing to calculate momenta");
611  }
612  }
613  if( finishedTwoClu )
614  {
615  Rotate(vec, vecLen);
616  return;
617  }
618 
619  //
620  // PNBlackTrackEnergy is the kinetic energy available for
621  // proton/neutron black track particles [was enp(1) in fortran code]
622  // DTABlackTrackEnergy is the kinetic energy available for
623  // deuteron/triton/alpha particles [was enp(3) in fortran code]
624  // the atomic weight of the target nucleus is >= 1.5 AND
625  // neither the incident nor the target particles have changed AND
626  // there is no kinetic energy available for either proton/neutron
627  // or for deuteron/triton/alpha black track particles
628  // For diffraction scattering on heavy nuclei use elastic routines instead
629 
630  theReactionDynamics.TwoBody( vec, vecLen,
631  modifiedOriginal, originalTarget,
632  currentParticle, targetParticle,
633  targetNucleus, targetHasChanged );
634 }
635 
637  const G4ReactionProduct &currentParticle,
638  const G4ReactionProduct &targetParticle,
639  G4ReactionProduct &leadParticle )
640 {
641  // the following was in GenerateXandPt and TwoCluster
642  // add a parameter to the GenerateXandPt function telling it about the strange particle
643  //
644  // assumes that the original particle was a strange particle
645  //
646  G4bool lead = false;
647  if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
648  (currentParticle.GetDefinition() != G4Proton::Proton()) &&
649  (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
650  {
651  lead = true;
652  leadParticle = currentParticle; // set lead to the incident particle
653  }
654  else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
655  (targetParticle.GetDefinition() != G4Proton::Proton()) &&
656  (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
657  {
658  lead = true;
659  leadParticle = targetParticle; // set lead to the target particle
660  }
661  return lead;
662 }
663 
664 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
665  G4int &vecLen)
666 {
667  G4double rotation = 2.*pi*G4UniformRand();
668  cache = rotation;
669  G4int i;
670  for( i=0; i<vecLen; ++i )
671  {
672  G4ThreeVector momentum = vec[i]->GetMomentum();
673  momentum = momentum.rotate(rotation, what);
674  vec[i]->SetMomentum(momentum);
675  }
676 }
677 
678 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
679 {
680  G4int nsec = aParticleChange->GetNumberOfSecondaries();
681  if (nsec==0) return nullptr;
682  int i = 0;
683  G4bool found = false;
684  while (i!=nsec && !found){
685  // G4cout<<"Checking "<<aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleName()<<G4endl;
686  // if (aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleType()=="rhadron") found = true;
687  if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=nullptr) found = true;
688  i++;
689  }
690  i--;
691  if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
692  return nullptr;
693 }
694 
const double GeV
Definition: MathUtil.h:16
G4bool MarkLeadingStrangeParticle(const G4ReactionProduct &currentParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
G4bool ApplicabilityTester(const G4ParticleDefinition &aPart)
void TwoBody(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4DynamicParticle *originalTarget, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &targetHasChanged)
G4ParticleDefinition * GetCloud()
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
susybsm::HSCParticle pa
Definition: classes.h:8
static bool s_isRMeson(int pdg)
virtual G4double GetMicroscopicCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement, G4double aTemp)
Definition: Electron.h:6
std::vector< G4int > ReactionProduct
const Double_t pi
const double MeV
static bool s_isMesonino(int pdg)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ReactionProduct GetFinalState(const G4Track &aTrack, G4ParticleDefinition *&aTarget)
G4bool TwoCluster(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
G4double GetInclusiveCrossSection(const G4DynamicParticle *aParticle, const G4Element *anElement)
G4ParticleDefinition * GetSpectator()
G4bool GenerateXandPt(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
void CalculateMomenta(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
FullModelHadronicProcess(G4ProcessHelper *aHelper, const G4String &processName="FullModelHadronicProcess")
double p1[4]
Definition: TauolaWrapper.h:89
void SuppressChargedPions(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen, const G4ReactionProduct &modifiedOriginal, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged)
void Rotate(G4FastVector< G4ReactionProduct, MYGHADLISTSIZE > &vec, G4int &vecLen)
const G4DynamicParticle * FindRhadron(G4ParticleChange *)
dbl *** dir
Definition: mlp_gen.cc:35
G4bool IsApplicable(const G4ParticleDefinition &aP) override