Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4INCLXXInterface.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // INCL++ intra-nuclear cascade model
27 // Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28 // Davide Mancusi, CEA
29 // Alain Boudard, CEA
30 // Sylvie Leray, CEA
31 // Joseph Cugnon, University of Liege
32 //
33 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
37 #include <cmath>
38 
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4INCLXXInterface.hh"
42 #include "G4GenericIon.hh"
43 #include "G4INCLCascade.hh"
45 #include "G4ReactionProduct.hh"
47 #include "G4String.hh"
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
51 #include "G4INCLVersion.hh"
52 
55  theINCLModel(NULL),
56  thePreCompoundModel(aPreCompound),
57  theInterfaceStore(G4INCLXXInterfaceStore::GetInstance()),
58  complainedAboutBackupModel(false),
59  complainedAboutPreCompound(false),
60  theIonTable(G4IonTable::GetIonTable())
61 {
62  if(!thePreCompoundModel) {
65  thePreCompoundModel = static_cast<G4VPreCompoundModel*>(p);
66  if(!thePreCompoundModel) { thePreCompoundModel = new G4PreCompoundModel; }
67  }
68 
69  // Use the environment variable G4INCLXX_NO_DE_EXCITATION to disable de-excitation
70  if(getenv("G4INCLXX_NO_DE_EXCITATION")) {
71  G4String message = "de-excitation is completely disabled!";
72  theInterfaceStore->EmitWarning(message);
73  theDeExcitation = 0;
74  } else {
77  theDeExcitation = static_cast<G4VPreCompoundModel*>(p);
79  }
80 
81  // use the envvar G4INCLXX_DUMP_REMNANT to dump information about the
82  // remnants on stdout
83  if(getenv("G4INCLXX_DUMP_REMNANT"))
84  dumpRemnantInfo = true;
85  else
86  dumpRemnantInfo = false;
87 
88  theBackupModel = new G4BinaryLightIonReaction;
89  theBackupModelNucleon = new G4BinaryCascade;
90 }
91 
93 {
94 }
95 
96 G4bool G4INCLXXInterface::AccurateProjectile(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) const {
97  // Use direct kinematics if the projectile is a nucleon or a pion
98  const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
99  if(projectileDef == G4Proton::Proton()
100  || projectileDef == G4Neutron::Neutron()
101  || projectileDef == G4PionPlus::PionPlus()
102  || projectileDef == G4PionZero::PionZero()
103  || projectileDef == G4PionMinus::PionMinus())
104  return false;
105 
106  // Here all projectiles should be nuclei
107  const G4int pA = projectileDef->GetAtomicMass();
108  if(pA<=0) {
109  std::stringstream ss;
110  ss << "the model does not know how to handle a collision between a "
111  << projectileDef->GetParticleName() << " projectile and a Z="
112  << theNucleus.GetZ_asInt() << ", A=" << theNucleus.GetA_asInt();
113  theInterfaceStore->EmitWarning(ss.str());
114  return true;
115  }
116 
117  // If either nucleus is a LCP (A<=4), run the collision as light on heavy
118  const G4int tA = theNucleus.GetA_asInt();
119  if(tA<=4 || pA<=4) {
120  if(pA<tA)
121  return false;
122  else
123  return true;
124  }
125 
126  // If one of the nuclei is heavier than theMaxProjMassINCL, run the collision
127  // as light on heavy.
128  // Note that here we are sure that either the projectile or the target is
129  // smaller than theMaxProjMass; otherwise theBackupModel would have been
130  // called.
131  const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
132  if(pA > theMaxProjMassINCL)
133  return true;
134  else if(tA > theMaxProjMassINCL)
135  return false;
136  else
137  // In all other cases, use the global setting
138  return theInterfaceStore->GetAccurateProjectile();
139 }
140 
142 {
143  // For reactions on nucleons, use the backup model (without complaining)
144  if(aTrack.GetDefinition()->GetAtomicMass()<=1 && theNucleus.GetA_asInt()<=1) {
145  return theBackupModelNucleon->ApplyYourself(aTrack, theNucleus);
146  }
147 
148  // For systems heavier than theMaxProjMassINCL, use another model (typically
149  // BIC)
150  const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
151  if(aTrack.GetDefinition()->GetAtomicMass() > theMaxProjMassINCL
152  && theNucleus.GetA_asInt() > theMaxProjMassINCL) {
153  if(!complainedAboutBackupModel) {
154  complainedAboutBackupModel = true;
155  std::stringstream ss;
156  ss << "INCL++ refuses to handle reactions between nuclei with A>"
157  << theMaxProjMassINCL
158  << ". A backup model ("
159  << theBackupModel->GetModelName()
160  << ") will be used instead.";
161  theInterfaceStore->EmitBigWarning(ss.str());
162  }
163  return theBackupModel->ApplyYourself(aTrack, theNucleus);
164  }
165 
166  // For energies lower than cascadeMinEnergyPerNucleon, use PreCompound
167  const G4double cascadeMinEnergyPerNucleon = theInterfaceStore->GetCascadeMinEnergyPerNucleon();
168  const G4double trackKinE = aTrack.GetKineticEnergy();
169  const G4ParticleDefinition *trackDefinition = aTrack.GetDefinition();
170  if((trackDefinition==G4Neutron::NeutronDefinition() || trackDefinition==G4Proton::ProtonDefinition())
171  && trackKinE < cascadeMinEnergyPerNucleon) {
172  if(!complainedAboutPreCompound) {
173  complainedAboutPreCompound = true;
174  std::stringstream ss;
175  ss << "INCL++ refuses to handle nucleon-induced reactions below "
176  << cascadeMinEnergyPerNucleon / MeV
177  << " MeV. A PreCoumpound model ("
178  << thePreCompoundModel->GetModelName()
179  << ") will be used instead.";
180  theInterfaceStore->EmitBigWarning(ss.str());
181  }
182  return thePreCompoundModel->ApplyYourself(aTrack, theNucleus);
183  }
184 
185  // Calculate the total four-momentum in the entrance channel
186  const G4int nucleusA = theNucleus.GetA_asInt();
187  const G4int nucleusZ = theNucleus.GetZ_asInt();
188  const G4double theNucleusMass = theIonTable->GetIonMass(nucleusZ, nucleusA);
189  const G4double theTrackMass = trackDefinition->GetPDGMass();
190  const G4double theTrackEnergy = trackKinE + theTrackMass;
191  const G4double theTrackMomentumAbs2 = theTrackEnergy*theTrackEnergy - theTrackMass*theTrackMass;
192  const G4double theTrackMomentumAbs = ((theTrackMomentumAbs2>0.0) ? std::sqrt(theTrackMomentumAbs2) : 0.0);
193  const G4ThreeVector theTrackMomentum = aTrack.Get4Momentum().getV().unit() * theTrackMomentumAbs;
194 
195  G4LorentzVector goodTrack4Momentum(theTrackMomentum, theTrackEnergy);
196  G4LorentzVector fourMomentumIn;
197  fourMomentumIn.setE(theTrackEnergy + theNucleusMass);
198  fourMomentumIn.setVect(theTrackMomentum);
199 
200  // Check if inverse kinematics should be used
201  const G4bool inverseKinematics = AccurateProjectile(aTrack, theNucleus);
202 
203  // If we are running in inverse kinematics, redefine aTrack and theNucleus
204  G4LorentzRotation *toInverseKinematics = NULL;
205  G4LorentzRotation *toDirectKinematics = NULL;
206  G4HadProjectile const *aProjectileTrack = &aTrack;
207  G4Nucleus *theTargetNucleus = &theNucleus;
208  if(inverseKinematics) {
209  G4ParticleDefinition *oldTargetDef = theIonTable->GetIon(nucleusZ, nucleusA, 0);
210  const G4ParticleDefinition *oldProjectileDef = aTrack.GetDefinition();
211 
212  if(oldProjectileDef != 0 && oldTargetDef != 0) {
213  const G4int newTargetA = oldProjectileDef->GetAtomicMass();
214  const G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
215 
216  if(newTargetA > 0 && newTargetZ > 0) {
217  // This should give us the same energy per nucleon
218  theTargetNucleus = new G4Nucleus(newTargetA, newTargetZ);
219  toInverseKinematics = new G4LorentzRotation(goodTrack4Momentum.boostVector());
220  G4LorentzVector theProjectile4Momentum(0.0, 0.0, 0.0, theNucleusMass);
221  G4DynamicParticle swappedProjectileParticle(oldTargetDef, (*toInverseKinematics) * theProjectile4Momentum);
222  aProjectileTrack = new G4HadProjectile(swappedProjectileParticle);
223  } else {
224  G4String message = "badly defined target after swapping. Falling back to normal (non-swapped) mode.";
225  theInterfaceStore->EmitWarning(message);
226  toInverseKinematics = new G4LorentzRotation;
227  }
228  } else {
229  G4String message = "oldProjectileDef or oldTargetDef was null";
230  theInterfaceStore->EmitWarning(message);
231  toInverseKinematics = new G4LorentzRotation;
232  }
233  }
234 
235  // INCL assumes the projectile particle is going in the direction of
236  // the Z-axis. Here we construct proper rotation to convert the
237  // momentum vectors of the outcoming particles to the original
238  // coordinate system.
239  G4LorentzVector projectileMomentum = aProjectileTrack->Get4Momentum();
240 
241  // INCL++ assumes that the projectile is going in the direction of
242  // the z-axis. In principle, if the coordinate system used by G4
243  // hadronic framework is defined differently we need a rotation to
244  // transform the INCL++ reaction products to the appropriate
245  // frame. Please note that it isn't necessary to apply this
246  // transform to the projectile because when creating the INCL++
247  // projectile particle; \see{toINCLKineticEnergy} needs to use only the
248  // projectile energy (direction is simply assumed to be along z-axis).
249  G4RotationMatrix toZ;
250  toZ.rotateZ(-projectileMomentum.phi());
251  toZ.rotateY(-projectileMomentum.theta());
252  G4RotationMatrix toLabFrame3 = toZ.inverse();
253  G4LorentzRotation toLabFrame(toLabFrame3);
254  // However, it turns out that the projectile given to us by G4
255  // hadronic framework is already going in the direction of the
256  // z-axis so this rotation is actually unnecessary. Both toZ and
257  // toLabFrame turn out to be unit matrices as can be seen by
258  // uncommenting the folowing two lines:
259  // G4cout <<"toZ = " << toZ << G4endl;
260  // G4cout <<"toLabFrame = " << toLabFrame << G4endl;
261 
262  theResult.Clear(); // Make sure the output data structure is clean.
263  theResult.SetStatusChange(stopAndKill);
264 
265  std::list<G4Fragment> remnants;
266 
267  const G4int maxTries = 200;
268  G4int nTries = 0;
269  // INCL can generate transparent events. However, this is meaningful
270  // only in the standalone code. In Geant4 we must "force" INCL to
271  // produce a valid cascade.
272  G4bool eventIsOK = false;
273  do {
274  const G4INCL::ParticleSpecies theSpecies = toINCLParticleSpecies(*aProjectileTrack);
275  const G4double kineticEnergy = toINCLKineticEnergy(*aProjectileTrack);
276 
277  // The INCL model will be created at the first use
279 
280  const G4INCL::EventInfo eventInfo = theINCLModel->processEvent(theSpecies, kineticEnergy, theTargetNucleus->GetA_asInt(), theTargetNucleus->GetZ_asInt());
281  // eventIsOK = !eventInfo.transparent && nTries < maxTries;
282  eventIsOK = !eventInfo.transparent;
283  if(eventIsOK) {
284 
285  // If the collision was run in inverse kinematics, prepare to boost back
286  // all the reaction products
287  if(inverseKinematics) {
288  toDirectKinematics = new G4LorentzRotation(toInverseKinematics->inverse());
289  }
290 
291  G4LorentzVector fourMomentumOut;
292 
293  for(G4int i = 0; i < eventInfo.nParticles; i++) {
294  G4int A = eventInfo.A[i];
295  G4int Z = eventInfo.Z[i];
296  // G4cout <<"INCL particle A = " << A << " Z = " << Z << G4endl;
297  G4double kinE = eventInfo.EKin[i];
298  G4double px = eventInfo.px[i];
299  G4double py = eventInfo.py[i];
300  G4double pz = eventInfo.pz[i];
301  G4DynamicParticle *p = toG4Particle(A, Z , kinE, px, py, pz);
302  if(p != 0) {
303  G4LorentzVector momentum = p->Get4Momentum();
304 
305  // Apply the toLabFrame rotation
306  momentum *= toLabFrame;
307  // Apply the inverse-kinematics boost
308  if(inverseKinematics) {
309  momentum *= *toDirectKinematics;
310  momentum.setVect(-momentum.vect());
311  }
312 
313  // Set the four-momentum of the reaction products
314  p->Set4Momentum(momentum);
315  fourMomentumOut += momentum;
316  theResult.AddSecondary(p);
317 
318  } else {
319  G4String message = "the model produced a particle that couldn't be converted to Geant4 particle.";
320  theInterfaceStore->EmitWarning(message);
321  }
322  }
323 
324  for(G4int i = 0; i < eventInfo.nRemnants; i++) {
325  const G4int A = eventInfo.ARem[i];
326  const G4int Z = eventInfo.ZRem[i];
327  // G4cout <<"INCL particle A = " << A << " Z = " << Z << G4endl;
328  const G4double kinE = eventInfo.EKinRem[i];
329  const G4double px = eventInfo.pxRem[i];
330  const G4double py = eventInfo.pyRem[i];
331  const G4double pz = eventInfo.pzRem[i];
332  G4ThreeVector spin(
333  eventInfo.jxRem[i]*hbar_Planck,
334  eventInfo.jyRem[i]*hbar_Planck,
335  eventInfo.jzRem[i]*hbar_Planck
336  );
337  const G4double excitationE = eventInfo.EStarRem[i];
338  const G4double nuclearMass = G4NucleiProperties::GetNuclearMass(A, Z) + excitationE;
339  const G4double scaling = remnant4MomentumScaling(nuclearMass,
340  kinE,
341  px, py, pz);
342  G4LorentzVector fourMomentum(scaling * px, scaling * py, scaling * pz,
343  nuclearMass + kinE);
344  if(std::abs(scaling - 1.0) > 0.01) {
345  std::stringstream ss;
346  ss << "momentum scaling = " << scaling
347  << "\n Lorentz vector = " << fourMomentum
348  << ")\n A = " << A << ", Z = " << Z
349  << "\n E* = " << excitationE << ", nuclearMass = " << nuclearMass
350  << "\n remnant i=" << i << ", nRemnants=" << eventInfo.nRemnants
351  << "\n Reaction was: " << aTrack.GetKineticEnergy()/MeV
352  << "-MeV " << aTrack.GetDefinition()->GetParticleName() << " + "
353  << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
354  << ", in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
355  theInterfaceStore->EmitWarning(ss.str());
356  }
357 
358  // Apply the toLabFrame rotation
359  fourMomentum *= toLabFrame;
360  spin *= toLabFrame3;
361  // Apply the inverse-kinematics boost
362  if(inverseKinematics) {
363  fourMomentum *= *toDirectKinematics;
364  fourMomentum.setVect(-fourMomentum.vect());
365  }
366 
367  fourMomentumOut += fourMomentum;
368  G4Fragment remnant(A, Z, fourMomentum);
369  remnant.SetAngularMomentum(spin);
370  if(dumpRemnantInfo) {
371  G4cerr << "G4INCLXX_DUMP_REMNANT: " << remnant << " spin: " << spin << G4endl;
372  }
373  remnants.push_back(remnant);
374  }
375 
376  // Check four-momentum conservation
377  const G4LorentzVector violation4Momentum = fourMomentumOut - fourMomentumIn;
378  const G4double energyViolation = std::abs(violation4Momentum.e());
379  const G4double momentumViolation = violation4Momentum.rho();
380  if(energyViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
381  std::stringstream ss;
382  ss << "energy conservation violated by " << energyViolation/MeV << " MeV in "
383  << aTrack.GetKineticEnergy()/MeV << "-MeV " << aTrack.GetDefinition()->GetParticleName()
384  << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
385  << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
386  theInterfaceStore->EmitWarning(ss.str());
387  eventIsOK = false;
388  const G4int nSecondaries = theResult.GetNumberOfSecondaries();
389  for(G4int j=0; j<nSecondaries; ++j)
390  delete theResult.GetSecondary(j)->GetParticle();
391  theResult.Clear();
392  theResult.SetStatusChange(stopAndKill);
393  remnants.clear();
394  } else if(momentumViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
395  std::stringstream ss;
396  ss << "momentum conservation violated by " << momentumViolation/MeV << " MeV in "
397  << aTrack.GetKineticEnergy()/MeV << "-MeV " << aTrack.GetDefinition()->GetParticleName()
398  << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
399  << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
400  theInterfaceStore->EmitWarning(ss.str());
401  eventIsOK = false;
402  const G4int nSecondaries = theResult.GetNumberOfSecondaries();
403  for(G4int j=0; j<nSecondaries; ++j)
404  delete theResult.GetSecondary(j)->GetParticle();
405  theResult.Clear();
406  theResult.SetStatusChange(stopAndKill);
407  remnants.clear();
408  }
409  }
410  nTries++;
411  } while(!eventIsOK && nTries < maxTries);
412 
413  // Clean up the objects that we created for the inverse kinematics
414  if(inverseKinematics) {
415  delete aProjectileTrack;
416  delete theTargetNucleus;
417  delete toInverseKinematics;
418  delete toDirectKinematics;
419  }
420 
421  if(!eventIsOK) {
422  std::stringstream ss;
423  ss << "maximum number of tries exceeded for the proposed "
424  << aTrack.GetKineticEnergy()/MeV << "-MeV " << aTrack.GetDefinition()->GetParticleName()
425  << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
426  << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
427  theInterfaceStore->EmitWarning(ss.str());
428  theResult.SetStatusChange(isAlive);
429  theResult.SetEnergyChange(aTrack.GetKineticEnergy());
430  theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
431  return &theResult;
432  }
433 
434  // De-excitation:
435 
436  if(theDeExcitation != 0) {
437  for(std::list<G4Fragment>::iterator i = remnants.begin();
438  i != remnants.end(); i++) {
439  G4ReactionProductVector *deExcitationResult = theDeExcitation->DeExcite((*i));
440 
441  for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
442  fragment != deExcitationResult->end(); ++fragment) {
443  G4ParticleDefinition *def = (*fragment)->GetDefinition();
444  if(def != 0) {
445  G4DynamicParticle *theFragment = new G4DynamicParticle(def, (*fragment)->GetMomentum());
446  theResult.AddSecondary(theFragment);
447  }
448  }
449 
450  for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
451  fragment != deExcitationResult->end(); ++fragment) {
452  delete (*fragment);
453  }
454  deExcitationResult->clear();
455  delete deExcitationResult;
456  }
457  }
458 
459  remnants.clear();
460 
461  return &theResult;
462 }
463 
465  return 0;
466 }
467 
468 G4INCL::ParticleType G4INCLXXInterface::toINCLParticleType(G4ParticleDefinition const * const pdef) const {
469  if( pdef == G4Proton::Proton()) return G4INCL::Proton;
470  else if(pdef == G4Neutron::Neutron()) return G4INCL::Neutron;
471  else if(pdef == G4PionPlus::PionPlus()) return G4INCL::PiPlus;
472  else if(pdef == G4PionMinus::PionMinus()) return G4INCL::PiMinus;
473  else if(pdef == G4PionZero::PionZero()) return G4INCL::PiZero;
474  else if(pdef == G4Deuteron::Deuteron()) return G4INCL::Composite;
475  else if(pdef == G4Triton::Triton()) return G4INCL::Composite;
476  else if(pdef == G4He3::He3()) return G4INCL::Composite;
477  else if(pdef == G4Alpha::Alpha()) return G4INCL::Composite;
479  else return G4INCL::UnknownParticle;
480 }
481 
482 G4INCL::ParticleSpecies G4INCLXXInterface::toINCLParticleSpecies(G4HadProjectile const &aTrack) const {
483  const G4ParticleDefinition *pdef = aTrack.GetDefinition();
484  const G4INCL::ParticleType theType = toINCLParticleType(pdef);
485  if(theType!=G4INCL::Composite)
486  return G4INCL::ParticleSpecies(theType);
487  else {
488  G4INCL::ParticleSpecies theSpecies;
489  theSpecies.theType=theType;
490  theSpecies.theA=(G4int) pdef->GetBaryonNumber();
491  theSpecies.theZ=(G4int) pdef->GetPDGCharge();
492  return theSpecies;
493  }
494 }
495 
496 G4double G4INCLXXInterface::toINCLKineticEnergy(G4HadProjectile const &aTrack) const {
497  return aTrack.GetKineticEnergy();
498 }
499 
500 G4ParticleDefinition *G4INCLXXInterface::toG4ParticleDefinition(G4int A, G4int Z) const {
501  if (A == 1 && Z == 1) return G4Proton::Proton();
502  else if(A == 1 && Z == 0) return G4Neutron::Neutron();
503  else if(A == 0 && Z == 1) return G4PionPlus::PionPlus();
504  else if(A == 0 && Z == -1) return G4PionMinus::PionMinus();
505  else if(A == 0 && Z == 0) return G4PionZero::PionZero();
506  else if(A == 2 && Z == 1) return G4Deuteron::Deuteron();
507  else if(A == 3 && Z == 1) return G4Triton::Triton();
508  else if(A == 3 && Z == 2) return G4He3::He3();
509  else if(A == 4 && Z == 2) return G4Alpha::Alpha();
510  else if(A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition
511  return theIonTable->GetIon(Z, A, 0);
512  } else { // Error, unrecognized particle
513  return 0;
514  }
515 }
516 
517 G4DynamicParticle *G4INCLXXInterface::toG4Particle(G4int A, G4int Z,
518  G4double kinE,
519  G4double px,
520  G4double py, G4double pz) const {
521  const G4ParticleDefinition *def = toG4ParticleDefinition(A, Z);
522  if(def == 0) { // Check if we have a valid particle definition
523  return 0;
524  }
525  const G4double energy = kinE * MeV;
526  const G4ThreeVector momentum(px, py, pz);
527  const G4ThreeVector momentumDirection = momentum.unit();
528  G4DynamicParticle *p = new G4DynamicParticle(def, momentumDirection, energy);
529  return p;
530 }
531 
532 G4double G4INCLXXInterface::remnant4MomentumScaling(G4double mass,
533  G4double kineticE,
534  G4double px, G4double py,
535  G4double pz) const {
536  const G4double p2 = px*px + py*py + pz*pz;
537  if(p2 > 0.0) {
538  const G4double pnew2 = kineticE*kineticE + 2.0*kineticE*mass;
539  return std::sqrt(pnew2)/std::sqrt(p2);
540  } else {
541  return 1.0;
542  }
543 }
544 
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
Hep3Vector boostVector() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetMaxProjMassINCL() const
Getter for theMaxProjMassINCL.
G4double GetCascadeMinEnergyPerNucleon() const
Getter for cascadeMinEnergyPerNucleon.
G4HadSecondary * GetSecondary(size_t i)
CLHEP::HepLorentzRotation G4LorentzRotation
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
void SetAngularMomentum(const G4ThreeVector &value)
Definition: G4Fragment.hh:287
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
G4INCLXXInterface(G4VPreCompoundModel *const aPreCompound=0)
const G4String & GetModelName() const
void EmitBigWarning(const G4String &message) const
Emit a BIG warning to G4cout.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
int G4int
Definition: G4Types.hh:78
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Header file for the G4INCLXXInterfaceStore class.
const G4String & GetParticleName() const
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
const G4String & GetIonName(G4int Z, G4int A, G4int lvl=0) const
Definition: G4IonTable.cc:992
HepRotation inverse() const
G4int GetAtomicNumber() const
double phi() const
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
void SetStatusChange(G4HadFinalStateStatus aS)
Short_t ARem[maxSizeRemnants]
Remnant mass number.
std::vector< G4ReactionProduct * > G4ReactionProductVector
static G4INCLXXInterfaceStore * GetInstance()
Get the singleton instance.
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
Singleton class for configuring the INCL++ Geant4 interface.
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
Hep3Vector vect() const
Short_t Z[maxSizeParticles]
Particle charge number.
const G4ParticleDefinition * GetDefinition() const
double theta() const
Short_t nParticles
Number of particles in the final state.
G4INCL::INCL * GetINCLModel()
Get the cached INCL model engine.
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
bool G4bool
Definition: G4Types.hh:79
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1232
const EventInfo & processEvent()
G4double GetKineticEnergy() const
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
double rho() const
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
const G4String & GetParticleType() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4LorentzVector & Get4Momentum() const
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4int GetAtomicMass() const
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
G4LorentzVector Get4Momentum() const
Int_t nRemnants
Number of remnants.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
G4HadronicInteraction * FindModel(const G4String &name)
Short_t A[maxSizeParticles]
Particle mass number.
void Set4Momentum(const G4LorentzVector &momentum)
void SetEnergyChange(G4double anEnergy)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)=0
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
G4double GetPDGMass() const
static G4HadronicInteractionRegistry * Instance()
Hep3Vector unit() const
G4DynamicParticle * GetParticle()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4bool GetAccurateProjectile() const
Getter for accurateProjectile.
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
HepLorentzRotation inverse() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
void setVect(const Hep3Vector &)
Bool_t transparent
True if the event is transparent.
void EmitWarning(const G4String &message)
Emit a warning to G4cout.
double G4double
Definition: G4Types.hh:76
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
G4double GetPDGCharge() const
static G4He3 * He3()
Definition: G4He3.cc:94
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
void SetMomentumChange(const G4ThreeVector &aV)
G4int GetNumberOfSecondaries() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
void AddSecondary(G4DynamicParticle *aP)
G4GLOB_DLL std::ostream G4cerr
Hep3Vector getV() const