Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
38#include <cmath>
39
41#include "G4SystemOfUnits.hh"
42#include "G4INCLXXInterface.hh"
43#include "G4GenericIon.hh"
44#include "G4INCLCascade.hh"
46#include "G4ReactionProduct.hh"
47#include "G4HadSecondary.hh"
48#include "G4ParticleTable.hh"
51#include "G4String.hh"
53#include "G4SystemOfUnits.hh"
55#include "G4INCLVersion.hh"
56#include "G4VEvaporation.hh"
61#include "G4INCLConfig.hh"
62#include "G4INCLRandom.hh"
63
65#include "G4HyperTriton.hh"
66#include "G4HyperH4.hh"
67#include "G4HyperAlpha.hh"
68#include "G4DoubleHyperH4.hh"
70#include "G4HyperHe5.hh"
71
73 G4VIntraNuclearTransportModel(G4INCLXXInterfaceStore::GetInstance()->getINCLXXVersionName()),
74 theINCLModel(NULL),
75 thePreCompoundModel(aPreCompound),
76 theInterfaceStore(G4INCLXXInterfaceStore::GetInstance()),
77 theTally(NULL),
78 complainedAboutBackupModel(false),
79 complainedAboutPreCompound(false),
80 theIonTable(G4IonTable::GetIonTable()),
81 theINCLXXLevelDensity(NULL),
82 theINCLXXFissionProbability(NULL),
83 secID(-1)
84{
85 if(!thePreCompoundModel) {
88 thePreCompoundModel = static_cast<G4VPreCompoundModel*>(p);
89 if(!thePreCompoundModel) { thePreCompoundModel = new G4PreCompoundModel; }
90 }
91
92 // Use the environment variable G4INCLXX_NO_DE_EXCITATION to disable de-excitation
93 if(std::getenv("G4INCLXX_NO_DE_EXCITATION")) {
94 G4String message = "de-excitation is completely disabled!";
95 theInterfaceStore->EmitWarning(message);
97 } else {
100 theDeExcitation = static_cast<G4VPreCompoundModel*>(p);
102
103 // set the fission parameters for G4ExcitationHandler
104 G4VEvaporationChannel * const theFissionChannel =
105 theDeExcitation->GetExcitationHandler()->GetEvaporation()->GetFissionChannel();
106 G4CompetitiveFission * const theFissionChannelCast = dynamic_cast<G4CompetitiveFission *>(theFissionChannel);
107 if(theFissionChannelCast) {
108 theINCLXXLevelDensity = new G4FissionLevelDensityParameterINCLXX;
109 theFissionChannelCast->SetLevelDensityParameter(theINCLXXLevelDensity);
110 theINCLXXFissionProbability = new G4FissionProbability;
111 theINCLXXFissionProbability->SetFissionLevelDensityParameter(theINCLXXLevelDensity);
112 theFissionChannelCast->SetEmissionStrategy(theINCLXXFissionProbability);
113 theInterfaceStore->EmitBigWarning("INCL++/G4ExcitationHandler uses its own level-density parameter for fission");
114 } else {
115 theInterfaceStore->EmitBigWarning("INCL++/G4ExcitationHandler could not use its own level-density parameter for fission");
116 }
117 }
118
119 // use the envvar G4INCLXX_DUMP_REMNANT to dump information about the
120 // remnants on stdout
121 if(std::getenv("G4INCLXX_DUMP_REMNANT"))
122 dumpRemnantInfo = true;
123 else
124 dumpRemnantInfo = false;
125
126 theBackupModel = new G4BinaryLightIonReaction;
127 theBackupModelNucleon = new G4BinaryCascade;
128 secID = G4PhysicsModelCatalog::GetModelID( "model_INCLXXCascade" );
129}
130
132{
133 delete theINCLXXLevelDensity;
134 delete theINCLXXFissionProbability;
135}
136
137G4bool G4INCLXXInterface::AccurateProjectile(const G4HadProjectile &aTrack, const G4Nucleus &theNucleus) const {
138 // Use direct kinematics if the projectile is a nucleon or a pion
139 const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
140 if(std::abs(projectileDef->GetBaryonNumber()) < 2) // Every non-composite particle. abs()-> in case of anti-nucleus projectile
141 return false;
142
143 // Here all projectiles should be nuclei
144 const G4int pA = projectileDef->GetAtomicMass();
145 if(pA<=0) {
146 std::stringstream ss;
147 ss << "the model does not know how to handle a collision between a "
148 << projectileDef->GetParticleName() << " projectile and a Z="
149 << theNucleus.GetZ_asInt() << ", A=" << theNucleus.GetA_asInt();
150 theInterfaceStore->EmitBigWarning(ss.str());
151 return true;
152 }
153
154 // If either nucleus is a LCP (A<=4), run the collision as light on heavy
155 const G4int tA = theNucleus.GetA_asInt();
156 if(tA<=4 || pA<=4) {
157 if(pA<tA)
158 return false;
159 else
160 return true;
161 }
162
163 // If one of the nuclei is heavier than theMaxProjMassINCL, run the collision
164 // as light on heavy.
165 // Note that here we are sure that either the projectile or the target is
166 // smaller than theMaxProjMass; otherwise theBackupModel would have been
167 // called.
168 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
169 if(pA > theMaxProjMassINCL)
170 return true;
171 else if(tA > theMaxProjMassINCL)
172 return false;
173 else
174 // In all other cases, use the global setting
175 return theInterfaceStore->GetAccurateProjectile();
176}
177
179{
180 G4ParticleDefinition const * const trackDefinition = aTrack.GetDefinition();
181 const G4bool isIonTrack = trackDefinition->GetParticleType()==G4GenericIon::GenericIon()->GetParticleType();
182 const G4int trackA = trackDefinition->GetAtomicMass();
183 const G4int trackZ = (G4int) trackDefinition->GetPDGCharge();
184 const G4int trackPDG = (G4int) trackDefinition->GetPDGEncoding();
185 const G4int trackL = trackDefinition->GetNumberOfLambdasInHypernucleus();
186 const G4int nucleusA = theNucleus.GetA_asInt();
187 const G4int nucleusZ = theNucleus.GetZ_asInt();
188
189 // For reactions induced by weird projectiles (e.g. He2), bail out
190 if((isIonTrack && ((trackZ<=0 && trackL==0) || trackA<=trackZ)) ||
191 (nucleusA>1 && (nucleusZ<=0 || nucleusA<=nucleusZ))) {
192 theResult.Clear();
193 theResult.SetStatusChange(isAlive);
194 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
195 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
196 return &theResult;
197 }
198
199 // For reactions on nucleons, use the backup model (without complaining),
200 // except for anti_proton and anti_neutron projectile (in this case, INCLXX is used).
201 if(trackA<=1 && nucleusA<=1 && (trackPDG!=-2212 && trackPDG!=-2112)) {
202 return theBackupModelNucleon->ApplyYourself(aTrack, theNucleus);
203 }
204
205 // For systems heavier than theMaxProjMassINCL, use another model (typically
206 // BIC)
207 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
208 if(trackA > theMaxProjMassINCL && nucleusA > theMaxProjMassINCL) {
209 if(!complainedAboutBackupModel) {
210 complainedAboutBackupModel = true;
211 std::stringstream ss;
212 ss << "INCL++ refuses to handle reactions between nuclei with A>"
213 << theMaxProjMassINCL
214 << ". A backup model ("
215 << theBackupModel->GetModelName()
216 << ") will be used instead.";
217 theInterfaceStore->EmitBigWarning(ss.str());
218 }
219 return theBackupModel->ApplyYourself(aTrack, theNucleus);
220 }
221
222 // For energies lower than cascadeMinEnergyPerNucleon, use PreCompound
223 const G4double cascadeMinEnergyPerNucleon = theInterfaceStore->GetCascadeMinEnergyPerNucleon();
224 const G4double trackKinE = aTrack.GetKineticEnergy();
225 if((trackDefinition==G4Neutron::NeutronDefinition() || trackDefinition==G4Proton::ProtonDefinition())
226 && trackKinE < cascadeMinEnergyPerNucleon) {
227 if(!complainedAboutPreCompound) {
228 complainedAboutPreCompound = true;
229 std::stringstream ss;
230 ss << "INCL++ refuses to handle nucleon-induced reactions below "
231 << cascadeMinEnergyPerNucleon / MeV
232 << " MeV. A PreCoumpound model ("
233 << thePreCompoundModel->GetModelName()
234 << ") will be used instead.";
235 theInterfaceStore->EmitBigWarning(ss.str());
236 }
237 return thePreCompoundModel->ApplyYourself(aTrack, theNucleus);
238 }
239
240 // Calculate the total four-momentum in the entrance channel
241 const G4double theNucleusMass = theIonTable->GetIonMass(nucleusZ, nucleusA);
242 const G4double theTrackMass = trackDefinition->GetPDGMass();
243 const G4double theTrackEnergy = trackKinE + theTrackMass;
244 const G4double theTrackMomentumAbs2 = theTrackEnergy*theTrackEnergy - theTrackMass*theTrackMass;
245 const G4double theTrackMomentumAbs = ((theTrackMomentumAbs2>0.0) ? std::sqrt(theTrackMomentumAbs2) : 0.0);
246 const G4ThreeVector theTrackMomentum = aTrack.Get4Momentum().getV().unit() * theTrackMomentumAbs;
247
248 G4LorentzVector goodTrack4Momentum(theTrackMomentum, theTrackEnergy);
249 G4LorentzVector fourMomentumIn;
250 fourMomentumIn.setE(theTrackEnergy + theNucleusMass);
251 fourMomentumIn.setVect(theTrackMomentum);
252
253 // Check if inverse kinematics should be used
254 const G4bool inverseKinematics = AccurateProjectile(aTrack, theNucleus);
255
256 // If we are running in inverse kinematics, redefine aTrack and theNucleus
257 G4LorentzRotation *toInverseKinematics = NULL;
258 G4LorentzRotation *toDirectKinematics = NULL;
259 G4HadProjectile const *aProjectileTrack = &aTrack;
260 G4Nucleus *theTargetNucleus = &theNucleus;
261 if(inverseKinematics) {
262 G4ParticleDefinition *oldTargetDef = theIonTable->GetIon(nucleusZ, nucleusA, 0);
263 const G4ParticleDefinition *oldProjectileDef = trackDefinition;
264
265 if(oldProjectileDef != 0 && oldTargetDef != 0) {
266 const G4int newTargetA = oldProjectileDef->GetAtomicMass();
267 const G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
268 const G4int newTargetL = oldProjectileDef->GetNumberOfLambdasInHypernucleus();
269 if(newTargetA > 0 && newTargetZ > 0) {
270 // This should give us the same energy per nucleon
271 theTargetNucleus = new G4Nucleus(newTargetA, newTargetZ, newTargetL);
272 toInverseKinematics = new G4LorentzRotation(goodTrack4Momentum.boostVector());
273 G4LorentzVector theProjectile4Momentum(0.0, 0.0, 0.0, theNucleusMass);
274 G4DynamicParticle swappedProjectileParticle(oldTargetDef, (*toInverseKinematics) * theProjectile4Momentum);
275 aProjectileTrack = new G4HadProjectile(swappedProjectileParticle);
276 } else {
277 G4String message = "badly defined target after swapping. Falling back to normal (non-swapped) mode.";
278 theInterfaceStore->EmitWarning(message);
279 toInverseKinematics = new G4LorentzRotation;
280 }
281 } else {
282 G4String message = "oldProjectileDef or oldTargetDef was null";
283 theInterfaceStore->EmitWarning(message);
284 toInverseKinematics = new G4LorentzRotation;
285 }
286 }
287
288 // INCL assumes the projectile particle is going in the direction of
289 // the Z-axis. Here we construct proper rotation to convert the
290 // momentum vectors of the outcoming particles to the original
291 // coordinate system.
292 G4LorentzVector projectileMomentum = aProjectileTrack->Get4Momentum();
293
294 // INCL++ assumes that the projectile is going in the direction of
295 // the z-axis. In principle, if the coordinate system used by G4
296 // hadronic framework is defined differently we need a rotation to
297 // transform the INCL++ reaction products to the appropriate
298 // frame. Please note that it isn't necessary to apply this
299 // transform to the projectile because when creating the INCL++
300 // projectile particle; \see{toINCLKineticEnergy} needs to use only the
301 // projectile energy (direction is simply assumed to be along z-axis).
303 toZ.rotateZ(-projectileMomentum.phi());
304 toZ.rotateY(-projectileMomentum.theta());
305 G4RotationMatrix toLabFrame3 = toZ.inverse();
306 G4LorentzRotation toLabFrame(toLabFrame3);
307 // However, it turns out that the projectile given to us by G4
308 // hadronic framework is already going in the direction of the
309 // z-axis so this rotation is actually unnecessary. Both toZ and
310 // toLabFrame turn out to be unit matrices as can be seen by
311 // uncommenting the folowing two lines:
312 // G4cout <<"toZ = " << toZ << G4endl;
313 // G4cout <<"toLabFrame = " << toLabFrame << G4endl;
314
315 theResult.Clear(); // Make sure the output data structure is clean.
316 theResult.SetStatusChange(stopAndKill);
317
318 std::list<G4Fragment> remnants;
319
320 const G4int maxTries = 200;
321 G4int nTries = 0;
322 // INCL can generate transparent events. However, this is meaningful
323 // only in the standalone code. In Geant4 we must "force" INCL to
324 // produce a valid cascade.
325 G4bool eventIsOK = false;
326 do {
327 const G4INCL::ParticleSpecies theSpecies = toINCLParticleSpecies(*aProjectileTrack);
328 const G4double kineticEnergy = toINCLKineticEnergy(*aProjectileTrack);
329
330 // The INCL model will be created at the first use
332
333 const G4INCL::EventInfo eventInfo = theINCLModel->processEvent(theSpecies, kineticEnergy,
334 theTargetNucleus->GetA_asInt(),
335 theTargetNucleus->GetZ_asInt(),
336 -theTargetNucleus->GetL()); // Strangeness has opposite sign
337 // eventIsOK = !eventInfo.transparent && nTries < maxTries; // of the number of Lambdas
338 eventIsOK = !eventInfo.transparent;
339 if(eventIsOK) {
340
341 // If the collision was run in inverse kinematics, prepare to boost back
342 // all the reaction products
343 if(inverseKinematics) {
344 toDirectKinematics = new G4LorentzRotation(toInverseKinematics->inverse());
345 }
346
347 G4LorentzVector fourMomentumOut;
348
349 for(G4int i = 0; i < eventInfo.nParticles; ++i) {
350 G4int A = eventInfo.A[i];
351 G4int Z = eventInfo.Z[i];
352 G4int S = eventInfo.S[i]; // Strangeness
353 G4int PDGCode = eventInfo.PDGCode[i];
354 G4double kinE = eventInfo.EKin[i];
355 G4double px = eventInfo.px[i];
356 G4double py = eventInfo.py[i];
357 G4double pz = eventInfo.pz[i];
358 G4DynamicParticle *p = toG4Particle(A, Z, S, PDGCode, kinE, px, py, pz);
359 if(p != 0) {
360 G4LorentzVector momentum = p->Get4Momentum();
361
362 // Apply the toLabFrame rotation
363 momentum *= toLabFrame;
364 // Apply the inverse-kinematics boost
365 if(inverseKinematics) {
366 momentum *= *toDirectKinematics;
367 momentum.setVect(-momentum.vect());
368 }
369
370 // Set the four-momentum of the reaction products
371 p->Set4Momentum(momentum);
372 fourMomentumOut += momentum;
373
374 // Propagate the particle's parent resonance information
375 G4HadSecondary secondary(p, 1.0, secID);
376 G4ParticleDefinition* parentResonanceDef = nullptr;
377 if ( eventInfo.parentResonancePDGCode[i] != 0 ) {
378 parentResonanceDef = G4ParticleTable::GetParticleTable()->FindParticle(eventInfo.parentResonancePDGCode[i]);
379 }
380 secondary.SetParentResonanceDef(parentResonanceDef);
381 secondary.SetParentResonanceID(eventInfo.parentResonanceID[i]);
382
383 theResult.AddSecondary(secondary);
384
385 } else {
386 G4String message = "the model produced a particle that couldn't be converted to Geant4 particle.";
387 theInterfaceStore->EmitWarning(message);
388 }
389 }
390
391 for(G4int i = 0; i < eventInfo.nRemnants; ++i) {
392 const G4int A = eventInfo.ARem[i];
393 const G4int Z = eventInfo.ZRem[i];
394 const G4int S = eventInfo.SRem[i];
395 // Check that the remnant is a physical bound state: if not, resample the collision.
396 if(( Z == 0 && S == 0 && A > 1 ) || // No bound states for nn, nnn, nnnn, ...
397 ( Z == 0 && S != 0 && A < 4 ) || // No bound states for nl, ll, nnl, nll, lll
398 ( Z != 0 && S != 0 && A == Z + std::abs(S) )) { // No bound states for pl, ppl, pll, ...
399 std::stringstream ss;
400 ss << "unphysical residual fragment : Z=" << Z << " S=" << S << " A=" << A
401 << " skipping it and resampling the collision";
402 theInterfaceStore->EmitWarning(ss.str());
403 eventIsOK = false;
404 continue;
405 }
406 const G4double kinE = eventInfo.EKinRem[i];
407 const G4double px = eventInfo.pxRem[i];
408 const G4double py = eventInfo.pyRem[i];
409 const G4double pz = eventInfo.pzRem[i];
410 G4ThreeVector spin(
411 eventInfo.jxRem[i]*hbar_Planck,
412 eventInfo.jyRem[i]*hbar_Planck,
413 eventInfo.jzRem[i]*hbar_Planck
414 );
415 const G4double excitationE = eventInfo.EStarRem[i];
416 G4double nuclearMass = excitationE;
417
418 if ( S == 0 ) {
419 nuclearMass += G4NucleiProperties::GetNuclearMass(A, Z);
420 } else {
421 // Assumed that the opposite of the strangeness of the remnant gives the number of Lambdas inside it
422 nuclearMass += G4HyperNucleiProperties::GetNuclearMass(A, Z, std::abs(S));
423 }
424 const G4double scaling = remnant4MomentumScaling(nuclearMass, kinE, px, py, pz);
425 G4LorentzVector fourMomentum(scaling * px, scaling * py, scaling * pz,
426 nuclearMass + kinE);
427 if(std::abs(scaling - 1.0) > 0.01) {
428 std::stringstream ss;
429 ss << "momentum scaling = " << scaling
430 << "\n Lorentz vector = " << fourMomentum
431 << ")\n A = " << A << ", Z = " << Z << ", S = " << S
432 << "\n E* = " << excitationE << ", nuclearMass = " << nuclearMass
433 << "\n remnant i=" << i << ", nRemnants=" << eventInfo.nRemnants
434 << "\n Reaction was: " << aTrack.GetKineticEnergy()/MeV
435 << "-MeV " << trackDefinition->GetParticleName() << " + "
436 << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
437 << ", in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
438 theInterfaceStore->EmitWarning(ss.str());
439 }
440
441 // Apply the toLabFrame rotation
442 fourMomentum *= toLabFrame;
443 spin *= toLabFrame3;
444 // Apply the inverse-kinematics boost
445 if(inverseKinematics) {
446 fourMomentum *= *toDirectKinematics;
447 fourMomentum.setVect(-fourMomentum.vect());
448 }
449
450 fourMomentumOut += fourMomentum;
451 G4Fragment remnant(A, Z, std::abs(S), fourMomentum); // Assumed that -strangeness gives the number of Lambdas
452 remnant.SetAngularMomentum(spin);
453 remnant.SetCreatorModelID(secID);
454 if(dumpRemnantInfo) {
455 G4cerr << "G4INCLXX_DUMP_REMNANT: " << remnant << " spin: " << spin << G4endl;
456 }
457 remnants.push_back(remnant);
458 }
459
460 // Give up is the event is not ok (e.g. unphysical residual)
461 if(!eventIsOK) {
462 const G4int nSecondaries = (G4int)theResult.GetNumberOfSecondaries();
463 for(G4int j=0; j<nSecondaries; ++j) delete theResult.GetSecondary(j)->GetParticle();
464 theResult.Clear();
465 theResult.SetStatusChange(stopAndKill);
466 remnants.clear();
467 } else {
468 /* // Check four-momentum conservation
469 G4INCL::Config *theConfig;
470 theConfig=&theInterfaceStore->GetINCLConfig();
471 G4double nbatrestThreshold=theConfig->getnbAtrestThreshold();
472 G4double pbatrestThreshold=theConfig->getAtrestThreshold();
473 if (((trackDefinition->GetParticleName() == "anti_neutron") && // in INCL antinucleon at rest is considered with 0 kinetic energy
474 (aTrack.GetKineticEnergy() <= nbatrestThreshold)) ||
475 ((trackDefinition->GetParticleName() == "anti_proton") &&
476 (aTrack.GetKineticEnergy() <= pbatrestThreshold)) ||
477 ((trackDefinition->GetParticleName() == "anti_neutron") &&
478 ((theTargetNucleus->GetA_asInt()==1 || theTargetNucleus->GetA_asInt()==2) && theTargetNucleus->GetZ_asInt()==1)) ||
479 ((trackDefinition->GetParticleName() == "anti_proton") &&
480 ((theTargetNucleus->GetA_asInt()==1 || theTargetNucleus->GetA_asInt()==2) && theTargetNucleus->GetZ_asInt()==1)))
481 {
482 fourMomentumIn.setE(theNucleusMass + theTrackMass);
483 fourMomentumIn.setVect(theTrackMomentum-theTrackMomentum);
484 }*/
485 const G4LorentzVector violation4Momentum = fourMomentumOut - fourMomentumIn;
486 const G4double energyViolation = std::abs(violation4Momentum.e());
487 const G4double momentumViolation = violation4Momentum.rho();
488 if(energyViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
489 std::stringstream ss;
490 ss << "energy conservation violated by " << energyViolation/MeV << " MeV in "
491 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
492 << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
493 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
494 theInterfaceStore->EmitWarning(ss.str());
495 eventIsOK = false;
496 const G4int nSecondaries = (G4int)theResult.GetNumberOfSecondaries();
497 for(G4int j=0; j<nSecondaries; ++j) delete theResult.GetSecondary(j)->GetParticle();
498 theResult.Clear();
499 theResult.SetStatusChange(stopAndKill);
500 remnants.clear();
501 } else if(momentumViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
502 std::stringstream ss;
503 ss << "momentum conservation violated by " << momentumViolation/MeV << " MeV in "
504 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
505 << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
506 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
507 theInterfaceStore->EmitWarning(ss.str());
508 eventIsOK = false;
509 const G4int nSecondaries = (G4int)theResult.GetNumberOfSecondaries();
510 for(G4int j=0; j<nSecondaries; ++j) delete theResult.GetSecondary(j)->GetParticle();
511 theResult.Clear();
512 theResult.SetStatusChange(stopAndKill);
513 remnants.clear();
514 }
515 }
516 }
517 nTries++;
518 } while(!eventIsOK && nTries < maxTries); /* Loop checking, 10.07.2015, D.Mancusi */
519
520 // Clean up the objects that we created for the inverse kinematics
521 if(inverseKinematics) {
522 delete aProjectileTrack;
523 delete theTargetNucleus;
524 delete toInverseKinematics;
525 delete toDirectKinematics;
526 }
527
528 if(!eventIsOK) {
529 std::stringstream ss;
530 ss << "maximum number of tries exceeded for the proposed "
531 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
532 << " + " << theIonTable->GetIonName(nucleusZ, nucleusA, 0)
533 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
534 theInterfaceStore->EmitWarning(ss.str());
535 theResult.SetStatusChange(isAlive);
536 theResult.SetEnergyChange(aTrack.GetKineticEnergy());
537 theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
538 return &theResult;
539 }
540
541 // De-excitation:
542
543 if(theDeExcitation != 0) {
544 for(std::list<G4Fragment>::iterator i = remnants.begin();
545 i != remnants.end(); ++i) {
546 G4ReactionProductVector *deExcitationResult = theDeExcitation->DeExcite((*i));
547
548 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
549 fragment != deExcitationResult->end(); ++fragment) {
550 const G4ParticleDefinition *def = (*fragment)->GetDefinition();
551 if(def != 0) {
552 G4DynamicParticle *theFragment = new G4DynamicParticle(def, (*fragment)->GetMomentum());
553 theResult.AddSecondary(theFragment, (*fragment)->GetCreatorModelID());
554 }
555 }
556
557 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
558 fragment != deExcitationResult->end(); ++fragment) {
559 delete (*fragment);
560 }
561 deExcitationResult->clear();
562 delete deExcitationResult;
563 }
564 }
565
566 remnants.clear();
567
568 if((theTally = theInterfaceStore->GetTally()))
569 theTally->Tally(aTrack, theNucleus, theResult);
570
571 return &theResult;
572}
573
577
578G4INCL::ParticleType G4INCLXXInterface::toINCLParticleType(G4ParticleDefinition const * const pdef) const {
579 if( pdef == G4Proton::Proton()) return G4INCL::Proton;
580 else if(pdef == G4Neutron::Neutron()) return G4INCL::Neutron;
581 else if(pdef == G4PionPlus::PionPlus()) return G4INCL::PiPlus;
582 else if(pdef == G4PionMinus::PionMinus()) return G4INCL::PiMinus;
583 else if(pdef == G4PionZero::PionZero()) return G4INCL::PiZero;
584 else if(pdef == G4KaonPlus::KaonPlus()) return G4INCL::KPlus;
585 else if(pdef == G4KaonZero::KaonZero()) return G4INCL::KZero;
586 else if(pdef == G4KaonMinus::KaonMinus()) return G4INCL::KMinus;
587 else if(pdef == G4AntiKaonZero::AntiKaonZero()) return G4INCL::KZeroBar;
588 // For K0L & K0S we do not take into account K0/K0B oscillations
589 else if(pdef == G4KaonZeroLong::KaonZeroLong()) return G4UniformRand() < 0.5 ? G4INCL::KZeroBar : G4INCL::KZero;
590 else if(pdef == G4KaonZeroShort::KaonZeroShort()) return G4UniformRand() < 0.5 ? G4INCL::KZeroBar : G4INCL::KZero;
591 else if(pdef == G4Deuteron::Deuteron()) return G4INCL::Composite;
592 else if(pdef == G4Triton::Triton()) return G4INCL::Composite;
593 else if(pdef == G4He3::He3()) return G4INCL::Composite;
594 else if(pdef == G4Alpha::Alpha()) return G4INCL::Composite;
595 else if(pdef == G4AntiProton::AntiProton()) return G4INCL::antiProton;
596 else if(pdef == G4AntiNeutron::AntiNeutron()) return G4INCL::antiNeutron;
597 else if(pdef->GetParticleType() == G4GenericIon::GenericIon()->GetParticleType()) return G4INCL::Composite;
598 else return G4INCL::UnknownParticle;
599}
600
601G4INCL::ParticleSpecies G4INCLXXInterface::toINCLParticleSpecies(G4HadProjectile const &aTrack) const {
602 const G4ParticleDefinition *pdef = aTrack.GetDefinition();
603 const G4INCL::ParticleType theType = toINCLParticleType(pdef);
604 if(theType!=G4INCL::Composite)
605 return G4INCL::ParticleSpecies(theType);
606 else {
607 G4INCL::ParticleSpecies theSpecies;
608 theSpecies.theType=theType;
609 theSpecies.theA=pdef->GetAtomicMass();
610 theSpecies.theZ=pdef->GetAtomicNumber();
611 return theSpecies;
612 }
613}
614
615G4double G4INCLXXInterface::toINCLKineticEnergy(G4HadProjectile const &aTrack) const {
616 return aTrack.GetKineticEnergy();
617}
618
619G4ParticleDefinition *G4INCLXXInterface::toG4ParticleDefinition(G4int A, G4int Z, G4int S, G4int PDGCode) const {
620 if (PDGCode == 2212) { return G4Proton::Proton();
621 } else if(PDGCode == 2112) { return G4Neutron::Neutron();
622 } else if(PDGCode == 211) { return G4PionPlus::PionPlus();
623 } else if(PDGCode == 111) { return G4PionZero::PionZero();
624 } else if(PDGCode == -211) { return G4PionMinus::PionMinus();
625
626 } else if(PDGCode == 221) { return G4Eta::Eta();
627 } else if(PDGCode == 22) { return G4Gamma::Gamma();
628
629 } else if(PDGCode == 3122) { return G4Lambda::Lambda();
630 } else if(PDGCode == 3222) { return G4SigmaPlus::SigmaPlus();
631 } else if(PDGCode == 3212) { return G4SigmaZero::SigmaZero();
632 } else if(PDGCode == 3112) { return G4SigmaMinus::SigmaMinus();
633 } else if(PDGCode == 321) { return G4KaonPlus::KaonPlus();
634 } else if(PDGCode == -321) { return G4KaonMinus::KaonMinus();
635 } else if(PDGCode == 130) { return G4KaonZeroLong::KaonZeroLong();
636 } else if(PDGCode == 310) { return G4KaonZeroShort::KaonZeroShort();
637 } else if(PDGCode == 311 || PDGCode == -311) {
639 else return G4KaonZeroLong::KaonZeroLong();
640
641 } else if(PDGCode == 1002) { return G4Deuteron::Deuteron();
642 } else if(PDGCode == 1003) { return G4Triton::Triton();
643 } else if(PDGCode == 2003) { return G4He3::He3();
644 } else if(PDGCode == 2004) { return G4Alpha::Alpha();
645
646 } else if(PDGCode == -2212) { return G4AntiProton::AntiProton();
647 } else if(PDGCode == -2112) { return G4AntiNeutron::AntiNeutron();
648 } else if(S != 0) { // Assumed that -S gives the number of Lambdas
649 if (A == 3 && Z == 1 && S == -1 ) return G4HyperTriton::Definition();
650 if (A == 4 && Z == 1 && S == -1 ) return G4HyperH4::Definition();
651 if (A == 4 && Z == 2 && S == -1 ) return G4HyperAlpha::Definition();
652 if (A == 4 && Z == 1 && S == -2 ) return G4DoubleHyperH4::Definition();
653 if (A == 4 && Z == 0 && S == -2 ) return G4DoubleHyperDoubleNeutron::Definition();
654 if (A == 5 && Z == 2 && S == -1 ) return G4HyperHe5::Definition();
655 } else if(A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition.
656 return theIonTable->GetIon(Z, A, 0);
657 }
658 return 0; // Error, unrecognized particle
659}
660
661G4DynamicParticle *G4INCLXXInterface::toG4Particle(G4int A, G4int Z, G4int S, G4int PDGCode,
662 G4double kinE, G4double px,
663 G4double py, G4double pz) const {
664 const G4ParticleDefinition *def = toG4ParticleDefinition(A, Z, S, PDGCode);
665 if(def == 0) { // Check if we have a valid particle definition
666 return 0;
667 }
668 const G4double energy = kinE * MeV;
669 const G4ThreeVector momentum(px, py, pz);
670 const G4ThreeVector momentumDirection = momentum.unit();
671 G4DynamicParticle *p = new G4DynamicParticle(def, momentumDirection, energy);
672 return p;
673}
674
675G4double G4INCLXXInterface::remnant4MomentumScaling(G4double mass,
676 G4double kineticE,
677 G4double px, G4double py,
678 G4double pz) const {
679 const G4double p2 = px*px + py*py + pz*pz;
680 if(p2 > 0.0) {
681 const G4double pnew2 = kineticE*kineticE + 2.0*kineticE*mass;
682 return std::sqrt(pnew2)/std::sqrt(p2);
683 } else {
684 return 1.0;
685 }
686}
687
688void G4INCLXXInterface::ModelDescription(std::ostream& outFile) const {
689 outFile
690 << "The Liège Intranuclear Cascade (INCL++) is a model for reactions induced\n"
691 << "by nucleons, antinucleons, pions, kaons, Lambda, Sigma and light ion on any nucleus.\n"
692 << "The reaction is described as an avalanche of binary nucleon-nucleon collisions,\n"
693 << "which can lead to the emission of energetic particles and to the formation of an\n"
694 << "excited thermalised nucleus (remnant). The de-excitation of the remnant is\n"
695 << "outside the scope of INCL++ and is typically described by another model.\n\n"
696 << "INCL++ has been reasonably well tested for nucleon (~50 MeV to ~15 GeV),\n"
697 << "pion (idem) and light-ion projectiles (up to A=18, ~10A MeV to 1A GeV).\n"
698 << "Most tests involved target nuclei close to the stability valley, with\n"
699 << "numbers between 4 and 300.\n\n"
700 << "Reference: D. Mancusi et al., Phys. Rev. C90 (2014) 054602\n\n";
701}
702
704 return theDeExcitation->GetModelName();
705}
G4double S(G4double temp)
@ isAlive
@ stopAndKill
Header file for the G4INCLXXInterfaceStore class.
CLHEP::HepLorentzRotation G4LorentzRotation
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::HepRotation G4RotationMatrix
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
HepLorentzRotation inverse() const
double theta() const
Hep3Vector boostVector() const
Hep3Vector getV() const
Hep3Vector vect() const
void setVect(const Hep3Vector &)
HepRotation inverse() const
HepRotation & rotateZ(double delta)
Definition Rotation.cc:87
HepRotation & rotateY(double delta)
Definition Rotation.cc:74
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4AntiKaonZero * AntiKaonZero()
static G4AntiNeutron * AntiNeutron()
static G4AntiProton * AntiProton()
void SetEmissionStrategy(G4VEmissionProbability *aFissionProb)
void SetLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4DoubleHyperDoubleNeutron * Definition()
static G4DoubleHyperH4 * Definition()
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Eta * Eta()
Definition G4Eta.cc:108
Revised level-density parameter for fission after INCL++.
void SetFissionLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
void SetCreatorModelID(G4int value)
void SetAngularMomentum(const G4ThreeVector &)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void SetParentResonanceDef(const G4ParticleDefinition *parentDef)
void SetParentResonanceID(const G4int parentID)
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4He3 * He3()
Definition G4He3.cc:90
static G4HyperAlpha * Definition()
static G4HyperH4 * Definition()
Definition G4HyperH4.cc:43
static G4HyperHe5 * Definition()
Definition G4HyperHe5.cc:43
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4HyperTriton * Definition()
Singleton class for configuring the INCL++ Geant4 interface.
static G4INCLXXInterfaceStore * GetInstance()
Get the singleton instance.
G4INCL::INCL * GetINCLModel()
Get the cached INCL model engine.
void EmitBigWarning(const G4String &message) const
Emit a BIG warning to G4cout.
virtual void ModelDescription(std::ostream &outFile) const
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
G4INCLXXInterface(G4VPreCompoundModel *const aPreCompound=0)
G4String const & GetDeExcitationModelName() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4KaonZero * KaonZero()
static G4Lambda * Lambda()
Definition G4Lambda.cc:105
static G4Neutron * NeutronDefinition()
Definition G4Neutron.cc:96
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
G4int GetL() const
Definition G4Nucleus.hh:87
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
G4int GetNumberOfLambdasInHypernucleus() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4PionZero * PionZero()
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
static G4SigmaZero * SigmaZero()
static G4Triton * Triton()
Definition G4Triton.cc:90
G4VIntraNuclearTransportModel(const G4String &mName="CascadeModel", G4VPreCompoundModel *ptr=nullptr)
G4double energy(const ThreeVector &p, const G4double m)
G4double shoot()
Short_t S[maxSizeParticles]
Particle strangeness number.
Int_t parentResonancePDGCode[maxSizeParticles]
Particle's parent resonance PDG code.
Short_t nRemnants
Number of remnants.
Bool_t transparent
True if the event is transparent.
Short_t A[maxSizeParticles]
Particle mass number.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Int_t parentResonanceID[maxSizeParticles]
Particle's parent resonance unique ID identifier.
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
Short_t Z[maxSizeParticles]
Particle charge number.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Short_t nParticles
Number of particles in the final state.
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Short_t SRem[maxSizeRemnants]
Remnant strangeness number.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Int_t PDGCode[maxSizeParticles]
PDG numbering of the particles.
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.