Geant4-11
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"
49#include "G4String.hh"
51#include "G4SystemOfUnits.hh"
53#include "G4INCLVersion.hh"
54#include "G4VEvaporation.hh"
59
61#include "G4HyperTriton.hh"
62#include "G4HyperH4.hh"
63#include "G4HyperAlpha.hh"
64#include "G4DoubleHyperH4.hh"
66#include "G4HyperHe5.hh"
67
70 theINCLModel(NULL),
71 thePreCompoundModel(aPreCompound),
72 theInterfaceStore(G4INCLXXInterfaceStore::GetInstance()),
73 theTally(NULL),
74 complainedAboutBackupModel(false),
75 complainedAboutPreCompound(false),
76 theIonTable(G4IonTable::GetIonTable()),
77 theINCLXXLevelDensity(NULL),
78 theINCLXXFissionProbability(NULL),
79 secID(-1)
80{
86 }
87
88 // Use the environment variable G4INCLXX_NO_DE_EXCITATION to disable de-excitation
89 if(std::getenv("G4INCLXX_NO_DE_EXCITATION")) {
90 G4String message = "de-excitation is completely disabled!";
93 } else {
96 theDeExcitation = static_cast<G4VPreCompoundModel*>(p);
98
99 // set the fission parameters for G4ExcitationHandler
100 G4VEvaporationChannel * const theFissionChannel =
102 G4CompetitiveFission * const theFissionChannelCast = dynamic_cast<G4CompetitiveFission *>(theFissionChannel);
103 if(theFissionChannelCast) {
105 theFissionChannelCast->SetLevelDensityParameter(theINCLXXLevelDensity);
108 theFissionChannelCast->SetEmissionStrategy(theINCLXXFissionProbability);
109 theInterfaceStore->EmitBigWarning("INCL++/G4ExcitationHandler uses its own level-density parameter for fission");
110 } else {
111 theInterfaceStore->EmitBigWarning("INCL++/G4ExcitationHandler could not use its own level-density parameter for fission");
112 }
113 }
114
115 // use the envvar G4INCLXX_DUMP_REMNANT to dump information about the
116 // remnants on stdout
117 if(std::getenv("G4INCLXX_DUMP_REMNANT"))
118 dumpRemnantInfo = true;
119 else
120 dumpRemnantInfo = false;
121
124 secID = G4PhysicsModelCatalog::GetModelID( "model_INCLXXCascade" );
125}
126
128{
131}
132
134 // Use direct kinematics if the projectile is a nucleon or a pion
135 const G4ParticleDefinition *projectileDef = aTrack.GetDefinition();
136 if(std::abs(projectileDef->GetBaryonNumber()) < 2) // Every non-composite particle. abs()-> in case of anti-nucleus projectile
137 return false;
138
139 // Here all projectiles should be nuclei
140 const G4int pA = projectileDef->GetAtomicMass();
141 if(pA<=0) {
142 std::stringstream ss;
143 ss << "the model does not know how to handle a collision between a "
144 << projectileDef->GetParticleName() << " projectile and a Z="
145 << theNucleus.GetZ_asInt() << ", A=" << theNucleus.GetA_asInt();
147 return true;
148 }
149
150 // If either nucleus is a LCP (A<=4), run the collision as light on heavy
151 const G4int tA = theNucleus.GetA_asInt();
152 if(tA<=4 || pA<=4) {
153 if(pA<tA)
154 return false;
155 else
156 return true;
157 }
158
159 // If one of the nuclei is heavier than theMaxProjMassINCL, run the collision
160 // as light on heavy.
161 // Note that here we are sure that either the projectile or the target is
162 // smaller than theMaxProjMass; otherwise theBackupModel would have been
163 // called.
164 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
165 if(pA > theMaxProjMassINCL)
166 return true;
167 else if(tA > theMaxProjMassINCL)
168 return false;
169 else
170 // In all other cases, use the global setting
172}
173
175{
176 G4ParticleDefinition const * const trackDefinition = aTrack.GetDefinition();
177 const G4bool isIonTrack = trackDefinition->GetParticleType()==G4GenericIon::GenericIon()->GetParticleType();
178 const G4int trackA = trackDefinition->GetAtomicMass();
179 const G4int trackZ = (G4int) trackDefinition->GetPDGCharge();
180 const G4int trackL = trackDefinition->GetNumberOfLambdasInHypernucleus();
181 const G4int nucleusA = theNucleus.GetA_asInt();
182 const G4int nucleusZ = theNucleus.GetZ_asInt();
183
184 // For reactions induced by weird projectiles (e.g. He2), bail out
185 if((isIonTrack && ((trackZ<=0 && trackL==0) || trackA<=trackZ)) ||
186 (nucleusA>1 && (nucleusZ<=0 || nucleusA<=nucleusZ))) {
191 return &theResult;
192 }
193
194 // For reactions on nucleons, use the backup model (without complaining)
195 if(trackA<=1 && nucleusA<=1) {
196 return theBackupModelNucleon->ApplyYourself(aTrack, theNucleus);
197 }
198
199 // For systems heavier than theMaxProjMassINCL, use another model (typically
200 // BIC)
201 const G4int theMaxProjMassINCL = theInterfaceStore->GetMaxProjMassINCL();
202 if(trackA > theMaxProjMassINCL && nucleusA > theMaxProjMassINCL) {
205 std::stringstream ss;
206 ss << "INCL++ refuses to handle reactions between nuclei with A>"
207 << theMaxProjMassINCL
208 << ". A backup model ("
210 << ") will be used instead.";
212 }
213 return theBackupModel->ApplyYourself(aTrack, theNucleus);
214 }
215
216 // For energies lower than cascadeMinEnergyPerNucleon, use PreCompound
217 const G4double cascadeMinEnergyPerNucleon = theInterfaceStore->GetCascadeMinEnergyPerNucleon();
218 const G4double trackKinE = aTrack.GetKineticEnergy();
219 if((trackDefinition==G4Neutron::NeutronDefinition() || trackDefinition==G4Proton::ProtonDefinition())
220 && trackKinE < cascadeMinEnergyPerNucleon) {
223 std::stringstream ss;
224 ss << "INCL++ refuses to handle nucleon-induced reactions below "
225 << cascadeMinEnergyPerNucleon / MeV
226 << " MeV. A PreCoumpound model ("
228 << ") will be used instead.";
230 }
231 return thePreCompoundModel->ApplyYourself(aTrack, theNucleus);
232 }
233
234 // Calculate the total four-momentum in the entrance channel
235 const G4double theNucleusMass = theIonTable->GetIonMass(nucleusZ, nucleusA);
236 const G4double theTrackMass = trackDefinition->GetPDGMass();
237 const G4double theTrackEnergy = trackKinE + theTrackMass;
238 const G4double theTrackMomentumAbs2 = theTrackEnergy*theTrackEnergy - theTrackMass*theTrackMass;
239 const G4double theTrackMomentumAbs = ((theTrackMomentumAbs2>0.0) ? std::sqrt(theTrackMomentumAbs2) : 0.0);
240 const G4ThreeVector theTrackMomentum = aTrack.Get4Momentum().getV().unit() * theTrackMomentumAbs;
241
242 G4LorentzVector goodTrack4Momentum(theTrackMomentum, theTrackEnergy);
243 G4LorentzVector fourMomentumIn;
244 fourMomentumIn.setE(theTrackEnergy + theNucleusMass);
245 fourMomentumIn.setVect(theTrackMomentum);
246
247 // Check if inverse kinematics should be used
248 const G4bool inverseKinematics = AccurateProjectile(aTrack, theNucleus);
249
250 // If we are running in inverse kinematics, redefine aTrack and theNucleus
251 G4LorentzRotation *toInverseKinematics = NULL;
252 G4LorentzRotation *toDirectKinematics = NULL;
253 G4HadProjectile const *aProjectileTrack = &aTrack;
254 G4Nucleus *theTargetNucleus = &theNucleus;
255 if(inverseKinematics) {
256 G4ParticleDefinition *oldTargetDef = theIonTable->GetIon(nucleusZ, nucleusA, 0);
257 const G4ParticleDefinition *oldProjectileDef = trackDefinition;
258
259 if(oldProjectileDef != 0 && oldTargetDef != 0) {
260 const G4int newTargetA = oldProjectileDef->GetAtomicMass();
261 const G4int newTargetZ = oldProjectileDef->GetAtomicNumber();
262 const G4int newTargetL = oldProjectileDef->GetNumberOfLambdasInHypernucleus();
263 if(newTargetA > 0 && newTargetZ > 0) {
264 // This should give us the same energy per nucleon
265 theTargetNucleus = new G4Nucleus(newTargetA, newTargetZ, newTargetL);
266 toInverseKinematics = new G4LorentzRotation(goodTrack4Momentum.boostVector());
267 G4LorentzVector theProjectile4Momentum(0.0, 0.0, 0.0, theNucleusMass);
268 G4DynamicParticle swappedProjectileParticle(oldTargetDef, (*toInverseKinematics) * theProjectile4Momentum);
269 aProjectileTrack = new G4HadProjectile(swappedProjectileParticle);
270 } else {
271 G4String message = "badly defined target after swapping. Falling back to normal (non-swapped) mode.";
273 toInverseKinematics = new G4LorentzRotation;
274 }
275 } else {
276 G4String message = "oldProjectileDef or oldTargetDef was null";
278 toInverseKinematics = new G4LorentzRotation;
279 }
280 }
281
282 // INCL assumes the projectile particle is going in the direction of
283 // the Z-axis. Here we construct proper rotation to convert the
284 // momentum vectors of the outcoming particles to the original
285 // coordinate system.
286 G4LorentzVector projectileMomentum = aProjectileTrack->Get4Momentum();
287
288 // INCL++ assumes that the projectile is going in the direction of
289 // the z-axis. In principle, if the coordinate system used by G4
290 // hadronic framework is defined differently we need a rotation to
291 // transform the INCL++ reaction products to the appropriate
292 // frame. Please note that it isn't necessary to apply this
293 // transform to the projectile because when creating the INCL++
294 // projectile particle; \see{toINCLKineticEnergy} needs to use only the
295 // projectile energy (direction is simply assumed to be along z-axis).
297 toZ.rotateZ(-projectileMomentum.phi());
298 toZ.rotateY(-projectileMomentum.theta());
299 G4RotationMatrix toLabFrame3 = toZ.inverse();
300 G4LorentzRotation toLabFrame(toLabFrame3);
301 // However, it turns out that the projectile given to us by G4
302 // hadronic framework is already going in the direction of the
303 // z-axis so this rotation is actually unnecessary. Both toZ and
304 // toLabFrame turn out to be unit matrices as can be seen by
305 // uncommenting the folowing two lines:
306 // G4cout <<"toZ = " << toZ << G4endl;
307 // G4cout <<"toLabFrame = " << toLabFrame << G4endl;
308
309 theResult.Clear(); // Make sure the output data structure is clean.
311
312 std::list<G4Fragment> remnants;
313
314 const G4int maxTries = 200;
315 G4int nTries = 0;
316 // INCL can generate transparent events. However, this is meaningful
317 // only in the standalone code. In Geant4 we must "force" INCL to
318 // produce a valid cascade.
319 G4bool eventIsOK = false;
320 do {
321 const G4INCL::ParticleSpecies theSpecies = toINCLParticleSpecies(*aProjectileTrack);
322 const G4double kineticEnergy = toINCLKineticEnergy(*aProjectileTrack);
323
324 // The INCL model will be created at the first use
326
327 const G4INCL::EventInfo eventInfo = theINCLModel->processEvent(theSpecies, kineticEnergy,
328 theTargetNucleus->GetA_asInt(),
329 theTargetNucleus->GetZ_asInt(),
330 -theTargetNucleus->GetL()); // Strangeness has opposite sign
331 // eventIsOK = !eventInfo.transparent && nTries < maxTries; // of the number of Lambdas
332 eventIsOK = !eventInfo.transparent;
333 if(eventIsOK) {
334
335 // If the collision was run in inverse kinematics, prepare to boost back
336 // all the reaction products
337 if(inverseKinematics) {
338 toDirectKinematics = new G4LorentzRotation(toInverseKinematics->inverse());
339 }
340
341 G4LorentzVector fourMomentumOut;
342
343 for(G4int i = 0; i < eventInfo.nParticles; ++i) {
344 G4int A = eventInfo.A[i];
345 G4int Z = eventInfo.Z[i];
346 G4int S = eventInfo.S[i]; // Strangeness
347 G4int PDGCode = eventInfo.PDGCode[i];
348 // G4cout <<"INCL particle A = " << A << " Z = " << Z << " S = " << S << G4endl;
349 G4double kinE = eventInfo.EKin[i];
350 G4double px = eventInfo.px[i];
351 G4double py = eventInfo.py[i];
352 G4double pz = eventInfo.pz[i];
353 G4DynamicParticle *p = toG4Particle(A, Z, S, PDGCode, kinE, px, py, pz);
354 if(p != 0) {
355 G4LorentzVector momentum = p->Get4Momentum();
356
357 // Apply the toLabFrame rotation
358 momentum *= toLabFrame;
359 // Apply the inverse-kinematics boost
360 if(inverseKinematics) {
361 momentum *= *toDirectKinematics;
362 momentum.setVect(-momentum.vect());
363 }
364
365 // Set the four-momentum of the reaction products
366 p->Set4Momentum(momentum);
367 fourMomentumOut += momentum;
369
370 } else {
371 G4String message = "the model produced a particle that couldn't be converted to Geant4 particle.";
373 }
374 }
375
376 for(G4int i = 0; i < eventInfo.nRemnants; ++i) {
377 const G4int A = eventInfo.ARem[i];
378 const G4int Z = eventInfo.ZRem[i];
379 const G4int S = eventInfo.SRem[i];
380 // G4cout <<"INCL particle A = " << A << " Z = " << Z << " S= " << S << G4endl;
381 const G4double kinE = eventInfo.EKinRem[i];
382 const G4double px = eventInfo.pxRem[i];
383 const G4double py = eventInfo.pyRem[i];
384 const G4double pz = eventInfo.pzRem[i];
385 G4ThreeVector spin(
386 eventInfo.jxRem[i]*hbar_Planck,
387 eventInfo.jyRem[i]*hbar_Planck,
388 eventInfo.jzRem[i]*hbar_Planck
389 );
390 const G4double excitationE = eventInfo.EStarRem[i];
391 G4double nuclearMass = excitationE;
392 if ( S == 0 ) {
393 nuclearMass += G4NucleiProperties::GetNuclearMass(A, Z);
394 } else {
395 // Assumed that the opposite of the strangeness of the remnant gives the number of Lambdas inside it
396 nuclearMass += G4HyperNucleiProperties::GetNuclearMass(A, Z, std::abs(S));
397 }
398 const G4double scaling = remnant4MomentumScaling(nuclearMass, kinE, px, py, pz);
399 G4LorentzVector fourMomentum(scaling * px, scaling * py, scaling * pz,
400 nuclearMass + kinE);
401 if(std::abs(scaling - 1.0) > 0.01) {
402 std::stringstream ss;
403 ss << "momentum scaling = " << scaling
404 << "\n Lorentz vector = " << fourMomentum
405 << ")\n A = " << A << ", Z = " << Z << ", S = " << S
406 << "\n E* = " << excitationE << ", nuclearMass = " << nuclearMass
407 << "\n remnant i=" << i << ", nRemnants=" << eventInfo.nRemnants
408 << "\n Reaction was: " << aTrack.GetKineticEnergy()/MeV
409 << "-MeV " << trackDefinition->GetParticleName() << " + "
410 << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
411 << ", in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
413 }
414
415 // Apply the toLabFrame rotation
416 fourMomentum *= toLabFrame;
417 spin *= toLabFrame3;
418 // Apply the inverse-kinematics boost
419 if(inverseKinematics) {
420 fourMomentum *= *toDirectKinematics;
421 fourMomentum.setVect(-fourMomentum.vect());
422 }
423
424 fourMomentumOut += fourMomentum;
425 G4Fragment remnant(A, Z, std::abs(S), fourMomentum); // Assumed that -strangeness gives the number of Lambdas
426 remnant.SetAngularMomentum(spin);
427 remnant.SetCreatorModelID(secID);
428 if(dumpRemnantInfo) {
429 G4cerr << "G4INCLXX_DUMP_REMNANT: " << remnant << " spin: " << spin << G4endl;
430 }
431 remnants.push_back(remnant);
432 }
433
434 // Check four-momentum conservation
435 const G4LorentzVector violation4Momentum = fourMomentumOut - fourMomentumIn;
436 const G4double energyViolation = std::abs(violation4Momentum.e());
437 const G4double momentumViolation = violation4Momentum.rho();
438 if(energyViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
439 std::stringstream ss;
440 ss << "energy conservation violated by " << energyViolation/MeV << " MeV in "
441 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
442 << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
443 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
445 eventIsOK = false;
446 const G4int nSecondaries = theResult.GetNumberOfSecondaries();
447 for(G4int j=0; j<nSecondaries; ++j)
451 remnants.clear();
452 } else if(momentumViolation > G4INCLXXInterfaceStore::GetInstance()->GetConservationTolerance()) {
453 std::stringstream ss;
454 ss << "momentum conservation violated by " << momentumViolation/MeV << " MeV in "
455 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
456 << " + " << theIonTable->GetIonName(theNucleus.GetZ_asInt(), theNucleus.GetA_asInt(), 0)
457 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics. Will resample.";
459 eventIsOK = false;
460 const G4int nSecondaries = theResult.GetNumberOfSecondaries();
461 for(G4int j=0; j<nSecondaries; ++j)
465 remnants.clear();
466 }
467 }
468 nTries++;
469 } while(!eventIsOK && nTries < maxTries); /* Loop checking, 10.07.2015, D.Mancusi */
470
471 // Clean up the objects that we created for the inverse kinematics
472 if(inverseKinematics) {
473 delete aProjectileTrack;
474 delete theTargetNucleus;
475 delete toInverseKinematics;
476 delete toDirectKinematics;
477 }
478
479 if(!eventIsOK) {
480 std::stringstream ss;
481 ss << "maximum number of tries exceeded for the proposed "
482 << aTrack.GetKineticEnergy()/MeV << "-MeV " << trackDefinition->GetParticleName()
483 << " + " << theIonTable->GetIonName(nucleusZ, nucleusA, 0)
484 << " inelastic reaction, in " << (inverseKinematics ? "inverse" : "direct") << " kinematics.";
489 return &theResult;
490 }
491
492 // De-excitation:
493
494 if(theDeExcitation != 0) {
495 for(std::list<G4Fragment>::iterator i = remnants.begin();
496 i != remnants.end(); ++i) {
497 G4ReactionProductVector *deExcitationResult = theDeExcitation->DeExcite((*i));
498
499 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
500 fragment != deExcitationResult->end(); ++fragment) {
501 const G4ParticleDefinition *def = (*fragment)->GetDefinition();
502 if(def != 0) {
503 G4DynamicParticle *theFragment = new G4DynamicParticle(def, (*fragment)->GetMomentum());
504 theResult.AddSecondary(theFragment, (*fragment)->GetCreatorModelID());
505 }
506 }
507
508 for(G4ReactionProductVector::iterator fragment = deExcitationResult->begin();
509 fragment != deExcitationResult->end(); ++fragment) {
510 delete (*fragment);
511 }
512 deExcitationResult->clear();
513 delete deExcitationResult;
514 }
515 }
516
517 remnants.clear();
518
520 theTally->Tally(aTrack, theNucleus, theResult);
521
522 return &theResult;
523}
524
526 return 0;
527}
528
530 if( pdef == G4Proton::Proton()) return G4INCL::Proton;
531 else if(pdef == G4Neutron::Neutron()) return G4INCL::Neutron;
532 else if(pdef == G4PionPlus::PionPlus()) return G4INCL::PiPlus;
533 else if(pdef == G4PionMinus::PionMinus()) return G4INCL::PiMinus;
534 else if(pdef == G4PionZero::PionZero()) return G4INCL::PiZero;
535 else if(pdef == G4KaonPlus::KaonPlus()) return G4INCL::KPlus;
536 else if(pdef == G4KaonZero::KaonZero()) return G4INCL::KZero;
537 else if(pdef == G4KaonMinus::KaonMinus()) return G4INCL::KMinus;
538 else if(pdef == G4AntiKaonZero::AntiKaonZero()) return G4INCL::KZeroBar;
539 // For K0L & K0S we do not take into account K0/K0B oscillations
540 else if(pdef == G4KaonZeroLong::KaonZeroLong()) return G4UniformRand() < 0.5 ? G4INCL::KZeroBar : G4INCL::KZero;
541 else if(pdef == G4KaonZeroShort::KaonZeroShort()) return G4UniformRand() < 0.5 ? G4INCL::KZeroBar : G4INCL::KZero;
542 else if(pdef == G4Deuteron::Deuteron()) return G4INCL::Composite;
543 else if(pdef == G4Triton::Triton()) return G4INCL::Composite;
544 else if(pdef == G4He3::He3()) return G4INCL::Composite;
545 else if(pdef == G4Alpha::Alpha()) return G4INCL::Composite;
547 else return G4INCL::UnknownParticle;
548}
549
551 const G4ParticleDefinition *pdef = aTrack.GetDefinition();
552 const G4INCL::ParticleType theType = toINCLParticleType(pdef);
553 if(theType!=G4INCL::Composite)
554 return G4INCL::ParticleSpecies(theType);
555 else {
556 G4INCL::ParticleSpecies theSpecies;
557 theSpecies.theType=theType;
558 theSpecies.theA=pdef->GetAtomicMass();
559 theSpecies.theZ=pdef->GetAtomicNumber();
560 return theSpecies;
561 }
562}
563
565 return aTrack.GetKineticEnergy();
566}
567
569 if (PDGCode == 2212) { return G4Proton::Proton();
570 } else if(PDGCode == 2112) { return G4Neutron::Neutron();
571 } else if(PDGCode == 211) { return G4PionPlus::PionPlus();
572 } else if(PDGCode == 111) { return G4PionZero::PionZero();
573 } else if(PDGCode == -211) { return G4PionMinus::PionMinus();
574
575 } else if(PDGCode == 221) { return G4Eta::Eta();
576 } else if(PDGCode == 22) { return G4Gamma::Gamma();
577
578 } else if(PDGCode == 3122) { return G4Lambda::Lambda();
579 } else if(PDGCode == 3222) { return G4SigmaPlus::SigmaPlus();
580 } else if(PDGCode == 3212) { return G4SigmaZero::SigmaZero();
581 } else if(PDGCode == 3112) { return G4SigmaMinus::SigmaMinus();
582 } else if(PDGCode == 321) { return G4KaonPlus::KaonPlus();
583 } else if(PDGCode == -321) { return G4KaonMinus::KaonMinus();
584 } else if(PDGCode == 130) { return G4KaonZeroLong::KaonZeroLong();
585 } else if(PDGCode == 310) { return G4KaonZeroShort::KaonZeroShort();
586
587 } else if(PDGCode == 1002) { return G4Deuteron::Deuteron();
588 } else if(PDGCode == 1003) { return G4Triton::Triton();
589 } else if(PDGCode == 2003) { return G4He3::He3();
590 } else if(PDGCode == 2004) { return G4Alpha::Alpha();
591 } else if(S != 0) { // Assumed that -S gives the number of Lambdas
592 if (A == 3 && Z == 1 && S == -1 ) return G4HyperTriton::Definition();
593 if (A == 4 && Z == 1 && S == -1 ) return G4HyperH4::Definition();
594 if (A == 4 && Z == 2 && S == -1 ) return G4HyperAlpha::Definition();
595 if (A == 4 && Z == 1 && S == -2 ) return G4DoubleHyperH4::Definition();
596 if (A == 4 && Z == 0 && S == -2 ) return G4DoubleHyperDoubleNeutron::Definition();
597 if (A == 5 && Z == 2 && S == -1 ) return G4HyperHe5::Definition();
598 } else if(A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition.
599 return theIonTable->GetIon(Z, A, 0);
600 }
601 return 0; // Error, unrecognized particle
602}
603
605 G4double kinE, G4double px,
606 G4double py, G4double pz) const {
607 const G4ParticleDefinition *def = toG4ParticleDefinition(A, Z, S, PDGCode);
608 if(def == 0) { // Check if we have a valid particle definition
609 return 0;
610 }
611 const G4double energy = kinE * MeV;
612 const G4ThreeVector momentum(px, py, pz);
613 const G4ThreeVector momentumDirection = momentum.unit();
614 G4DynamicParticle *p = new G4DynamicParticle(def, momentumDirection, energy);
615 return p;
616}
617
619 G4double kineticE,
620 G4double px, G4double py,
621 G4double pz) const {
622 const G4double p2 = px*px + py*py + pz*pz;
623 if(p2 > 0.0) {
624 const G4double pnew2 = kineticE*kineticE + 2.0*kineticE*mass;
625 return std::sqrt(pnew2)/std::sqrt(p2);
626 } else {
627 return 1.0;
628 }
629}
630
631void G4INCLXXInterface::ModelDescription(std::ostream& outFile) const {
632 outFile
633 << "The Li�ge Intranuclear Cascade (INCL++) is a model for reactions induced\n"
634 << "by nucleons, pions and light ion on any nucleus. The reaction is\n"
635 << "described as an avalanche of binary nucleon-nucleon collisions, which can\n"
636 << "lead to the emission of energetic particles and to the formation of an\n"
637 << "excited thermalised nucleus (remnant). The de-excitation of the remnant is\n"
638 << "outside the scope of INCL++ and is typically described by another model.\n\n"
639 << "INCL++ has been reasonably well tested for nucleon (~50 MeV to ~15 GeV),\n"
640 << "pion (idem) and light-ion projectiles (up to A=18, ~10A MeV to 1A GeV).\n"
641 << "Most tests involved target nuclei close to the stability valley, with\n"
642 << "numbers between 4 and 250.\n\n"
643 << "Reference: D. Mancusi et al., Phys. Rev. C90 (2014) 054602\n\n";
644}
645
648}
G4double S(G4double temp)
@ isAlive
@ stopAndKill
Header file for the G4INCLXXInterfaceStore class.
CLHEP::HepLorentzRotation G4LorentzRotation
std::vector< G4ReactionProduct * > G4ReactionProductVector
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
#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:88
static G4AntiKaonZero * AntiKaonZero()
void SetEmissionStrategy(G4VEmissionProbability *aFissionProb)
void SetLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4DoubleHyperDoubleNeutron * Definition()
static G4DoubleHyperH4 * Definition()
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Eta * Eta()
Definition: G4Eta.cc:108
G4VEvaporation * GetEvaporation()
Revised level-density parameter for fission after INCL++.
void SetFissionLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
void SetAngularMomentum(const G4ThreeVector &)
Definition: G4Fragment.cc:281
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
void SetEnergyChange(G4double anEnergy)
G4HadSecondary * GetSecondary(size_t i)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4DynamicParticle * GetParticle()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
static G4He3 * He3()
Definition: G4He3.cc:93
static G4HyperAlpha * Definition()
Definition: G4HyperAlpha.cc:45
static G4HyperH4 * Definition()
Definition: G4HyperH4.cc:45
static G4HyperHe5 * Definition()
Definition: G4HyperHe5.cc:45
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4HyperTriton * Definition()
Singleton class for configuring the INCL++ Geant4 interface.
void EmitWarning(const G4String &message)
Emit a warning to G4cout.
static G4INCLXXInterfaceStore * GetInstance()
Get the singleton instance.
G4int GetMaxProjMassINCL() const
Getter for theMaxProjMassINCL.
G4INCL::INCL * GetINCLModel()
Get the cached INCL model engine.
G4double GetConservationTolerance() const
Getter for conservationTolerance.
void EmitBigWarning(const G4String &message) const
Emit a BIG warning to G4cout.
G4double GetCascadeMinEnergyPerNucleon() const
Getter for cascadeMinEnergyPerNucleon.
G4INCLXXVInterfaceTally * GetTally() const
Getter for the interface tally.
G4bool GetAccurateProjectile() const
Getter for accurateProjectile.
virtual void ModelDescription(std::ostream &outFile) const
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
G4INCLXXInterfaceStore *const theInterfaceStore
G4HadFinalState theResult
G4INCL::ParticleType toINCLParticleType(G4ParticleDefinition const *const) const
Convert G4ParticleDefinition to corresponding INCL particle type.
G4ParticleDefinition * toG4ParticleDefinition(G4int A, G4int Z, G4int S, G4int PDGCode) const
Convert A, Z and S to a G4ParticleDefinition.
G4double remnant4MomentumScaling(G4double mass, G4double kineticE, G4double px, G4double py, G4double pz) const
Rescale remnant momentum if necessary.
G4INCL::ParticleSpecies toINCLParticleSpecies(G4HadProjectile const &) const
Convert G4HadProjectile to corresponding INCL particle species.
G4double toINCLKineticEnergy(G4HadProjectile const &) const
Convert G4HadProjectile to corresponding INCL particle kinetic energy.
G4INCLXXInterface(G4VPreCompoundModel *const aPreCompound=0)
G4FissionProbability * theINCLXXFissionProbability
G4DynamicParticle * toG4Particle(G4int A, G4int Z, G4int S, G4int PDGCode, G4double kinE, G4double px, G4double py, G4double pz) const
Convert an INCL particle to a G4DynamicParticle.
G4HadronicInteraction * theBackupModelNucleon
G4VLevelDensityParameter * theINCLXXLevelDensity
G4String const & GetDeExcitationModelName() const
G4HadronicInteraction * theBackupModel
G4bool AccurateProjectile(const G4HadProjectile &aTrack, const G4Nucleus &theTargetNucleus) const
G4INCLXXVInterfaceTally * theTally
G4IonTable *const theIonTable
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4VPreCompoundModel * thePreCompoundModel
G4INCL::INCL * theINCLModel
virtual void Tally(G4HadProjectile const &aTrack, G4Nucleus const &theNucleus, G4HadFinalState const &result)=0
const EventInfo & processEvent()
const G4String & GetIonName(G4int Z, G4int A, G4int lvl=0) const
Definition: G4IonTable.cc:1229
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
Definition: G4IonTable.cc:1517
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4KaonZero * KaonZero()
Definition: G4KaonZero.cc:103
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4int GetL() const
Definition: G4Nucleus.hh:108
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
G4double GetPDGCharge() const
G4int GetNumberOfLambdasInHypernucleus() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:107
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:101
static G4Triton * Triton()
Definition: G4Triton.cc:93
G4VEvaporationChannel * GetFissionChannel()
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0
G4ExcitationHandler * GetExcitationHandler() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4double energy(const ThreeVector &p, const G4double m)
float hbar_Planck
Definition: hepunit.py:263
Short_t S[maxSizeParticles]
Particle strangeness number.
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].
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.