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"
8 #include "SimG4Core/CustomPhysics/interface/FullModelReactionDynamics.hh"
11 
12 FullModelHadronicProcess::FullModelHadronicProcess(G4ProcessHelper * aHelper, const G4String& processName) :
13  G4VDiscreteProcess(processName), theHelper(aHelper)
14 {}
15 
16 
17 FullModelHadronicProcess::~FullModelHadronicProcess(){}
18 
19 G4bool FullModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
20 {
21  return theHelper->ApplicabilityTester(aP);
22 }
23 
24 G4double FullModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *aParticle,
25  const G4Element *anElement,
26  G4double aTemp)
27 {
28  //Get the cross section for this particle/element combination from the ProcessHelper
29  G4double InclXsec = theHelper->GetInclusiveCrossSection(aParticle,anElement);
30  // G4cout<<"Returned cross section from helper was: "<<InclXsec/millibarn<<" millibarn"<<G4endl;
31 
32  // Need to provide Set-methods for these in time
33  G4double HighestEnergyLimit = 10 * TeV ;
34  G4double LowestEnergyLimit = 1 * eV;
35 
36  G4double ParticleEnergy = aParticle->GetKineticEnergy();
37 
38 
39  if (ParticleEnergy > HighestEnergyLimit || ParticleEnergy < LowestEnergyLimit){
40  return 0;
41  } else {
42  // G4cout << "Microscopic Cross Section: "<<InclXsec / millibarn<<" millibarn"<<G4endl;
43  return InclXsec;
44  }
45 
46 }
47 
48 G4double FullModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*)
49 {
50  G4Material *aMaterial = aTrack.GetMaterial();
51  const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
52  G4double sigma = 0.0;
53 
54  G4int nElements = aMaterial->GetNumberOfElements();
55 
56  const G4double *theAtomicNumDensityVector =
57  aMaterial->GetAtomicNumDensityVector();
58  G4double aTemp = aMaterial->GetTemperature();
59 
60  for( G4int i=0; i<nElements; ++i )
61  {
62  G4double xSection =
63  GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
64  sigma += theAtomicNumDensityVector[i] * xSection;
65  }
66 
67  return 1.0/sigma;
68 
69 }
70 
71 G4VParticleChange* FullModelHadronicProcess::PostStepDoIt(const G4Track& aTrack,
72  const G4Step& aStep)
73 {
74  // G4cout<<"***************** Entering FullModelHadronicProcess::PostStepDoIt **********************"<<G4endl;
75 
76  const G4TouchableHandle thisTouchable(aTrack.GetTouchableHandle());
77 
78  // A little setting up
79  aParticleChange.Initialize(aTrack);
80  // G4DynamicParticle* OrgPart = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
81  G4DynamicParticle* IncidentRhadron = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
82  CustomParticle* CustomIncident = static_cast<CustomParticle*>(IncidentRhadron->GetDefinition());
83  const G4ThreeVector aPosition = aTrack.GetPosition();
84  // std::cout<<"G: "<<aStep.GetStepLength()/cm<<std::endl;
85  const G4int theIncidentPDG = IncidentRhadron->GetDefinition()->GetPDGEncoding();
86  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
87  std::vector<G4ParticleDefinition*> theParticleDefinitions;
88  // std::vector<G4DynamicParticle*> *theDynamicParticles;//These are probably redundant, but they can easily be removed :-)
89  G4bool IncidentSurvives = false;
90  G4bool TargetSurvives = false;
91  G4Nucleus targetNucleus(aTrack.GetMaterial());
92  G4ParticleDefinition* outgoingRhadron=0;
93  G4ParticleDefinition* outgoingCloud=0;
94  G4ParticleDefinition* outgoingTarget=0;
95 // double gamma = IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass();
96 
97 
98  G4ThreeVector p_0 = IncidentRhadron->GetMomentum();
99 
100  // static int n=0;
101 
102  G4double e_kin_0 = IncidentRhadron->GetKineticEnergy();
103  // G4cout<<e_kin_0/GeV<<G4endl;
104 
105  G4DynamicParticle* cloudParticle = new G4DynamicParticle();
106  /*
107  if(CustomPDGParser::s_isRMeson(theIncidentPDG))
108  std::cout<<"Rmeson"<<std::endl;
109  if(CustomPDGParser::s_isRBaryon(theIncidentPDG))
110  std::cout<<"Rbaryon"<<std::endl;
111  */
112  /*
113  if(CustomPDGParser::s_isRBaryon(theIncidentPDG))
114  cloudParticle->SetDefinition(theParticleTable->FindParticle("rhadronbaryoncloud"));
115  if(CustomPDGParser::s_isRMeson(theIncidentPDG) || CustomPDGParser::s_isRGlueball(theIncidentPDG) )
116  cloudParticle->SetDefinition(theParticleTable->FindParticle("rhadronmesoncloud"));
117  */
118  cloudParticle->SetDefinition(CustomIncident->GetCloud());
119 
120  if(cloudParticle->GetDefinition() == 0)
121  {
122  std::cout << "FullModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
123  }
124  /*
125  std::cout<<"Incoming particle was "<<IncidentRhadron->GetDefinition()->GetParticleName()<<". Corresponding cloud is "<<cloudParticle->GetDefinition()->GetParticleName()<<std::endl;
126  G4cout<<"Kinetic energy was: "<<IncidentRhadron->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
127  */
128  double scale=cloudParticle->GetDefinition()->GetPDGMass()/IncidentRhadron->GetDefinition()->GetPDGMass();
129  // std::cout<<"Mass ratio: "<<scale<<std::endl;
130  G4LorentzVector cloudMomentum;
131  cloudMomentum.setVectM(IncidentRhadron->GetMomentum()*scale,cloudParticle->GetDefinition()->GetPDGMass());
132  G4LorentzVector gluinoMomentum;
133  // gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-scale),theParticleTable->FindParticle("~g")->GetPDGMass());
134  gluinoMomentum.setVectM(IncidentRhadron->GetMomentum()*(1.-scale),CustomIncident->GetSpectator()->GetPDGMass());
135  /*
136  G4cout <<"Are these the same?"<<G4endl;
137  G4cout<<gluinoMomentum<<G4endl;
138  G4cout<<(1.-scale) * IncidentRhadron->Get4Momentum()<<G4endl;
139  */
140  //These two for getting CMS transforms later (histogramming purposes...)
141  G4LorentzVector FullRhadron4Momentum = IncidentRhadron->Get4Momentum();
142  G4LorentzVector Cloud4Momentum = cloudMomentum;
143 
144  cloudParticle->Set4Momentum(cloudMomentum);
145 
146  G4DynamicParticle* OrgPart = cloudParticle;
147 
148  /*
149  std::cout<<"Original momentum: "<<IncidentRhadron->Get4Momentum().v().mag()/GeV<<" GeV, corresponding to gamma: "
150  <<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<std::endl;
151 
152  std::cout<<"Cloud momentum: "<<cloudParticle->Get4Momentum().v().mag()/GeV<<" GeV, corresponding to gamma: "
153  <<cloudParticle->GetTotalEnergy()/cloudParticle->GetDefinition()->GetPDGMass()<<std::endl;
154  */
155 
156  double E_0 = IncidentRhadron->GetKineticEnergy();
157 
158 
159  G4double ek = OrgPart->GetKineticEnergy();
160  G4double amas = OrgPart->GetDefinition()->GetPDGMass();
161 
162  G4double tkin = targetNucleus.Cinema( ek );
163  ek += tkin;
164  OrgPart->SetKineticEnergy( ek );
165  G4double et = ek + amas;
166  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
167  G4double pp = OrgPart->GetMomentum().mag();
168  if( pp > 0.0 )
169  {
170  G4ThreeVector momentum = OrgPart->GetMomentum();
171  OrgPart->SetMomentum( momentum * (p/pp) );
172  }
173 
174  // calculate black track energies
175  if(ek > 0.) { tkin = targetNucleus.EvaporationEffects( ek ); ek -= tkin; } // AR_NEWCODE_IMPORT
176  if(ek+gluinoMomentum.e()-gluinoMomentum.m()<=0.1*MeV||ek<=0.) {
177  //Very rare event...
178  G4cout<<"Kinetic energy is sick"<<G4endl;
179  G4cout<<"Full R-hadron: "<<(ek+gluinoMomentum.e()-gluinoMomentum.m())/MeV<<" MeV" <<G4endl;
180  G4cout<<"Quark system: "<<ek/MeV<<" MeV"<<G4endl;
181 // aParticleChange.ProposeTrackStatus( fStopAndKill ); // AR_NEWCODE_IMPORT
182  aParticleChange.ProposeTrackStatus( fStopButAlive ); // AR_NEWCODE_IMPORT
183  return &aParticleChange;
184  }
185  OrgPart->SetKineticEnergy( ek );
186  et = ek + amas;
187  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
188  pp = OrgPart->GetMomentum().mag();
189 
190  if( pp > 0.0 )
191  {
192  G4ThreeVector momentum = OrgPart->GetMomentum();
193  OrgPart->SetMomentum( momentum * (p/pp) );
194  }
195 
196 
197 
198  //Get the final state particles
199 
200  G4ParticleDefinition* aTarget;
201  ReactionProduct rp = theHelper->GetFinalState(aTrack,aTarget);
202  G4bool force2to2 = false;
203  // G4cout<<"Trying to get final state..."<<G4endl;
204  while(rp.size()!=2 && force2to2){
205  rp = theHelper->GetFinalState(aTrack,aTarget);
206  }
207  G4int NumberOfSecondaries = rp.size();
208  // G4cout<<"Multiplicity of selected final state: "<<rp.size()<<G4endl;
209 
210  //Getting CMS transforms. Boosting is done at histogram filling
211  G4LorentzVector Target4Momentum;
212  Target4Momentum.setVectM(0.,aTarget->GetPDGMass());
213  // Target4Momentum.setVectM(0.,targetNucleus.GetN()*GeV);
214  G4LorentzVector psum_full,psum_cloud;
215  psum_full = FullRhadron4Momentum + Target4Momentum;
216  psum_cloud = Cloud4Momentum + Target4Momentum;
217  G4ThreeVector trafo_full_cms = (-1)*psum_full.boostVector();
218  G4ThreeVector trafo_cloud_cms = (-1)*psum_cloud.boostVector();
219  /*
220  psum_full.boost(trafo_full_cms);
221  psum_cloud.boost(trafo_cloud_cms);
222  std::cout<<"Checking that the momenta are in deed zero:"<<psum_full.vect()<<std::endl;
223  */
224 
225  // OK Let's make some particles :-)
226  // We're not using them for anything yet, I know, but let's make sure the machinery is there
227 
228  for(ReactionProduct::iterator it = rp.begin();
229  it != rp.end();
230  it++)
231  {
232  G4ParticleDefinition* tempDef = theParticleTable->FindParticle(*it);
233  CustomParticle* tempCust = dynamic_cast<CustomParticle*>(tempDef);
234  if (tempDef==aTarget) TargetSurvives = true;
235 
236  // if (tempDef->GetParticleType()=="rhadron")
237  if (tempCust!=0)
238  {
239  outgoingRhadron = tempDef;
240  //Setting outgoing cloud definition
241  /*
242  if(CustomPDGParser::s_isRBaryon(*it))
243  outgoingCloud=theParticleTable->FindParticle("rhadronbaryoncloud");
244  if(CustomPDGParser::s_isRMeson(*it) || CustomPDGParser::s_isRGlueball(*it) )
245  outgoingCloud=theParticleTable->FindParticle("rhadronmesoncloud");
246  */
247  outgoingCloud=tempCust->GetCloud();
248  if(outgoingCloud == 0)
249  {
250  std::cout << "FullModelHadronicProcess::PostStepDoIt Definition of outgoing particle cloud not available!!" << std::endl;
251  }
252  /*
253  std::cout<<"Outgoing Rhadron is: "<<outgoingRhadron->GetParticleName()<<std::endl;
254  std::cout<<"Outgoing cloud is: "<<outgoingCloud->GetParticleName()<<std::endl;
255  */
256  }
257 
258  if (tempDef==G4Proton::Proton()||tempDef==G4Neutron::Neutron()) outgoingTarget = tempDef;
259  // if (tempDef->GetParticleType()!="rhadron"&&rp.size()==2) outgoingTarget = tempDef;
260  if (tempCust==0&&rp.size()==2) outgoingTarget = tempDef;
261  if (tempDef->GetPDGEncoding()==theIncidentPDG) {
262  IncidentSurvives = true;
263  } else {
264  theParticleDefinitions.push_back(tempDef);
265  /*
266  G4DynamicParticle* tempDyn = new G4DynamicParticle();
267  tempDyn->SetDefinition(tempDef);
268  theDynamicParticles->push_back(tempDyn);
269  */
270  }
271  }
272 
273  //Not using this, so...
274  // delete theDynamicParticles;
275 
276  if (outgoingTarget==0) outgoingTarget = theParticleTable->FindParticle(rp[1]);
277 
278  // A little debug information
279  /*
280  G4cout<<"The particles coming out of this reaction will be: ";
281  for (std::vector<G4DynamicParticle*>::iterator it = theDynamicParticles.begin();
282  it != theDynamicParticles.end();
283  it++){
284  G4cout<< (*it)->GetDefinition()->GetParticleName()<<" ";
285  }
286  G4cout<<G4endl;
287  */
288  // If the incident particle survives it is not a "secondary", although
289  // it would probably be easier to fStopAndKill the track every time.
290  if(IncidentSurvives) NumberOfSecondaries--;
291  aParticleChange.SetNumberOfSecondaries(NumberOfSecondaries);
292 
293 
294  // ADAPTED FROM G4LEPionMinusInelastic::ApplyYourself
295  // It is bloody ugly, but such is the way of cut 'n' paste
296 
297 
298  // Set up the incident
299  const G4HadProjectile* originalIncident = new G4HadProjectile(*OrgPart);//This is where rotation to z-axis is done (Aarrggh!)
300 
301 
302  //Maybe I need to calculate trafo from R-hadron... Bollocks! Labframe does not move. CMS does.
303  G4HadProjectile* org2 = new G4HadProjectile(*OrgPart);
304  G4LorentzRotation trans = org2->GetTrafoToLab();
305  delete org2;
306 
307  // if (originalIncident->GetKineticEnergy()<= 0.1*MeV) { //Needs rescaling. The kinetic energy of the quarksystem is the relevant quantity
308 
309  /*
310  G4cout<<"Kinetic energies: "<<G4endl;
311  G4cout<<"True kinetic energy: "<<originalIncident->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
312  G4cout<<"Mass: "<<originalIncident->GetDefinition()->GetPDGMass()/GeV<<" GeV"<<G4endl;
313 
314  G4double e_kin_rescaled = targetNucleus.EvaporationEffects(originalIncident->GetTotalEnergy()-originalIncident->GetDefinition()->GetPDGMass());
315 
316  G4cout<<"Rescaled kinetic energy: "<<e_kin_rescaled<<G4endl;
317 
318  const G4double cutOff = 0.1*MeV;
319 
320  if ( e_kin_rescaled < cutOff )
321  {
322  aParticleChange.ProposeTrackStatus( fStopAndKill );//If the dice decides not to cascade I stop the particle
323  return &aParticleChange;
324  }
325  */
326  // create the target particle
327 
328  G4DynamicParticle *originalTarget = new G4DynamicParticle;
329  originalTarget->SetDefinition(aTarget);
330 
331  G4ReactionProduct targetParticle(aTarget);
332 
333 
334  G4ReactionProduct currentParticle(const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
335  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
336  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
337 
338  /*
339  G4cout<<"After creation:"<<G4endl;
340  G4cout<<"currentParticle: "<<currentParticle.GetMomentum()/GeV<<" GeV vs. "<<OrgPart->Get4Momentum()/GeV<<" GeV"<<G4endl;
341  G4cout<<"targetParticle: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
342  G4cout<<"Fourmomentum from originalIncident: "<<originalIncident->Get4Momentum()<<" vs "<<OrgPart->Get4Momentum()<<G4endl;
343  */
344 
345 
346  G4ReactionProduct modifiedOriginal = currentParticle;
347 
348  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
349  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
350  G4bool incidentHasChanged = false;
351  if (!IncidentSurvives) incidentHasChanged = true; //I wonder if I am supposed to do this...
352  G4bool targetHasChanged = false;
353  if (!TargetSurvives) targetHasChanged = true; //Ditto here
354  G4bool quasiElastic = false;
355  if (rp.size()==2) quasiElastic = true; //Oh well...
356  G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
357  G4int vecLen = 0;
358  vec.Initialize( 0 );
359 
360 
361 
362  // I hope my understanding of "secondary" is correct here
363  // I think that it entails only what is not a surviving incident of target
364 
365  for (G4int i = 0; i!=NumberOfSecondaries;i++){
366  if(theParticleDefinitions[i]!=aTarget
367  && theParticleDefinitions[i]!=originalIncident->GetDefinition()
368  && theParticleDefinitions[i]!=outgoingRhadron
369  && theParticleDefinitions[i]!=outgoingTarget)
370  {
371  G4ReactionProduct* pa = new G4ReactionProduct;
372  pa->SetDefinition( theParticleDefinitions[i] );
373  (G4UniformRand() < 0.5) ? pa->SetSide( -1 ) : pa->SetSide( 1 );
374  vec.SetElement( vecLen++, pa );
375  }
376  }
377 
378  // if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingRhadron );
379 
380  if (incidentHasChanged) currentParticle.SetDefinitionAndUpdateE( outgoingCloud );
381  if (incidentHasChanged) modifiedOriginal.SetDefinition( outgoingCloud );//Is this correct? It solves the "free energy" checking in ReactionDynamics
382  if (targetHasChanged) targetParticle.SetDefinitionAndUpdateE( outgoingTarget );
383 
384  // G4cout<<"Calling CalculateMomenta... "<<G4endl;
385  /*
386  G4cout<<"Current particle starts as: "<<currentParticle.GetDefinition()->GetParticleName()<<G4endl;
387  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
388  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
389  */
390  // G4double e_temp = currentParticle.GetKineticEnergy();
391 
392  CalculateMomenta( vec, vecLen,
393  originalIncident, originalTarget, modifiedOriginal,
394  targetNucleus, currentParticle, targetParticle,
395  incidentHasChanged, targetHasChanged, quasiElastic );
396 
397  // G4cout <<"Cloud loss: "<<(e_temp-currentParticle.GetKineticEnergy())/GeV<<" GeV"<<G4endl;
398 
399  G4String cPname = currentParticle.GetDefinition()->GetParticleName();
400 
401  // if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud")
402  // {
403  /*
404  G4cout<<"Current particle is now: "<<cPname <<G4endl;
405  G4cout<<"with momentum: "<<currentParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
406  G4cout<<"and kinetic energy: "<<currentParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
407  G4cout<<"May be killed?: "<<currentParticle.GetMayBeKilled()<<G4endl;
408  G4cout<<"Modified original is: "<<modifiedOriginal.GetDefinition()->GetParticleName()<<G4endl;
409  G4cout<<"with momentum: "<<modifiedOriginal.GetMomentum()/GeV<<" GeV"<<G4endl;
410  G4cout<<"and kinetic energy: "<<modifiedOriginal.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
411  G4cout<<"May be killed?: "<<modifiedOriginal.GetMayBeKilled()<<G4endl;
412  G4cout<<"Target particle is: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
413  G4cout<<"with momentum: "<<targetParticle.GetMomentum()/GeV<<" GeV"<<G4endl;
414  G4cout<<"and kinetic energy: "<<targetParticle.GetKineticEnergy()/GeV<<" GeV"<<G4endl;
415  G4cout<<"May be killed?: "<<targetParticle.GetMayBeKilled()<<G4endl;
416  G4cout<<"incidentHasChanged: "<<incidentHasChanged<<G4endl;
417  G4cout<<"targetHasChanged: "<<targetHasChanged<<G4endl;
418  G4cout<<"Particles in vec:"<<G4endl;
419  for(int i=0; i<vecLen; ++i )
420  {
421  G4cout<< vec[i]->GetDefinition()->GetParticleName()<<G4endl;
422  }
423  */
424 
425  // }
426 
427  // G4cout<<"Done!"<<G4endl;
428  // const G4LorentzRotation& trans(originalIncident->GetTrafoToLab());
429  // G4cout<<"Check aParticleChange.GetNumberOfSecondaries(): "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
430 
431  aParticleChange.SetNumberOfSecondaries(vecLen+NumberOfSecondaries);
432  G4double e_kin=0;
433  G4LorentzVector cloud_p4_new; //Cloud 4-momentum in lab after collision
434  // n++;
435  // G4cout << n << G4endl;
436  /*
437  if(cPname!="rhadronmesoncloud"&&cPname!="rhadronbaryoncloud") {
438  G4cout<<"Cloud deleted!!! AAARRRRGGGHHH!!!"<<G4endl;
439  G4cout<<"Cloud name: "<<cPname<<G4endl;
440  G4cout<<"E_kin_0: "<<e_kin_0/GeV<<" GeV"<<G4endl;
441  // G4cout<<"n: "<<n<<G4endl;
442  // n=0;
443  }
444  */
445  cloud_p4_new.setVectM(currentParticle.GetMomentum(),currentParticle.GetMass());
446  cloud_p4_new *= trans;
447 
448  G4LorentzVector cloud_p4_old_full = Cloud4Momentum; //quark system in CMS BEFORE collision
449  cloud_p4_old_full.boost(trafo_full_cms);
450  G4LorentzVector cloud_p4_old_cloud = Cloud4Momentum; //quark system in cloud CMS BEFORE collision
451  cloud_p4_old_cloud.boost(trafo_cloud_cms);
452  G4LorentzVector cloud_p4_full = cloud_p4_new; //quark system in CMS AFTER collision
453  cloud_p4_full.boost(trafo_full_cms);
454  G4LorentzVector cloud_p4_cloud = cloud_p4_new; //quark system in cloud CMS AFTER collision
455  cloud_p4_cloud.boost(trafo_cloud_cms);
456 
457  G4LorentzVector p_g_cms = gluinoMomentum; //gluino in CMS BEFORE collision
458  p_g_cms.boost(trafo_full_cms);
459 
460  G4LorentzVector p4_new;
461  // p4_new.setVectM(cloud_p4_full.v()+p_g_cms.v(),outgoingRhadron->GetPDGMass());
462  // p4_new.boost(-trafo_full_cms);
463  // p4_new = cloud_p4_new + gluinoMomentum;
464  p4_new.setVectM( cloud_p4_new.v() + gluinoMomentum.v(), outgoingRhadron->GetPDGMass() );
465  // G4cout<<"P4-diff: "<<(p4_new-cloud_p4_new-gluinoMomentum)/GeV<<", magnitude: "<<(p4_new-cloud_p4_new-gluinoMomentum).m()/MeV<<" MeV" <<G4endl;
466 
467  G4ThreeVector p_new;
468  p_new = p4_new.vect();
469 
470  aParticleChange.ProposeLocalEnergyDeposit((p4_new-cloud_p4_new-gluinoMomentum).m());
471 
472  if( incidentHasChanged )
473  {
474  G4DynamicParticle* p0 = new G4DynamicParticle;
475  p0->SetDefinition( outgoingRhadron );
476  p0->SetMomentum( p_new );
477 
478  // May need to run SetDefinitionAndUpdateE here...
479  G4Track* Track0 = new G4Track(p0,
480  aTrack.GetGlobalTime(),
481  aPosition);
482  Track0->SetTouchableHandle(thisTouchable);
483  aParticleChange.AddSecondary(Track0);
484  /*
485  G4cout<<"Adding a particle "<<p0->GetDefinition()->GetParticleName()<<G4endl;
486  G4cout<<"with momentum: "<<p0->GetMomentum()/GeV<<" GeV"<<G4endl;
487  G4cout<<"and kinetic energy: "<<p0->GetKineticEnergy()/GeV<<" GeV"<<G4endl;
488  */
489  if(p0->GetKineticEnergy()>e_kin_0) {
490  G4cout<<"ALAAAAARM!!! (incident changed from "
491  <<IncidentRhadron->GetDefinition()->GetParticleName()
492  <<" to "<<p0->GetDefinition()->GetParticleName()<<")"<<G4endl;
493  G4cout<<"Energy loss: "<<(e_kin_0-p0->GetKineticEnergy())/GeV<<" GeV (should be positive)"<<G4endl;
494  //Turns out problem is only in 2 -> 3 (Won't fix 2 -> 2 energy deposition)
495  if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
496  } /*else {
497  G4cout<<"NO ALAAAAARM!!!"<<G4endl;
498  }*/
499  if(std::abs(p0->GetKineticEnergy()-e_kin_0)>100*GeV) {
500  G4cout<<"Diff. too big"<<G4endl;
501  }
502 
503  aParticleChange.ProposeTrackStatus( fStopAndKill );
504  }
505  else
506  {
507 
508  G4double p = p_new.mag();
509  if( p > DBL_MIN )
510  aParticleChange.ProposeMomentumDirection( p_new.x()/p, p_new.y()/p, p_new.z()/p );
511  else
512  aParticleChange.ProposeMomentumDirection( 1.0, 0.0, 0.0 );
513 
514  G4double aE = sqrt(p*p+(outgoingRhadron->GetPDGMass()*outgoingRhadron->GetPDGMass()) );
515  e_kin = aE - outgoingRhadron->GetPDGMass();
516  /* AR_NEWCODE_IMPORT
517  if(e_kin>e_kin_0) {
518  G4cout<<"ALAAAAARM!!!"<<G4endl;
519  G4cout<<"Energy loss: "<<(e_kin_0-e_kin)/GeV<<" GeV (should be positive)"<<G4endl;
520  if(rp.size()!=3) G4cout<<"DOUBLE ALAAAAARM!!!"<<G4endl;
521  }
522  if(std::abs(e_kin-e_kin_0)>100*GeV) {
523  G4cout<<"Diff. too big"<<G4endl;
524  }
525 
526  if (std::fabs(aE)<.1*eV) aE=.1*eV;
527  aParticleChange.ProposeEnergy( aE- outgoingRhadron->GetPDGMass() ); //I do not know if this is correct...
528  if(std::abs(e_kin-e_kin_0)>100*GeV) {
529  G4cout<<"Diff. too big"<<G4endl;
530  }
531  */
532  }
533 
534  // return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep);
535  if( targetParticle.GetMass() > 0.0 ) // targetParticle can be eliminated in TwoBody
536  {
537  G4DynamicParticle *p1 = new G4DynamicParticle;
538  p1->SetDefinition( targetParticle.GetDefinition() );
539  // G4cout<<"Target secondary: "<<targetParticle.GetDefinition()->GetParticleName()<<G4endl;
540  G4ThreeVector momentum = targetParticle.GetMomentum();
541  momentum = momentum.rotate(cache,what);
542  p1->SetMomentum( momentum );
543  p1->SetMomentum( (trans*p1->Get4Momentum()).vect());
544  G4Track* Track1 = new G4Track(p1,
545  aTrack.GetGlobalTime(),
546  aPosition);
547  Track1->SetTouchableHandle(thisTouchable);
548  aParticleChange.AddSecondary(Track1);
549  }
550  G4DynamicParticle *pa;
551  /*
552  G4cout<<"vecLen: "<<vecLen<<G4endl;
553  G4cout<<"#sec's: "<<aParticleChange.GetNumberOfSecondaries()<<G4endl;
554  */
555 
556 
557 
558  for(int i=0; i<vecLen; ++i )
559  {
560  pa = new G4DynamicParticle();
561  pa->SetDefinition( vec[i]->GetDefinition() );
562  pa->SetMomentum( vec[i]->GetMomentum() );
563  pa->Set4Momentum(trans*(pa->Get4Momentum()));
564  G4ThreeVector pvec = pa->GetMomentum();
565  G4Track* Trackn = new G4Track(pa,
566  aTrack.GetGlobalTime(),
567  aPosition);
568  Trackn->SetTouchableHandle(thisTouchable);
569  aParticleChange.AddSecondary(Trackn);
570 
571  // debug
572 
573 // G4cerr << "FullModelHadronicProcess: New secondary " << i
574 // << " ID " << Trackn->GetTrackID()
575 // << " PDG " << Trackn->GetDefinition()->GetParticleName()
576 // << " position " << Trackn->GetPosition()
577 // << " volume " << Trackn->GetTouchable()
578 // << " handle " << Trackn->GetTouchableHandle() << G4endl;
579 
580  delete vec[i];
581  }
582 
583  // Histogram filling
584  const G4DynamicParticle* theRhadron = FindRhadron(&aParticleChange);
585 
586 
587  if (theRhadron!=NULL||IncidentSurvives)
588  {
589 
590  double E_new;
591  if(IncidentSurvives)
592  {
593  // E_new = currentParticle.GetKineticEnergy();
594  E_new = e_kin;
595  } else {
596  E_new = theRhadron->GetKineticEnergy();
597  if(CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
598  !=CustomPDGParser::s_isRMeson(theIncidentPDG)
599  ||
600  CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
601  !=CustomPDGParser::s_isMesonino(theIncidentPDG)
602  ) {
603 
604  G4cout<<"Rm: "<<CustomPDGParser::s_isRMeson(theRhadron->GetDefinition()->GetPDGEncoding())
605  <<" vs: "<<CustomPDGParser::s_isRMeson(theIncidentPDG)<<G4endl;
606  G4cout<<"Sm: "<<CustomPDGParser::s_isMesonino(theRhadron->GetDefinition()->GetPDGEncoding())
607  <<" vs: "<<CustomPDGParser::s_isMesonino(theIncidentPDG)<<G4endl;
608 
609  }
610  }
611 
612  //Calculating relevant scattering angles.
613  G4LorentzVector p4_old_full = FullRhadron4Momentum; //R-hadron in CMS BEFORE collision
614  p4_old_full.boost(trafo_full_cms);
615  G4LorentzVector p4_old_cloud = FullRhadron4Momentum; //R-hadron in cloud CMS BEFORE collision
616  p4_old_cloud.boost(trafo_cloud_cms);
617  G4LorentzVector p4_full = p4_new; //R-hadron in CMS AFTER collision
618  // G4cout<<p4_full.v()/GeV<<G4endl;
619  p4_full=p4_full.boost(trafo_full_cms);
620  // G4cout<<p4_full.m()<<" / "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).m()<<G4endl;
621  G4LorentzVector p4_cloud = p4_new; //R-hadron in cloud CMS AFTER collision
622  p4_cloud.boost(trafo_cloud_cms);
623 
624 
625  /*
626  G4double dtheta_fullcms = p4_full.vect().angle(p4_old_full.vect());
627  G4double dtheta_cloudcms = p4_cloud.vect().angle(p4_old_cloud.vect());
628  G4double dtheta_lab = p_new.angle(p_0);//acos(p_0*p_new/(p_0.mag()*p_new.mag()));
629 
630  G4double cloud_dtheta_fullcms = cloud_p4_full.vect().angle(cloud_p4_old_full.vect());
631  G4double cloud_dtheta_cloudcms = cloud_p4_cloud.vect().angle(p4_old_cloud.vect());
632  G4double cloud_dtheta_lab = cloud_p4_new.vect().angle(p_0);
633  //Writing out momenta for manual check of boosts:
634  G4cout<<"******************************************"<<G4endl;
635 
636  G4cout<<"R-hadron, before: "<<G4endl;
637  G4cout<<"Lab: "<<FullRhadron4Momentum.v().mag()/GeV<<" GeV"<<G4endl;
638  G4cout<<"CMS: "<<p4_old_full.v().mag()/GeV<<" GeV"<<G4endl;
639  G4cout<<"after: "<<G4endl;
640  G4cout<<"Lab: "<<p4_new.v().mag()/GeV<<" GeV"<<G4endl;
641  G4cout<<"CMS: "<<p4_full.v().mag()/GeV<<" GeV"<<G4endl;
642  G4cout<<"cos(theta*): "<<cos(dtheta_fullcms)<<G4endl;
643  G4cout<<"Gluino: "<<G4endl;
644  G4cout<<"Lab: "<<gluinoMomentum.v().mag()/GeV<<" GeV"<<G4endl;
645  G4cout<<"CMS: "<<p_g_cms.v().mag()/GeV<<" GeV"<<G4endl;
646 
647  G4cout<<"Cloud, before: "<<G4endl;
648  G4cout<<"Lab: "<<Cloud4Momentum.v().mag()/GeV<<" GeV"<<G4endl;
649  G4cout<<"CMS: "<<cloud_p4_old_full.v().mag()/GeV<<" GeV"<<G4endl;
650  G4cout<<"after: "<<G4endl;
651  G4cout<<"CMS: "<<cloud_p4_full.v().mag()/GeV<<" GeV"<<G4endl;
652  G4cout<<"cloud cos(theta*): "<<cos(cloud_dtheta_fullcms)<<G4endl;
653  G4cout<<"Longitudinal: "<<cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()/GeV<<" GeV"<<G4endl;
654  G4cout<<"~Combined longitudinal: "<<(cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()+p4_old_full.v().mag())/GeV<<" GeV"<<G4endl;
655  if ((cos(cloud_dtheta_fullcms)*cloud_p4_full.v().mag()+p4_old_full.v().mag())<0.) G4cout<<"ALAAARM"<<G4endl;
656  G4cout<<"Psum cos(theta*): "<<cos((p_g_cms.v()+cloud_p4_full.v()).angle(p4_old_full.v()))<<G4endl;
657  G4cout<<"True R-hadron (CMS): "<<(p_g_cms.v()+cloud_p4_full.v()).mag()/GeV<<" GeV"<<G4endl;
658  G4cout<<"******************************************"<<G4endl;
659  G4cout<<"Fucking manual fucking calculation:"<<G4endl;
660  G4cout<<"Momenta:"<<G4endl;
661  G4cout<<"Cloud, lab: "<<cloud_p4_new.v()/GeV<<" + gluino: "<<gluinoMomentum.v()/GeV
662  <<" = "<<(cloud_p4_new.v()+gluinoMomentum.v())/GeV
663  <<". Boost to CMS: "<<(cloud_p4_new+gluinoMomentum).boost(trafo_full_cms).v()/GeV<<G4endl;
664  G4cout<<"Cloud, CMS: "<<cloud_p4_full.v()/GeV<<" + gluino: "<<p_g_cms.v()/GeV
665  <<" = "<<(cloud_p4_full.v()+p_g_cms.v())/GeV
666  <<". Boost to Lab: "<<(cloud_p4_full+p_g_cms).boost(-trafo_full_cms).v()/GeV<<G4endl;
667  G4cout<<"Ref: "<<p4_new.v()/GeV<<" / "<<p4_full.v()/GeV<<G4endl;
668  G4cout<<"******************************************"<<G4endl;
669  */
670 
671 
672  // std::cout<<"Lab, fullcms, cloudcms: "<<dtheta_lab<<", "<<dtheta_fullcms<<", "<<dtheta_cloudcms<<std::endl;
673  G4double AbsDeltaE = E_0-E_new;
674  // G4cout <<"Energy loss: "<<AbsDeltaE/GeV<<G4endl;
675  if (AbsDeltaE > 10*GeV) {
676  G4cout<<"Energy loss larger than 10 GeV..."<<G4endl;
677  G4cout<<"E_0: "<<E_0/GeV<<" GeV"<<G4endl;
678  G4cout<<"E_new: "<<E_new/GeV<<" GeV"<<G4endl;
679  G4cout<<"Gamma: "<<IncidentRhadron->GetTotalEnergy()/IncidentRhadron->GetDefinition()->GetPDGMass()<<G4endl;
680  G4cout<<"x: "<<aPosition.x()/cm<<" cm"<<G4endl;
681  }
682  }
683  delete originalIncident;
684  delete originalTarget;
685  // aParticleChange.DumpInfo();
686  // G4cout << "Exiting FullModelHadronicProcess::PostStepDoIt"<<G4endl;
687 
688  //clear interaction length
689  ClearNumberOfInteractionLengthLeft();
690 
691 
692  return &aParticleChange;
693 
694 }
695 
696 
697 void FullModelHadronicProcess::CalculateMomenta(
698  G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
699  G4int &vecLen,
700  const G4HadProjectile *originalIncident, // the original incident particle
701  const G4DynamicParticle *originalTarget,
702  G4ReactionProduct &modifiedOriginal, // Fermi motion and evap. effects included
703  G4Nucleus &targetNucleus,
704  G4ReactionProduct &currentParticle,
705  G4ReactionProduct &targetParticle,
706  G4bool &incidentHasChanged,
707  G4bool &targetHasChanged,
708  G4bool quasiElastic )
709 {
710  FullModelReactionDynamics theReactionDynamics;
711 
712  cache = 0;
713  what = originalIncident->Get4Momentum().vect();
714 
715  //Commented out like in casqmesmin.F where CALL STPAIR is commented out
716  /*
717  theReactionDynamics.ProduceStrangeParticlePairs( vec, vecLen,
718  modifiedOriginal, originalTarget,
719  currentParticle, targetParticle,
720  incidentHasChanged, targetHasChanged );
721  */
722 
723  if( quasiElastic )
724  {
725  // G4cout<<"We are calling TwoBody..."<<G4endl;
726  theReactionDynamics.TwoBody( vec, vecLen,
727  modifiedOriginal, originalTarget,
728  currentParticle, targetParticle,
729  targetNucleus, targetHasChanged );
730 
731  return;
732  }
733 
734  //If ProduceStrangeParticlePairs is commented out, let's cut this one as well
735  G4ReactionProduct leadingStrangeParticle;
736  G4bool leadFlag = MarkLeadingStrangeParticle( currentParticle,
737  targetParticle,
738  leadingStrangeParticle );
739 
740 
741 
742  //
743  // Note: the number of secondaries can be reduced in GenerateXandPt and TwoCluster
744  //
745  G4bool finishedGenXPt = false;
746  G4bool annihilation = false;
747  if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
748  currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
749  {
750  // original was an anti-particle and annihilation has taken place
751  annihilation = true;
752  G4double ekcor = 1.0;
753  G4double ek = originalIncident->GetKineticEnergy();
754  G4double ekOrg = ek;
755 
756  const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
757  if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
758  const G4double atomicWeight = targetNucleus.GetN();
759  ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
760 
761  G4double tkin = targetNucleus.Cinema( ek );
762  ek += tkin;
763  ekOrg += tkin;
764  modifiedOriginal.SetKineticEnergy( ekOrg );
765  }
766 
767  const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
768  G4double rand1 = G4UniformRand();
769  G4double rand2 = G4UniformRand();
770  if( (annihilation || (vecLen >= 6) || (modifiedOriginal.GetKineticEnergy()/GeV >= 1.0)) &&
771  (((originalIncident->GetDefinition() == G4KaonPlus::KaonPlus()) ||
772  (originalIncident->GetDefinition() == G4KaonMinus::KaonMinus()) ||
773  (originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong()) ||
774  (originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort())) &&
775  ((rand1 < 0.5) || (rand2 > twsup[vecLen]))))
776  finishedGenXPt =
777  theReactionDynamics.GenerateXandPt( vec, vecLen,
778  modifiedOriginal, originalIncident,
779  currentParticle, targetParticle,
780  targetNucleus, incidentHasChanged,
781  targetHasChanged, leadFlag,
782  leadingStrangeParticle );
783  if( finishedGenXPt )
784  {
785  Rotate(vec, vecLen);
786  return;
787  }
788 
789  G4bool finishedTwoClu = false;
790  if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
791  {
792 
793  for(G4int i=0; i<vecLen; i++) delete vec[i];
794  vecLen = 0;
795  }
796  else
797  {
798 
799  theReactionDynamics.SuppressChargedPions( vec, vecLen,
800  modifiedOriginal, currentParticle,
801  targetParticle, targetNucleus,
802  incidentHasChanged, targetHasChanged );
803 
804  try
805  {
806  finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
807  modifiedOriginal, originalIncident,
808  currentParticle, targetParticle,
809  targetNucleus, incidentHasChanged,
810  targetHasChanged, leadFlag,
811  leadingStrangeParticle );
812  }
813  catch(G4HadReentrentException aC)
814  {
815  aC.Report(G4cout);
816  throw G4HadReentrentException(__FILE__, __LINE__, "Failing to calculate momenta");
817  }
818  }
819  if( finishedTwoClu )
820  {
821  Rotate(vec, vecLen);
822  return;
823  }
824 
825  //
826  // PNBlackTrackEnergy is the kinetic energy available for
827  // proton/neutron black track particles [was enp(1) in fortran code]
828  // DTABlackTrackEnergy is the kinetic energy available for
829  // deuteron/triton/alpha particles [was enp(3) in fortran code]
830  //const G4double pnCutOff = 0.1;
831  //const G4double dtaCutOff = 0.1;
832  //if( (targetNucleus.GetN() >= 1.5)
833  // && !(incidentHasChanged || targetHasChanged)
834  // && (targetNucleus.GetPNBlackTrackEnergy()/MeV <= pnCutOff)
835  // && (targetNucleus.GetDTABlackTrackEnergy()/MeV <= dtaCutOff) )
836  //{
837  // the atomic weight of the target nucleus is >= 1.5 AND
838  // neither the incident nor the target particles have changed AND
839  // there is no kinetic energy available for either proton/neutron
840  // or for deuteron/triton/alpha black track particles
841  // For diffraction scattering on heavy nuclei use elastic routines instead
842  //G4cerr << "*** Error in G4InelasticInteraction::CalculateMomenta" << G4endl;
843  //G4cerr << "*** the elastic scattering would be better here ***" <<G4endl;
844  //}
845  theReactionDynamics.TwoBody( vec, vecLen,
846  modifiedOriginal, originalTarget,
847  currentParticle, targetParticle,
848  targetNucleus, targetHasChanged );
849 }
850 
851 
852 G4bool FullModelHadronicProcess::MarkLeadingStrangeParticle(
853  const G4ReactionProduct &currentParticle,
854  const G4ReactionProduct &targetParticle,
855  G4ReactionProduct &leadParticle )
856 {
857  // the following was in GenerateXandPt and TwoCluster
858  // add a parameter to the GenerateXandPt function telling it about the strange particle
859  //
860  // assumes that the original particle was a strange particle
861  //
862  G4bool lead = false;
863  if( (currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
864  (currentParticle.GetDefinition() != G4Proton::Proton()) &&
865  (currentParticle.GetDefinition() != G4Neutron::Neutron()) )
866  {
867  lead = true;
868  leadParticle = currentParticle; // set lead to the incident particle
869  }
870  else if( (targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
871  (targetParticle.GetDefinition() != G4Proton::Proton()) &&
872  (targetParticle.GetDefinition() != G4Neutron::Neutron()) )
873  {
874  lead = true;
875  leadParticle = targetParticle; // set lead to the target particle
876  }
877  return lead;
878 }
879 
880 void FullModelHadronicProcess::Rotate(G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, G4int &vecLen)
881 {
882  G4double rotation = 2.*pi*G4UniformRand();
883  cache = rotation;
884  G4int i;
885  for( i=0; i<vecLen; ++i )
886  {
887  G4ThreeVector momentum = vec[i]->GetMomentum();
888  momentum = momentum.rotate(rotation, what);
889  vec[i]->SetMomentum(momentum);
890  }
891 }
892 
893 const G4DynamicParticle* FullModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
894 {
895  G4int nsec = aParticleChange->GetNumberOfSecondaries();
896  if (nsec==0) return 0;
897  int i = 0;
898  G4bool found = false;
899  while (i!=nsec && !found){
900  // G4cout<<"Checking "<<aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleName()<<G4endl;
901  // if (aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition()->GetParticleType()=="rhadron") found = true;
902  if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found = true;
903  i++;
904  }
905  i--;
906  if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
907  return 0;
908 }
909 
int i
Definition: DBlmapReader.cc:9
tuple pp
Definition: createTree.py:15
#define abs(x)
Definition: mlp_lapack.h:159
G4ParticleDefinition * GetCloud()
#define NULL
Definition: scimark2.h:8
static bool s_isRMeson(int pdg)
static bool s_isMesonino(int pdg)
T sqrt(T t)
Definition: SSEVec.h:46
G4ParticleDefinition * GetSpectator()
double p1[4]
Definition: TauolaWrapper.h:89
tuple cout
Definition: gather_cfg.py:121
double pi