CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FullModelHadronicProcess.cc
Go to the documentation of this file.
1 #include "G4HadReentrentException.hh"
2 #include "G4ProcessManager.hh"
3 #include "G4ParticleTable.hh"
4 
5 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
6 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
10 
11 using namespace CLHEP;
12 
13 FullModelHadronicProcess::FullModelHadronicProcess(G4ProcessHelper * aHelper,
14  const G4String& processName) :
15  G4VDiscreteProcess(processName), theHelper(aHelper)
16 {}
17 
18 FullModelHadronicProcess::~FullModelHadronicProcess()
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=0;
79  G4ParticleDefinition* outgoingCloud=0;
80  G4ParticleDefinition* outgoingTarget=0;
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() == 0)
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  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!=0) {
182  outgoingRhadron = tempDef;
183  //Setting outgoing cloud definition
184  outgoingCloud=tempCust->GetCloud();
185  if(outgoingCloud == 0) {
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==0&&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==0) 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 
233  /*
234  G4cout<<"Kinetic energies: "<<G4endl;
235  G4cout<<"True kinetic energy: "<<originalIncident->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
236  G4cout<<"Mass: "<<originalIncident->GetDefinition()->GetPDGMass()/GeV<<" GeV"<<G4endl;
237 
238  G4double e_kin_rescaled = targetNucleus.EvaporationEffects(originalIncident->GetTotalEnergy()-originalIncident->GetDefinition()->GetPDGMass());
239 
240  G4cout<<"Rescaled kinetic energy: "<<e_kin_rescaled<<G4endl;
241 
242  const G4double cutOff = 0.1*MeV;
243 
244  if ( e_kin_rescaled < cutOff )
245  {
246  aParticleChange.ProposeTrackStatus( fStopAndKill );//If the dice decides not to cascade I stop the particle
247  return &aParticleChange;
248  }
249  */
250  // create the target particle
251 
252  G4DynamicParticle *originalTarget = new G4DynamicParticle;
253  originalTarget->SetDefinition(aTarget);
254 
255  G4ReactionProduct targetParticle(aTarget);
256 
257  G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
258  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
259  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
260 
261  /*
262  G4cout<<"After creation:"<<G4endl;
263  G4cout<<"currentParticle: "<<currentParticle.GetMomentum()/GeV<<" GeV vs. "<<OrgPart->Get4Momentum()/GeV<<" GeV"<<G4endl;
264  G4cout<<"targetParticle: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
265  G4cout<<"Fourmomentum from originalIncident: "<<originalIncident->Get4Momentum()<<" vs "<<OrgPart->Get4Momentum()<<G4endl;
266  */
267 
268  G4ReactionProduct modifiedOriginal = currentParticle;
269 
270  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
271  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
272  G4bool incidentHasChanged = false;
273  if (!IncidentSurvives) incidentHasChanged = true; //I wonder if I am supposed to do this...
274  G4bool targetHasChanged = false;
275  if (!TargetSurvives) targetHasChanged = true; //Ditto here
276  G4bool quasiElastic = false;
277  if (rp.size()==2) quasiElastic = true; //Oh well...
278  G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> vec; // vec will contain the secondary particles
279  G4int vecLen = 0;
280  vec.Initialize( 0 );
281 
282  // I hope my understanding of "secondary" is correct here
283  // I think that it entails only what is not a surviving incident of target
284 
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)
290  {
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 );
295  }
296  }
297 
298  if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
299  if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );//Is this correct? It solves the "free energy" checking in ReactionDynamics
300  if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
301 
302  // G4cout<<"Calling CalculateMomenta... "<<G4endl;
303  /*
304  G4cout<<"Current particle starts as: "<<currentParticle.GetDefinition()->GetParticleName()<<G4endl;
305  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
306  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
307  */
308 
309  CalculateMomenta( vec, vecLen,
310  originalIncident, originalTarget, modifiedOriginal,
311  targetNucleus, currentParticle, targetParticle,
312  incidentHasChanged, targetHasChanged, quasiElastic );
313 
314  // G4cout <<"Cloud loss: "<<(e_temp-currentParticle.GetKineticEnergy())/GeV<<" GeV"<<G4endl;
315 
316  G4String cPname = currentParticle.GetDefinition()->GetParticleName();
317 
318  // if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud")
319  // {
320  /*
321  G4cout<<"Current particle is now: "<<cPname <<G4endl;
322  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
323  G4cout<<"and kinetic energy: "<<currentParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
324  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
325  G4cout<<"Modified original is: "<<modifiedOriginal.GetDefinition()->GetParticleName()<<G4endl;
326  G4cout<<"with momentum: "<<modifiedOriginal.GetMomentum()/GeV<<" GeV"<<G4endl;
327  G4cout<<"and kinetic energy: "<<modifiedOriginal.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
328  G4cout<<"May be killed?: "<<modifiedOriginal.GetMayBeKilled()<<G4endl;
329  G4cout<<"Target particle is: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
330  G4cout<<"with momentum: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
331  G4cout<<"and kinetic energy: "<<targetParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
332  G4cout<<"May be killed?: "<<targetParticle.GetMayBeKilled()<<G4endl;
333  G4cout<<"incidentHasChanged: "<<incidentHasChanged<<G4endl;
334  G4cout<<"targetHasChanged: "<<targetHasChanged<<G4endl;
335  G4cout<<"Particles in vec:"<<G4endl;
336  for(int i=0; i<vecLen; ++i )
337  {
338  G4cout<< vec[i]->GetDefinition()->GetParticleName()<<G4endl;
339  }
340  */
341  // G4cout<<"Done!"<<G4endl;
342  // const G4LorentzRotation& trans(originalIncident->GetTrafoToLab());
343 
344  aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
345  G4double e_kin=0;
346  G4LorentzVector cloud_p4_new; //Cloud 4-momentum in lab after collision
347  // n++;
348  // G4cout << n << G4endl;
349  /*
350  if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud") {
351  G4cout<<"Cloud deleted!!! AAARRRRGGGHHH!!!"<<G4endl;
352  G4cout<<"Cloud name: "<<cPname<<G4endl;
353  G4cout<<"E_kin_0: "<<e_kin_0/GeV<<" GeV"<<G4endl;
354  // G4cout<<"n: "<<n<<G4endl;
355  // n=0;
356  }
357  */
358  cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
359  cloud_p4_new *= trans;
360 
361  G4LorentzVector cloud_p4_old_full = Cloud4Momentum; //quark system in CMS BEFORE collision
362  cloud_p4_old_full.boost(trafo_full_cms);
363  G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum; //quark system in cloud CMS BEFORE collision
364  cloud_p4_old_cloud.boost(trafo_cloud_cms);
365  G4LorentzVector cloud_p4_full = cloud_p4_new; //quark system in CMS AFTER collision
366  cloud_p4_full.boost(trafo_full_cms);
367  G4LorentzVector cloud_p4_cloud = cloud_p4_new; //quark system in cloud CMS AFTER collision
368  cloud_p4_cloud.boost(trafo_cloud_cms);
369 
370  G4LorentzVector p_g_cms = gluinoMomentum; //gluino in CMS BEFORE collision
371  p_g_cms.boost(trafo_full_cms);
372 
373  G4LorentzVector p4_new( cloud_p4_new.v() + gluinoMomentum.v(), outgoingRhadron->GetPDGMass() );
374  // G4cout<<"P4-diff: "<<(p4_new-cloud_p4_new-gluinoMomentum)/GeV<<", magnitude: "
375  // <<(p4_new-cloud_p4_new-gluinoMomentum).m()/MeV<<" MeV" <<G4endl;
376 
377  G4ThreeVector p_new = p4_new.vect();
378 
379  aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).m());
380 
381  if( incidentHasChanged )
382  {
383  G4DynamicParticle* p0 = new G4DynamicParticle;
384  p0->SetDefinition( outgoingRhadron );
385  p0->SetMomentum( p_new );
386 
387  // May need to run SetDefinitionAndUpdateE here...
388  G4Track* Track0 = new G4Track(p0,
389  aTrack.GetGlobalTime(),
390  aPosition);
391  Track0->SetTouchableHandle(thisTouchable);
392  aParticleChange.AddSecondary(Track0);
393  /*
394  G4cout<<"Adding a particle "<<p0->GetDefinition()->GetParticleName()<<G4endl;
395  G4cout<<"with momentum: "<<p0->GetMomentum()/GeV<<" GeV"<<G4endl;
396  G4cout<<"and kinetic energy: "<<p0->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
397  */
398  if(p0->GetKineticEnergy()>e_kin_0) {
399  G4cout<<"ALAAAAARM!!! (incident changed from "
400  <<IncidentRhadron->GetDefinition()->GetParticleName()
401  <<" to "<<p0->GetDefinition()->GetParticleName()<<")"<<G4endl;
402  G4cout<<"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/GeV
403  <<" GeV (should be positive)"<<G4endl;
404  //Turns out problem is only in 2 -> 3 (Won't fix 2 -> 2 energy deposition)
405  if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
406  } /*else {
407  G4cout<<"NO ALAAAAARM!!!"<<G4endl;
408  }*/
409  if(std::abs(p0->GetKineticEnergy()-e_kin_0)>100*GeV) {
410  G4cout<<"Diff. too big"<<G4endl;
411  }
412  aParticleChange.ProposeTrackStatus( fStopAndKill );
413  }
414  else
415  {
416  G4double p = p_new.mag();
417  if( p > DBL_MIN )
418  aParticleChange.ProposeMomentumDirection( p_new.x()/p, p_new.y()/p, p_new.z()/p );
419  else
420  aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
421 
422  G4double aE = sqrt(p*p+(outgoingRhadron->GetPDGMass()*outgoingRhadron->GetPDGMass()) );
423  e_kin = aE - outgoingRhadron->GetPDGMass();
424  }
425 
426  // return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
427  if( targetParticle.GetMass() > 0.0 ) // targetParticle can be eliminated in TwoBody
428  {
429  G4DynamicParticle *p1 = new G4DynamicParticle;
430  p1->SetDefinition( targetParticle.GetDefinition() );
431  //G4cout<<"Target secondary: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
432  G4ThreeVector momentum = targetParticle.GetMomentum();
433  momentum = momentum.rotate(cache,what);
434  p1->SetMomentum( momentum );
435  p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
436  G4Track* Track1 = new G4Track(p1,
437  aTrack.GetGlobalTime(),
438  aPosition);
439  Track1->SetTouchableHandle(thisTouchable);
440  aParticleChange.AddSecondary(Track1);
441  }
442  G4DynamicParticle *pa;
443  /*
444  G4cout<<"vecLen: "<<vecLen<<G4endl;
445  G4cout<<"#sec's: "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
446  */
447 
448  for(int i=0; i<vecLen; ++i )
449  {
450  pa = new G4DynamicParticle();
451  pa->SetDefinition( vec[i]->GetDefinition() );
452  pa->SetMomentum( vec[i]->GetMomentum() );
453  pa->Set4Momentum(trans*(pa->Get4Momentum()));
454  G4ThreeVector pvec = pa->GetMomentum();
455  G4Track* Trackn = new G4Track(pa,
456  aTrack.GetGlobalTime(),
457  aPosition);
458  Trackn->SetTouchableHandle(thisTouchable);
459  aParticleChange.AddSecondary(Trackn);
460 
461  delete vec[i];
462  }
463 
464  // Histogram filling
465  const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
466 
467  if (theRhadron!=NULL||IncidentSurvives) {
468  double E_new;
469  if(IncidentSurvives) {
470  E_new = e_kin;
471  } else {
472  E_new = theRhadron->GetKineticEnergy();
473  if(CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
474  !=CustomPDGParser::s_isRMeson(theIncidentPDG)
475  ||
476  CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
477  !=CustomPDGParser::s_isMesonino(theIncidentPDG)
478  ) {
479  G4cout<<"Rm: "<<CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
480  <<" vs: "<<CustomPDGParser::s_isRMeson(theIncidentPDG)<<G4endl;
481  G4cout<<"Sm: "<<CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
482  <<" vs: "<<CustomPDGParser::s_isMesonino(theIncidentPDG)<<G4endl;
483 
484  }
485  }
486 
487  //Calculating relevant scattering angles.
488  G4LorentzVector p4_old_full = FullRhadron4Momentum; //R-hadron in CMS BEFORE collision
489  p4_old_full.boost(trafo_full_cms);
490  G4LorentzVector p4_old_cloud = FullRhadron4Momentum; //R-hadron in cloud CMS BEFORE collision
491  p4_old_cloud.boost(trafo_cloud_cms);
492  G4LorentzVector p4_full = p4_new; //R-hadron in CMS AFTER collision
493  // G4cout<<p4_full.v()/GeV<<G4endl;
494  p4_full=p4_full.boost(trafo_full_cms);
495  // G4cout<<p4_full.m()<<" / "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).m()<<G4endl;
496  G4LorentzVector p4_cloud = p4_new; //R-hadron in cloud CMS AFTER collision
497  p4_cloud.boost(trafo_cloud_cms);
498 
499  G4double AbsDeltaE = E_0-E_new;
500  // G4cout <<"Energy loss: "<<AbsDeltaE/GeV<<G4endl;
501  if (AbsDeltaE > 10*GeV) {
502  G4cout<<"Energy loss larger than 10 GeV..."<<G4endl;
503  G4cout<<"E_0: "<<E_0/GeV<<" GeV"<<G4endl;
504  G4cout<<"E_new: "<<E_new/GeV<<" GeV"<<G4endl;
505  G4cout<<"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
506  G4cout<<"x: "<<aPosition.x()/cm<<" cm"<<G4endl;
507  }
508  }
509  delete originalIncident;
510  delete originalTarget;
511  // aParticleChange.DumpInfo();
512  // G4cout << "Exiting FullModelHadronicProcess::PostStepDoIt"<<G4endl;
513 
514  //clear interaction length
515  ClearNumberOfInteractionLengthLeft();
516 
517  return &aParticleChange;
518 }
519 
520 
521 void FullModelHadronicProcess::CalculateMomenta(
522  G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
523  G4int &vecLen,
524  const G4HadProjectile *originalIncident, // the original incident particle
525  const G4DynamicParticle *originalTarget,
526  G4ReactionProduct &modifiedOriginal, // Fermi motion and evap. effects included
527  G4Nucleus &targetNucleus,
528  G4ReactionProduct &currentParticle,
529  G4ReactionProduct &targetParticle,
530  G4bool &incidentHasChanged,
531  G4bool &targetHasChanged,
532  G4bool quasiElastic )
533 {
534  FullModelReactionDynamics theReactionDynamics;
535 
536  cache = 0;
537  what = originalIncident->Get4Momentum().vect();
538 
539  if( quasiElastic )
540  {
541  // G4cout<<"We are calling TwoBody..."<<G4endl;
542  theReactionDynamics.TwoBody( vec, vecLen,
543  modifiedOriginal, originalTarget,
544  currentParticle, targetParticle,
545  targetNucleus, targetHasChanged );
546 
547  return;
548  }
549 
550  //If ProduceStrangeParticlePairs is commented out, let's cut this one as well
551  G4ReactionProduct leadingStrangeParticle;
552  G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
553  targetParticle,
554  leadingStrangeParticle );
555 
556  //
557  // Note: the number of secondaries can be reduced in GenerateXandPt and TwoCluster
558  //
559  G4bool finishedGenXPt = false;
560  G4bool annihilation = false;
561  if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
562  currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
563  {
564  // original was an anti-particle and annihilation has taken place
565  annihilation = true;
566  G4double ekcor = 1.0;
567  G4double ek = originalIncident->GetKineticEnergy();
568  G4double ekOrg = ek;
569 
570  const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
571  if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
572  const G4double atomicWeight = G4double(targetNucleus.GetN_asInt());
573  ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
574 
575  G4double tkin = targetNucleus.Cinema( ek );
576  ek += tkin;
577  ekOrg += tkin;
578  modifiedOriginal.SetKineticEnergy( ekOrg );
579  }
580 
581  const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
582  G4double rand1 = G4UniformRand();
583  G4double rand2 = G4UniformRand();
584  if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/GeV >= 1.0)) &&
585  (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
586  (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
587  (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
588  (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
589  ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
590  finishedGenXPt =
591  theReactionDynamics.GenerateXandPt( vec, vecLen,
592  modifiedOriginal, originalIncident,
593  currentParticle, targetParticle,
594  targetNucleus, incidentHasChanged,
595  targetHasChanged, leadFlag,
596  leadingStrangeParticle );
597  if( finishedGenXPt )
598  {
599  Rotate(vec, vecLen);
600  return;
601  }
602 
603  G4bool finishedTwoClu = false;
604  if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
605  {
606 
607  for(G4int i=0; i<vecLen; i++) delete vec[i];
608  vecLen = 0;
609  }
610  else
611  {
612 
613  theReactionDynamics.SuppressChargedPions( vec, vecLen,
614  modifiedOriginal, currentParticle,
615  targetParticle, targetNucleus,
616  incidentHasChanged, targetHasChanged );
617 
618  try
619  {
620  finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
621  modifiedOriginal, originalIncident,
622  currentParticle, targetParticle,
623  targetNucleus, incidentHasChanged,
624  targetHasChanged, leadFlag,
625  leadingStrangeParticle );
626  }
627  catch(G4HadReentrentException aC)
628  {
629  aC.Report(G4cout);
630  throw G4HadReentrentException(__FILE__, __LINE__, "Failing to calculate momenta");
631  }
632  }
633  if( finishedTwoClu )
634  {
635  Rotate(vec, vecLen);
636  return;
637  }
638 
639  //
640  // PNBlackTrackEnergy is the kinetic energy available for
641  // proton/neutron black track particles [was enp(1) in fortran code]
642  // DTABlackTrackEnergy is the kinetic energy available for
643  // deuteron/triton/alpha particles [was enp(3) in fortran code]
644  // the atomic weight of the target nucleus is >= 1.5 AND
645  // neither the incident nor the target particles have changed AND
646  // there is no kinetic energy available for either proton/neutron
647  // or for deuteron/triton/alpha black track particles
648  // For diffraction scattering on heavy nuclei use elastic routines instead
649 
650  theReactionDynamics.TwoBody( vec, vecLen,
651  modifiedOriginal, originalTarget,
652  currentParticle, targetParticle,
653  targetNucleus, targetHasChanged );
654 }
655 
656 G4bool FullModelHadronicProcess::MarkLeadingStrangeParticle(
657  const G4ReactionProduct &currentParticle,
658  const G4ReactionProduct &targetParticle,
659  G4ReactionProduct &leadParticle )
660 {
661  // the following was in GenerateXandPt and TwoCluster
662  // add a parameter to the GenerateXandPt function telling it about the strange particle
663  //
664  // assumes that the original particle was a strange particle
665  //
666  G4bool lead = false;
667  if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
668  (currentParticle.GetDefinition() != G4Proton::Proton()) &&
669  (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
670  {
671  lead = true;
672  leadParticle = currentParticle; // set lead to the incident particle
673  }
674  else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
675  (targetParticle.GetDefinition() != G4Proton::Proton()) &&
676  (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
677  {
678  lead = true;
679  leadParticle = targetParticle; // set lead to the target particle
680  }
681  return lead;
682 }
683 
684 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,MYGHADLISTSIZE> &vec,
685  G4int &vecLen)
686 {
687  G4double rotation = 2.*pi*G4UniformRand();
688  cache = rotation;
689  G4int i;
690  for( i=0; i<vecLen; ++i )
691  {
692  G4ThreeVector momentum = vec[i]->GetMomentum();
693  momentum = momentum.rotate(rotation, what);
694  vec[i]->SetMomentum(momentum);
695  }
696 }
697 
698 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
699 {
700  G4int nsec = aParticleChange->GetNumberOfSecondaries();
701  if (nsec==0) return 0;
702  int i = 0;
703  G4bool found = false;
704  while (i!=nsec && !found){
705  // G4cout<<"Checking "<<aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleName()<<G4endl;
706  // if (aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleType()=="rhadron") found = true;
707  if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found = true;
708  i++;
709  }
710  i--;
711  if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
712  return 0;
713 }
714 
int i
Definition: DBlmapReader.cc:9
const double GeV
Definition: MathUtil.h:16
G4ParticleDefinition * GetCloud()
#define NULL
Definition: scimark2.h:8
susybsm::HSCParticle pa
Definition: classes.h:8
static bool s_isRMeson(int pdg)
const Double_t pi
const double MeV
static bool s_isMesonino(int pdg)
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
G4ParticleDefinition * GetSpectator()
double p1[4]
Definition: TauolaWrapper.h:89
dbl *** dir
Definition: mlp_gen.cc:35