Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// INCL++ revision: v5.1.8
34//
35#define INCLXX_IN_GEANT4_MODE 1
36
37#include "globals.hh"
38
39#include <cmath>
40
42#include "G4SystemOfUnits.hh"
43#include "G4INCLXXInterface.hh"
44#include "G4GenericIon.hh"
45#include "G4INCLCascade.hh"
47#include "G4ReactionProduct.hh"
49#include "G4String.hh"
50
53 theINCLModel(NULL),
54 theInterfaceStore(G4INCLXXInterfaceStore::GetInstance()),
55 complainedAboutBackupModel(false)
56{
57 // Use the environment variable G4INCLXX_NO_DE_EXCITATION to disable de-excitation
58 if(getenv("G4INCLXX_NO_DE_EXCITATION")) {
59 G4String message = "de-excitation is completely disabled!";
60 theInterfaceStore->EmitWarning(message);
61 theExcitationHandler = 0;
62 } else {
63 theExcitationHandler = new G4ExcitationHandler;
64 }
65
66 theBackupModel = new G4BinaryLightIonReaction;
67}
68
70{
71 delete theBackupModel;
72 delete theExcitationHandler;
73}
74
75G4bool G4INCLXXInterface::AccurateProjectile(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) const {
76 // Use direct kinematics if the projectile is a nucleon or a pion
77 const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
78 if(projectileDef == G4Proton::Proton()
79 || projectileDef == G4Neutron::Neutron()
80 || projectileDef == G4PionPlus::PionPlus()
81 || projectileDef == G4PionZero::PionZero()
82 || projectileDef == G4PionMinus::PionMinus())
83 return false;
84
85 // Here all projectiles should be nuclei
86 const G4int pA = projectileDef->GetAtomicMass();
87 if(pA<=0) {
88 std::stringstream ss;
89 ss << "the model does not know how to handle a collision between a "
90 << projectileDef->GetParticleName() << " projectile and a Z="
91 << theNucleus.GetZ_asInt() << ", A=" << theNucleus.GetA_asInt();
92 theInterfaceStore->EmitWarning(ss.str());
93 return true;
94 }
95
96 // If either nucleus is a LCP (A<=4), run the collision as light on heavy
97 const G4int tA = theNucleus.GetA_asInt();
98 if(tA<=4 || pA<=4) {
99 if(pA<tA)
100 return false;
101 else
102 return true;
103 }
104
105 // If one of the nuclei is heavier than theMaxProjMassINCL, run the collision
106 // as light on heavy.
107 // Note that here we are sure that either the projectile or the target is
108 // smaller than theMaxProjMass; otherwise theBackupModel would have been
109 // called.
110 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
111 if(pA > theMaxProjMassINCL)
112 return true;
113 else if(tA > theMaxProjMassINCL)
114 return false;
115 else
116 // In all other cases, use the global setting
117 return theInterfaceStore->GetAccurateProjectile();
118}
119
121{
122 // For systems heavier than theMaxProjMassINCL, use another model (typically
123 // BIC)
124 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
125 if(aTrack.GetDefinition()->GetAtomicMass() > theMaxProjMassINCL
126 && theNucleus.GetA_asInt() > theMaxProjMassINCL) {
127 if(!complainedAboutBackupModel) {
128 complainedAboutBackupModel = true;
129 std::stringstream ss;
130 ss << "INCL++ refuses to handle reactions between nuclei with A>"
131 << theMaxProjMassINCL
132 << ". A backup model ("
133 << theBackupModel->GetModelName()
134 << ") will be used instead.";
135 G4cout << "[INCL++] Warning: " << ss.str() << G4endl;
136 }
137 return theBackupModel->ApplyYourself(aTrack, theNucleus);
138 }
139
140 const G4int maxTries = 200;
141
142 // Check if inverse kinematics should be used
143 const G4bool inverseKinematics = AccurateProjectile(aTrack, theNucleus);
144
145 // If we are running in inverse kinematics, redefine aTrack and theNucleus
146 G4LorentzRotation *toInverseKinematics = NULL;
147 G4LorentzRotation *toDirectKinematics = NULL;
148 G4HadProjectile const *aProjectileTrack = &aTrack;
149 G4Nucleus *theTargetNucleus = &theNucleus;
150 if(inverseKinematics) {
151 G4ParticleTable * const theParticleTable = G4ParticleTable::GetParticleTable();
152 const G4int oldTargetA = theNucleus.GetA_asInt();
153 const G4int oldTargetZ = theNucleus.GetZ_asInt();
154 G4ParticleDefinition *oldTargetDef = theParticleTable->GetIon(oldTargetZ, oldTargetA, 0.0);
155 const G4ParticleDefinition *oldProjectileDef = aTrack.GetDefinition();
156
157 if(oldProjectileDef != 0 && oldTargetDef != 0) {
158 const G4int newTargetA = oldProjectileDef->GetAtomicMass();
159 const G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
160
161 if(newTargetA > 0 && newTargetZ > 0) {
162 // This should give us the same energy per nucleon
163 theTargetNucleus = new G4Nucleus(newTargetA, newTargetZ);
164 const G4double theProjectileMass = theParticleTable->GetIonTable()->GetIonMass(oldTargetZ, oldTargetA);
165 toInverseKinematics = new G4LorentzRotation(aTrack.Get4Momentum().boostVector());
166 G4LorentzVector theProjectile4Momentum(0.0, 0.0, 0.0, theProjectileMass);
167 G4DynamicParticle swappedProjectileParticle(oldTargetDef, (*toInverseKinematics) * theProjectile4Momentum);
168 aProjectileTrack = new G4HadProjectile(swappedProjectileParticle);
169 } else {
170 G4String message = "badly defined target after swapping. Falling back to normal (non-swapped) mode.";
171 theInterfaceStore->EmitWarning(message);
172 toInverseKinematics = new G4LorentzRotation;
173 }
174 } else {
175 G4String message = "oldProjectileDef or oldTargetDef was null";
176 theInterfaceStore->EmitWarning(message);
177 toInverseKinematics = new G4LorentzRotation;
178 }
179 }
180
181 // INCL assumes the projectile particle is going in the direction of
182 // the Z-axis. Here we construct proper rotation to convert the
183 // momentum vectors of the outcoming particles to the original
184 // coordinate system.
185 G4LorentzVector projectileMomentum = aProjectileTrack->Get4Momentum();
186
187 // INCL++ assumes that the projectile is going in the direction of
188 // the z-axis. In principle, if the coordinate system used by G4
189 // hadronic framework is defined differently we need a rotation to
190 // transform the INCL++ reaction products to the appropriate
191 // frame. Please note that it isn't necessary to apply this
192 // transform to the projectile because when creating the INCL++
193 // projectile particle; \see{toINCLKineticEnergy} needs to use only the
194 // projectile energy (direction is simply assumed to be along z-axis).
196 toZ.rotateZ(-projectileMomentum.phi());
197 toZ.rotateY(-projectileMomentum.theta());
198 G4RotationMatrix toLabFrame3 = toZ.inverse();
199 G4LorentzRotation toLabFrame(toLabFrame3);
200 // However, it turns out that the projectile given to us by G4
201 // hadronic framework is already going in the direction of the
202 // z-axis so this rotation is actually unnecessary. Both toZ and
203 // toLabFrame turn out to be unit matrices as can be seen by
204 // uncommenting the folowing two lines:
205 // G4cout <<"toZ = " << toZ << G4endl;
206 // G4cout <<"toLabFrame = " << toLabFrame << G4endl;
207
208 theResult.Clear(); // Make sure the output data structure is clean.
209 theResult.SetStatusChange(stopAndKill);
210
211 std::list<G4Fragment> remnants;
212
213 G4int nTries = 0;
214 // INCL can generate transparent events. However, this is meaningful
215 // only in the standalone code. In Geant4 we must "force" INCL to
216 // produce a valid cascade.
217 G4bool eventIsOK = false;
218 do {
219 const G4INCL::ParticleSpecies theSpecies = toINCLParticleSpecies(*aProjectileTrack);
220 const G4double kineticEnergy = toINCLKineticEnergy(*aProjectileTrack);
221
222 // The INCL model will be created at the first use
224
225 if(theInterfaceStore->GetDumpInput()) {
226 G4cout << theINCLModel->configToString() << G4endl;
227 }
228 const G4INCL::EventInfo eventInfo = theINCLModel->processEvent(theSpecies, kineticEnergy, theTargetNucleus->GetA_asInt(), theTargetNucleus->GetZ_asInt());
229 // eventIsOK = !eventInfo.transparent && nTries < maxTries;
230 eventIsOK = !eventInfo.transparent;
231 if(eventIsOK) {
232
233 // If the collision was run in inverse kinematics, prepare to boost back
234 // all the reaction products
235 if(inverseKinematics) {
236 toDirectKinematics = new G4LorentzRotation(toInverseKinematics->inverse());
237 }
238
239 for(G4int i = 0; i < eventInfo.nParticles; i++) {
240 G4int A = eventInfo.A[i];
241 G4int Z = eventInfo.Z[i];
242 // G4cout <<"INCL particle A = " << A << " Z = " << Z << G4endl;
243 G4double kinE = eventInfo.EKin[i];
244 G4double px = eventInfo.px[i];
245 G4double py = eventInfo.py[i];
246 G4double pz = eventInfo.pz[i];
247 G4DynamicParticle *p = toG4Particle(A, Z , kinE, px, py, pz);
248 if(p != 0) {
249 G4LorentzVector momentum = p->Get4Momentum();
250
251 // Apply the toLabFrame rotation
252 momentum *= toLabFrame;
253 // Apply the inverse-kinematics boost
254 if(inverseKinematics) {
255 momentum *= *toDirectKinematics;
256 momentum.setVect(-momentum.vect());
257 }
258
259 // Set the four-momentum of the reaction products
260 p->Set4Momentum(momentum);
261 theResult.AddSecondary(p);
262
263 } else {
264 G4String message = "the model produced a particle that couldn't be converted to Geant4 particle.";
265 theInterfaceStore->EmitWarning(message);
266 }
267 }
268
269 for(G4int i = 0; i < eventInfo.nRemnants; i++) {
270 const G4int A = eventInfo.ARem[i];
271 const G4int Z = eventInfo.ZRem[i];
272 // G4cout <<"INCL particle A = " << A << " Z = " << Z << G4endl;
273 const G4double kinE = eventInfo.EKinRem[i];
274 const G4double px = eventInfo.pxRem[i];
275 const G4double py = eventInfo.pyRem[i];
276 const G4double pz = eventInfo.pzRem[i];
277 G4ThreeVector spin(
278 eventInfo.jxRem[i]*hbar_Planck,
279 eventInfo.jyRem[i]*hbar_Planck,
280 eventInfo.jzRem[i]*hbar_Planck
281 );
282 const G4double excitationE = eventInfo.EStarRem[i];
283 const G4double nuclearMass = G4NucleiProperties::GetNuclearMass(A, Z) + excitationE;
284 const G4double scaling = remnant4MomentumScaling(nuclearMass,
285 kinE,
286 px, py, pz);
287 G4LorentzVector fourMomentum(scaling * px, scaling * py, scaling * pz,
288 nuclearMass + kinE);
289 if(std::abs(scaling - 1.0) > 0.01) {
290 std::stringstream ss;
291 ss << "momentum scaling = " << scaling
292 << "\n Lorentz vector = " << fourMomentum
293 << ")\n A = " << A << ", Z = " << Z
294 << "\n E* = " << excitationE << ", nuclearMass = " << nuclearMass
295 << "\n remnant i=" << i << ", nRemnants=" << eventInfo.nRemnants
296 << "\n Reaction was: " << aTrack.GetKineticEnergy()/MeV
297 << "-MeV " << aTrack.GetDefinition()->GetParticleName() << " + "
298 << G4ParticleTable::GetParticleTable()->GetIon(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0.0)->GetParticleName()
299 << ", in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
300 theInterfaceStore->EmitWarning(ss.str());
301 }
302
303 // Apply the toLabFrame rotation
304 fourMomentum *= toLabFrame;
305 spin *= toLabFrame3;
306 // Apply the inverse-kinematics boost
307 if(inverseKinematics) {
308 fourMomentum *= *toDirectKinematics;
309 fourMomentum.setVect(-fourMomentum.vect());
310 }
311
312 G4Fragment remnant(A, Z, fourMomentum);
313 remnant.SetAngularMomentum(spin);
314 remnants.push_back(remnant);
315 }
316 }
317 nTries++;
318 } while(!eventIsOK && nTries < maxTries);
319
320 // Clean up the objects that we created for the inverse kinematics
321 if(inverseKinematics) {
322 delete aProjectileTrack;
323 delete theTargetNucleus;
324 delete toInverseKinematics;
325 delete toDirectKinematics;
326 }
327
328 if(!eventIsOK) {
329 std::stringstream ss;
330 ss << "maximum number of tries exceeded for the proposed "
331 << aTrack.GetKineticEnergy()/MeV << "-MeV " << aTrack.GetDefinition()->GetParticleName()
332 << " + " << G4ParticleTable::GetParticleTable()->GetIon(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0.0)->GetParticleName()
333 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
334 theInterfaceStore->EmitWarning(ss.str());
335 theResult.SetStatusChange(isAlive);
336 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
337 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
338 return &theResult;
339 }
340
341 // De-excitation:
342
343 if(theExcitationHandler != 0) {
344 for(std::list<G4Fragment>::const_iterator i = remnants.begin();
345 i != remnants.end(); i++) {
346 G4ReactionProductVector *deExcitationResult = theExcitationHandler->BreakItUp((*i));
347
348 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
349 fragment != deExcitationResult->end(); ++fragment) {
350 G4ParticleDefinition *def = (*fragment)->GetDefinition();
351 if(def != 0) {
352 G4DynamicParticle *theFragment = new G4DynamicParticle(def, (*fragment)->GetMomentum());
353 theResult.AddSecondary(theFragment);
354 }
355 }
356
357 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
358 fragment != deExcitationResult->end(); ++fragment) {
359 delete (*fragment);
360 }
361 deExcitationResult->clear();
362 delete deExcitationResult;
363 }
364 }
365
366 remnants.clear();
367
368 return &theResult;
369}
370
372 return 0;
373}
374
375G4INCL::ParticleType G4INCLXXInterface::toINCLParticleType(G4ParticleDefinition const * const pdef) const {
376 if( pdef == G4Proton::Proton()) return G4INCL::Proton;
377 else if(pdef == G4Neutron::Neutron()) return G4INCL::Neutron;
378 else if(pdef == G4PionPlus::PionPlus()) return G4INCL::PiPlus;
379 else if(pdef == G4PionMinus::PionMinus()) return G4INCL::PiMinus;
380 else if(pdef == G4PionZero::PionZero()) return G4INCL::PiZero;
381 else if(pdef == G4Deuteron::Deuteron()) return G4INCL::Composite;
382 else if(pdef == G4Triton::Triton()) return G4INCL::Composite;
383 else if(pdef == G4He3::He3()) return G4INCL::Composite;
384 else if(pdef == G4Alpha::Alpha()) return G4INCL::Composite;
386 else return G4INCL::UnknownParticle;
387}
388
389G4INCL::ParticleSpecies G4INCLXXInterface::toINCLParticleSpecies(G4HadProjectile const &aTrack) const {
390 const G4ParticleDefinition *pdef = aTrack.GetDefinition();
391 const G4INCL::ParticleType theType = toINCLParticleType(pdef);
392 if(theType!=G4INCL::Composite)
393 return G4INCL::ParticleSpecies(theType);
394 else {
395 G4INCL::ParticleSpecies theSpecies;
396 theSpecies.theType=theType;
397 theSpecies.theA=(G4int) pdef->GetBaryonNumber();
398 theSpecies.theZ=(G4int) pdef->GetPDGCharge();
399 return theSpecies;
400 }
401}
402
403G4double G4INCLXXInterface::toINCLKineticEnergy(G4HadProjectile const &aTrack) const {
404 return aTrack.GetKineticEnergy();
405}
406
407G4ParticleDefinition *G4INCLXXInterface::toG4ParticleDefinition(G4int A, G4int Z) const {
408 if (A == 1 && Z == 1) return G4Proton::Proton();
409 else if(A == 1 && Z == 0) return G4Neutron::Neutron();
410 else if(A == 0 && Z == 1) return G4PionPlus::PionPlus();
411 else if(A == 0 && Z == -1) return G4PionMinus::PionMinus();
412 else if(A == 0 && Z == 0) return G4PionZero::PionZero();
413 else if(A == 2 && Z == 1) return G4Deuteron::Deuteron();
414 else if(A == 3 && Z == 1) return G4Triton::Triton();
415 else if(A == 3 && Z == 2) return G4He3::He3();
416 else if(A == 4 && Z == 2) return G4Alpha::Alpha();
417 else if(A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition
418 return G4ParticleTable::GetParticleTable()->GetIon(Z, A, 0.0);
419 } else { // Error, unrecognized particle
420 return 0;
421 }
422}
423
424G4DynamicParticle *G4INCLXXInterface::toG4Particle(G4int A, G4int Z,
425 G4double kinE,
426 G4double px,
427 G4double py, G4double pz) const {
428 const G4ParticleDefinition *def = toG4ParticleDefinition(A, Z);
429 if(def == 0) { // Check if we have a valid particle definition
430 return 0;
431 }
432 const G4double energy = kinE / MeV;
433 const G4ThreeVector momentum(px, py, pz);
434 const G4ThreeVector momentumDirection = momentum.unit();
435 G4DynamicParticle *p = new G4DynamicParticle(def, momentumDirection, energy);
436 return p;
437}
438
439G4double G4INCLXXInterface::remnant4MomentumScaling(G4double mass,
440 G4double kineticE,
441 G4double px, G4double py,
442 G4double pz) const {
443 const G4double p2 = px*px + py*py + pz*pz;
444 if(p2 > 0.0) {
445 const G4double pnew2 = kineticE*kineticE + 2.0*kineticE*mass;
446 return std::sqrt(pnew2)/std::sqrt(p2);
447 } else {
448 return 1.0;
449 }
450}
451
@ isAlive
@ stopAndKill
Header file for the G4INCLXXInterfaceStore class.
CLHEP::HepLorentzRotation G4LorentzRotation
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Hep3Vector unit() const
HepLorentzRotation inverse() const
double theta() const
Hep3Vector boostVector() const
Hep3Vector vect() const
void setVect(const Hep3Vector &)
HepRotation inverse() const
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState) const
void SetAngularMomentum(const G4ThreeVector &value)
Definition: G4Fragment.hh:267
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
const G4String & GetModelName() const
static G4He3 * He3()
Definition: G4He3.cc:94
Singleton class for configuring the INCL++ Geant4 interface.
void EmitWarning(const G4String &message)
Emit a warning to G4cout.
G4int GetMaxProjMassINCL() const
Getter for theMaxProjMassINCL.
G4INCL::INCL * GetINCLModel()
Get the cached INCL model engine.
static G4INCLXXInterfaceStore * GetInstance()
Get the singleton instance.
G4bool GetDumpInput() const
Getter for dumpInput.
G4bool GetAccurateProjectile() const
Getter for accurateProjectile.
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4INCLXXInterface(const G4String &name="INCL++ cascade with G4ExcitationHandler")
const EventInfo & processEvent()
std::string configToString()
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95
Bool_t transparent
True if the event is transparent.
Short_t A[maxSizeParticles]
Particle mass number.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [hbar].
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
Short_t Z[maxSizeParticles]
Particle charge number.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [hbar].
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Int_t nParticles
Total number of emitted particles.
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Int_t nRemnants
Number of remnants.
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [hbar].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.