Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AblaInterface.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// ABLAXX statistical de-excitation model
27// Jose Luis Rodriguez, UDC (translation from ABLA07 and contact person)
28// Pekka Kaitaniemi, HIP (initial translation of ablav3p)
29// Aleksandra Kelic, GSI (ABLA07 code)
30// Davide Mancusi, CEA (contact person INCL)
31// Aatos Heikkinen, HIP (project coordination)
32//
33
34#include "G4AblaInterface.hh"
35
37#include "G4DoubleHyperH4.hh"
38#include "G4DynamicParticle.hh"
40#include "G4HyperAlpha.hh"
41#include "G4HyperH4.hh"
42#include "G4HyperHe5.hh"
43#include "G4HyperTriton.hh"
44#include "G4IonTable.hh"
48#include "G4ReactionProduct.hh"
50#include "G4SystemOfUnits.hh"
51#include "globals.hh"
52
53#include <cmath>
54#include <iostream>
55
57 : G4VPreCompoundModel(ptr, "ABLAXX"),
58 ablaResult(new G4VarNtp),
59 theABLAModel(new G4Abla(ablaResult)),
60 eventNumber(0),
61 secID(-1),
62 isInitialised(false)
63{
65 // G4cout << "### NEW PrecompoundModel " << this << G4endl;
68 G4cout << G4endl << "G4AblaInterface::InitialiseModel() was right." << G4endl;
69}
70
72{
73 applyYourselfResult.Clear();
74 delete ablaResult;
75 delete theABLAModel;
76 delete GetExcitationHandler();
77}
78
83
85{
86 if (isInitialised) return;
87 isInitialised = true;
88 theABLAModel->initEvapora();
89 theABLAModel->SetParameters();
91}
92
94 G4Nucleus& theNucleus)
95{
96 // This method is adapted from G4PreCompoundModel::ApplyYourself,
97 // and it is used only by Binary Cascade (BIC) when the latter is coupled with
98 // Abla for nuclear de-excitation. This method allows BIC+ABLA to be used also
99 // for proton and neutron projectile with kinetic energies below 45 MeV, by
100 // creating a "compound" nucleus made by the system "target nucleus +
101 // projectile", before calling the DeExcite method.
102 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
103 if (primary != G4Neutron::Definition() && primary != G4Proton::Definition()) {
105 ed << "G4AblaModel is used for ";
106 if (primary) ed << primary->GetParticleName();
107 G4Exception("G4AblaInterface::ApplyYourself()", "had040", FatalException, ed, "");
108 return nullptr;
109 }
110
111 G4int Zp = 0;
112 G4int Ap = 1;
113 if (primary == G4Proton::Definition()) Zp = 1;
114 G4double timePrimary = thePrimary.GetGlobalTime();
115 G4int A = theNucleus.GetA_asInt();
116 G4int Z = theNucleus.GetZ_asInt();
117 G4LorentzVector p = thePrimary.Get4Momentum();
119 p += G4LorentzVector(0.0, 0.0, 0.0, mass);
120
121 G4Fragment anInitialState(A + Ap, Z + Zp, p);
122 anInitialState.SetNumberOfExcitedParticle(1, Zp);
123 anInitialState.SetNumberOfHoles(1, Zp);
124 anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
125 anInitialState.SetCreatorModelID(secID);
126
127 G4ReactionProductVector* deExciteResult = DeExcite(anInitialState);
128
129 applyYourselfResult.Clear();
130 applyYourselfResult.SetStatusChange(stopAndKill);
131 for (auto const& prod : *deExciteResult) {
132 G4DynamicParticle* aNewDP =
133 new G4DynamicParticle(prod->GetDefinition(), prod->GetTotalEnergy(), prod->GetMomentum());
134 G4HadSecondary aNew = G4HadSecondary(aNewDP);
135 G4double time = std::max(prod->GetFormationTime(), 0.0);
136 aNew.SetTime(timePrimary + time);
137 aNew.SetCreatorModelID(prod->GetCreatorModelID());
138 delete prod;
139 applyYourselfResult.AddSecondary(aNew);
140 }
141 delete deExciteResult;
142 return &applyYourselfResult;
143}
144
146{
147 if (!isInitialised) InitialiseModel();
148
149 ablaResult->clear();
150
151 const G4int ARem = aFragment.GetA_asInt();
152 const G4int ZRem = aFragment.GetZ_asInt();
153 const G4int SRem = -aFragment.GetNumberOfLambdas(); // Strangeness = - (Number of lambdas)
154 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
155 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
156 const G4LorentzVector& pRem = aFragment.GetMomentum();
157 const G4double pxRem = pRem.x() / MeV;
158 const G4double pyRem = pRem.y() / MeV;
159 const G4double pzRem = pRem.z() / MeV;
160
161 ++eventNumber;
162
163 theABLAModel->DeexcitationAblaxx(ARem, ZRem, eStarRem, jRem, pxRem, pyRem, pzRem,
164 (G4int)eventNumber, SRem);
165
167
168 for (G4int j = 0; j < ablaResult->ntrack; ++j) { // Copy ABLA result to the EventInfo
169 G4ReactionProduct* product =
170 toG4Particle(ablaResult->avv[j], ablaResult->zvv[j], ablaResult->svv[j], ablaResult->enerj[j],
171 ablaResult->pxlab[j], ablaResult->pylab[j], ablaResult->pzlab[j]);
172 if (product) {
173 product->SetCreatorModelID(secID);
174 result->push_back(product);
175 }
176 }
177 return result;
178}
179
180G4ParticleDefinition* G4AblaInterface::toG4ParticleDefinition(G4int A, G4int Z, G4int S) const
181{
182 if (A == 1 && Z == 1 && S == 0)
183 return G4Proton::Proton();
184 else if (A == 1 && Z == 0 && S == 0)
185 return G4Neutron::Neutron();
186 else if (A == 1 && Z == 0 && S == -1)
187 return G4Lambda::Lambda();
188 else if (A == -1 && Z == 1 && S == 0)
189 return G4PionPlus::PionPlus();
190 else if (A == -1 && Z == -1 && S == 0)
191 return G4PionMinus::PionMinus();
192 else if (A == -1 && Z == 0 && S == 0)
193 return G4PionZero::PionZero();
194 else if (A == 0 && Z == 0 && S == 0)
195 return G4Gamma::Gamma();
196 else if (A == 2 && Z == 1 && S == 0)
197 return G4Deuteron::Deuteron();
198 else if (A == 3 && Z == 1 && S == 0)
199 return G4Triton::Triton();
200 else if (A == 3 && Z == 2 && S == 0)
201 return G4He3::He3();
202 else if (A == 3 && Z == 1 && S == -1)
204 else if (A == 4 && Z == 2 && S == 0)
205 return G4Alpha::Alpha();
206 else if (A == 4 && Z == 1 && S == -1)
207 return G4HyperH4::Definition();
208 else if (A == 4 && Z == 2 && S == -1)
210 else if (A == 4 && Z == 1 && S == -2)
212 else if (A == 4 && Z == 0 && S == -2)
214 else if (A == 5 && Z == 2 && S == -1)
215 return G4HyperHe5::Definition();
216 else if (A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition.
217 auto ionfromtable =
218 G4IonTable::GetIonTable()->GetIon(Z, A, std::abs(S), 0); // S is the number of lambdas
219 if (ionfromtable)
220 return ionfromtable;
221 else {
222 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
223 << " to G4ParticleDefinition, trouble ahead" << G4endl;
224 return 0;
225 }
226 }
227 else { // Error, unrecognized particle
228 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
229 << " to G4ParticleDefinition, trouble ahead" << G4endl;
230 return 0;
231 }
232}
233
234G4ReactionProduct* G4AblaInterface::toG4Particle(G4int A, G4int Z, G4int S, G4double kinE,
235 G4double px, G4double py, G4double pz) const
236{
237 G4ParticleDefinition* def = toG4ParticleDefinition(A, Z, S);
238 if (def == 0) { // Check if we have a valid particle definition
239 return 0;
240 }
241
242 const G4double energy = kinE * MeV;
243 const G4ThreeVector momentum(px, py, pz);
244 const G4ThreeVector momentumDirection = momentum.unit();
245 G4DynamicParticle p(def, momentumDirection, energy);
246 G4ReactionProduct* r = new G4ReactionProduct(def);
247 (*r) = p;
248 return r;
249}
250
251void G4AblaInterface::ModelDescription(std::ostream& outFile) const
252{
253 outFile << "ABLA++ does not provide an implementation of the ApplyYourself "
254 "method!\n\n";
255}
256
257void G4AblaInterface::DeExciteModelDescription(std::ostream& outFile) const
258{
259 outFile << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
260 << "the gamma emission and the evaporation of neutrons, light charged\n"
261 << "particles and IMFs, as well as fission where applicable. The code\n"
262 << "included in Geant4 is a C++ translation of the original Fortran\n"
263 << "code ABLA07. Although the model has been recently extended to\n"
264 << "hypernuclei by including the evaporation of lambda particles.\n"
265 << "More details about the physics are available in the Geant4\n"
266 << "Physics Reference Manual and in the reference articles.\n\n"
267 << "References:\n"
268 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of Joint\n"
269 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
270 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. "
271 "Leray, Y. Yariv, A. Mengoni, A. Stanculescu, and G. Mank (IAEA "
272 "INDC(NDS)-530, Vienna, 2008), pp. 181–221.\n\n"
273 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, 021602R (2018)\n"
274 << "(3) J.L. Rodriguez-Sanchez et al., Phys. Rev. C 105, 014623 (2022)\n"
275 << "(4) J.L. Rodriguez-Sanchez et al., Phys. Rev. Lett. 130, 132501 (2023)\n\n";
276}
G4double S(G4double temp)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
double mag() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &) final
virtual void DeExciteModelDescription(std::ostream &outFile) const
virtual G4HadFinalState * ApplyYourself(G4HadProjectile const &, G4Nucleus &) final
virtual void InitialiseModel() final
virtual void ModelDescription(std::ostream &outFile) const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
G4AblaInterface(G4ExcitationHandler *ptr=nullptr)
virtual ~G4AblaInterface()
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4DoubleHyperDoubleNeutron * Definition()
static G4DoubleHyperH4 * Definition()
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
G4int GetZ_asInt() const
void SetCreationTime(G4double time)
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
G4int GetNumberOfLambdas() const
G4ThreeVector GetAngularMomentum() const
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
G4int GetA_asInt() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
const G4String & GetModelName() const
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 G4HyperTriton * Definition()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
static G4Lambda * Lambda()
Definition G4Lambda.cc:105
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4Neutron * Definition()
Definition G4Neutron.cc:50
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
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4PionZero * PionZero()
static G4Proton * Definition()
Definition G4Proton.cc:45
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4Triton * Triton()
Definition G4Triton.cc:90
G4ExcitationHandler * GetExcitationHandler() const
G4VPreCompoundModel(G4ExcitationHandler *ptr=nullptr, const G4String &modelName="PrecompoundModel")
void SetExcitationHandler(G4ExcitationHandler *ptr)
G4double energy(const ThreeVector &p, const G4double m)