Geant4 11.2.2
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, GSI (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 "globals.hh"
35#include <iostream>
36#include <cmath>
37
38#include "G4AblaInterface.hh"
41#include "G4ReactionProduct.hh"
42#include "G4DynamicParticle.hh"
43#include "G4IonTable.hh"
44#include "G4SystemOfUnits.hh"
48#include "G4HyperTriton.hh"
49#include "G4HyperH4.hh"
50#include "G4HyperAlpha.hh"
51#include "G4DoubleHyperH4.hh"
53#include "G4HyperHe5.hh"
54
56 G4VPreCompoundModel(ptr, "ABLAXX"),
57 ablaResult(new G4VarNtp),
58 volant(new G4Volant),
59 theABLAModel(new G4Abla(volant, 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 volant;
75 delete ablaResult;
76 delete theABLAModel;
77 delete GetExcitationHandler();
78}
79
84
86{
87 if (isInitialised) return;
88 isInitialised = true;
89 theABLAModel->initEvapora();
90 theABLAModel->SetParameters();
92}
93
94
96 G4Nucleus & theNucleus)
97{
98 // This method is adapted from G4PreCompoundModel::ApplyYourself,
99 // and it is used only by Binary Cascade (BIC) when the latter is coupled with Abla
100 // for nuclear de-excitation.
101 // This method allows BIC+ABLA to be used also for proton and neutron projectile
102 // with kinetic energies below 45 MeV, by creating a "compound" nucleus made
103 // by the system "target nucleus + projectile", before calling the DeExcite
104 // method.
105 const G4ParticleDefinition* primary = thePrimary.GetDefinition();
106 if ( primary != G4Neutron::Definition() && primary != G4Proton::Definition() ) {
108 ed << "G4AblaModel is used for ";
109 if ( primary ) ed << primary->GetParticleName();
110 G4Exception( "G4AblaInterface::ApplyYourself()", "had040", FatalException, ed, "" );
111 return nullptr;
112 }
113
114 G4int Zp = 0;
115 G4int Ap = 1;
116 if ( primary == G4Proton::Definition() ) Zp = 1;
117 G4double timePrimary = thePrimary.GetGlobalTime();
118 G4int A = theNucleus.GetA_asInt();
119 G4int Z = theNucleus.GetZ_asInt();
120 G4LorentzVector p = thePrimary.Get4Momentum();
122 p += G4LorentzVector( 0.0, 0.0, 0.0, mass );
123
124 G4Fragment anInitialState(A + Ap, Z + Zp, p);
125 anInitialState.SetNumberOfExcitedParticle(1, Zp);
126 anInitialState.SetNumberOfHoles(1, Zp);
127 anInitialState.SetCreationTime( thePrimary.GetGlobalTime() );
128 anInitialState.SetCreatorModelID( secID );
129
130 G4ReactionProductVector* deExciteResult = DeExcite( anInitialState );
131
132 applyYourselfResult.Clear();
133 applyYourselfResult.SetStatusChange( stopAndKill );
134 for ( auto const & prod : *deExciteResult ) {
135 G4DynamicParticle * aNewDP =
136 new G4DynamicParticle( prod->GetDefinition(), prod->GetTotalEnergy(), prod->GetMomentum() );
137 G4HadSecondary aNew = G4HadSecondary( aNewDP );
138 G4double time = std::max( prod->GetFormationTime(), 0.0 );
139 aNew.SetTime( timePrimary + time );
140 aNew.SetCreatorModelID( prod->GetCreatorModelID() );
141 delete prod;
142 applyYourselfResult.AddSecondary( aNew );
143 }
144 delete deExciteResult;
145 return &applyYourselfResult;
146}
147
148
150 if (!isInitialised) InitialiseModel();
151
152 volant->clear();
153 ablaResult->clear();
154
155 const G4int ARem = aFragment.GetA_asInt();
156 const G4int ZRem = aFragment.GetZ_asInt();
157 const G4int SRem = -aFragment.GetNumberOfLambdas(); // Strangeness = - (Number of lambdas)
158 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
159 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
160 const G4LorentzVector& pRem = aFragment.GetMomentum();
161 const G4double pxRem = pRem.x() / MeV;
162 const G4double pyRem = pRem.y() / MeV;
163 const G4double pzRem = pRem.z() / MeV;
164
165 ++eventNumber;
166
167 theABLAModel->DeexcitationAblaxx(ARem, ZRem, eStarRem, jRem, pxRem, pyRem,
168 pzRem, (G4int)eventNumber, SRem);
169
171
172 for(G4int j = 0; j < ablaResult->ntrack; ++j)
173 { // Copy ABLA result to the EventInfo
174 G4ReactionProduct* product =
175 toG4Particle(ablaResult->avv[j], ablaResult->zvv[j], ablaResult->svv[j],
176 ablaResult->enerj[j], ablaResult->pxlab[j],
177 ablaResult->pylab[j], ablaResult->pzlab[j]);
178 if(product)
179 {
180 product->SetCreatorModelID(secID);
181 result->push_back(product);
182 }
183 }
184 return result;
185}
186
187G4ParticleDefinition *G4AblaInterface::toG4ParticleDefinition(G4int A, G4int Z, G4int S) const {
188 if (A == 1 && Z == 1 && S == 0 ) return G4Proton::Proton();
189 else if(A == 1 && Z == 0 && S == 0 ) return G4Neutron::Neutron();
190 else if(A == 1 && Z == 0 && S == -1) return G4Lambda::Lambda();
191 else if(A == -1 && Z == 1 && S == 0 ) return G4PionPlus::PionPlus();
192 else if(A == -1 && Z == -1 && S == 0 ) return G4PionMinus::PionMinus();
193 else if(A == -1 && Z == 0 && S == 0 ) return G4PionZero::PionZero();
194 else if(A == 0 && Z == 0 && S == 0 ) return G4Gamma::Gamma();
195 else if(A == 2 && Z == 1 && S == 0 ) return G4Deuteron::Deuteron();
196 else if(A == 3 && Z == 1 && S == 0 ) return G4Triton::Triton();
197 else if(A == 3 && Z == 2 && S == 0 ) return G4He3::He3();
198 else if(A == 3 && Z == 1 && S == -1) return G4HyperTriton::Definition();
199 else if(A == 4 && Z == 2 && S == 0 ) return G4Alpha::Alpha();
200 else if(A == 4 && Z == 1 && S == -1) return G4HyperH4::Definition();
201 else if(A == 4 && Z == 2 && S == -1) return G4HyperAlpha::Definition();
202 else if(A == 4 && Z == 1 && S == -2) return G4DoubleHyperH4::Definition();
203 else if(A == 4 && Z == 0 && S == -2) return G4DoubleHyperDoubleNeutron::Definition();
204 else if(A == 5 && Z == 2 && S == -1) return G4HyperHe5::Definition();
205 else if(A > 0 && Z > 0 && A > Z )
206 { // Returns ground state ion definition.
207 auto ionfromtable = G4IonTable::GetIonTable()->GetIon(Z, A, std::abs(S), 0); // S is the number of lambdas
208 if(ionfromtable)
209 return ionfromtable;
210 else
211 {
212 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
213 << " to G4ParticleDefinition, trouble ahead" << G4endl;
214 return 0;
215 }
216 }
217 else
218 { // Error, unrecognized particle
219 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
220 << " to G4ParticleDefinition, trouble ahead" << G4endl;
221 return 0;
222 }
223}
224
225G4ReactionProduct* G4AblaInterface::toG4Particle(G4int A, G4int Z, G4int S,
226 G4double kinE, G4double px,
227 G4double py, G4double pz) const {
228 G4ParticleDefinition* def = toG4ParticleDefinition(A, Z, S);
229 if(def == 0)
230 { // Check if we have a valid particle definition
231 return 0;
232 }
233
234 const G4double energy = kinE * MeV;
235 const G4ThreeVector momentum(px, py, pz);
236 const G4ThreeVector momentumDirection = momentum.unit();
237 G4DynamicParticle p(def, momentumDirection, energy);
239 (*r) = p;
240 return r;
241}
242
243void G4AblaInterface::ModelDescription(std::ostream& outFile) const
244{
245 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
246}
247
248void G4AblaInterface::DeExciteModelDescription(std::ostream& outFile) const
249{
250 outFile
251 << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
252 << "the gamma emission and the evaporation of neutrons, light charged\n"
253 << "particles and IMFs, as well as fission where applicable. The code\n"
254 << "included in Geant4 is a C++ translation of the original Fortran\n"
255 << "code ABLA07. Although the model has been recently extended to\n"
256 << "hypernuclei by including the evaporation of lambda particles.\n"
257 << "More details about the physics are available in the Geant4\n"
258 << "Physics Reference Manual and in the reference articles.\n\n"
259 << "References:\n"
260 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of "
261 "Joint\n"
262 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
263 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. Leray, "
264 "Y. Yariv,\n"
265 << "A. Mengoni, A. Stanculescu, and G. Mank (IAEA INDC(NDS)-530, Vienna, "
266 "2008), pp. 181–221.\n\n"
267 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, "
268 "021602 (2018)\n\n";
269}
270
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
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()
void initEvapora()
Definition G4Abla.cc:2133
void SetParameters()
Definition G4Abla.cc:2320
void DeexcitationAblaxx(G4int nucleusA, G4int nucleusZ, G4double excitationEnergy, G4double angularMomentum, G4double momX, G4double momY, G4double momZ, G4int eventnumber)
Definition G4Abla.cc:96
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
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
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:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
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
void SetCreatorModelID(const G4int mod)
static G4Triton * Triton()
Definition G4Triton.cc:90
G4ExcitationHandler * GetExcitationHandler() const
void SetExcitationHandler(G4ExcitationHandler *ptr)
G4double enerj[VARNTPSIZE]
G4double pylab[VARNTPSIZE]
G4int svv[VARNTPSIZE]
G4int avv[VARNTPSIZE]
G4double pzlab[VARNTPSIZE]
G4int zvv[VARNTPSIZE]
G4double pxlab[VARNTPSIZE]
G4double energy(const ThreeVector &p, const G4double m)