Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 delete volant;
74 delete ablaResult;
75 delete theABLAModel;
76 delete GetExcitationHandler();
77}
78
80{
82}
83
85{
86 if (isInitialised) return;
87 isInitialised = true;
88 theABLAModel->initEvapora();
89 theABLAModel->SetParameters();
91}
92
94 if (!isInitialised) InitialiseModel();
95
96 volant->clear();
97 ablaResult->clear();
98
99 const G4int ARem = aFragment.GetA_asInt();
100 const G4int ZRem = aFragment.GetZ_asInt();
101 const G4int SRem = -aFragment.GetNumberOfLambdas(); // Strangeness = - (Number of lambdas)
102 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
103 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
104 const G4LorentzVector& pRem = aFragment.GetMomentum();
105 const G4double pxRem = pRem.x() / MeV;
106 const G4double pyRem = pRem.y() / MeV;
107 const G4double pzRem = pRem.z() / MeV;
108
109 ++eventNumber;
110
111 theABLAModel->DeexcitationAblaxx(ARem, ZRem, eStarRem, jRem, pxRem, pyRem,
112 pzRem, (G4int)eventNumber, SRem);
113
115
116 for(G4int j = 0; j < ablaResult->ntrack; ++j)
117 { // Copy ABLA result to the EventInfo
118 G4ReactionProduct* product =
119 toG4Particle(ablaResult->avv[j], ablaResult->zvv[j], ablaResult->svv[j],
120 ablaResult->enerj[j], ablaResult->pxlab[j],
121 ablaResult->pylab[j], ablaResult->pzlab[j]);
122 if(product)
123 {
124 product->SetCreatorModelID(secID);
125 result->push_back(product);
126 }
127 }
128 return result;
129}
130
131G4ParticleDefinition *G4AblaInterface::toG4ParticleDefinition(G4int A, G4int Z, G4int S) const {
132 if (A == 1 && Z == 1 && S == 0 ) return G4Proton::Proton();
133 else if(A == 1 && Z == 0 && S == 0 ) return G4Neutron::Neutron();
134 else if(A == 1 && Z == 0 && S == -1) return G4Lambda::Lambda();
135 else if(A == -1 && Z == 1 && S == 0 ) return G4PionPlus::PionPlus();
136 else if(A == -1 && Z == -1 && S == 0 ) return G4PionMinus::PionMinus();
137 else if(A == -1 && Z == 0 && S == 0 ) return G4PionZero::PionZero();
138 else if(A == 0 && Z == 0 && S == 0 ) return G4Gamma::Gamma();
139 else if(A == 2 && Z == 1 && S == 0 ) return G4Deuteron::Deuteron();
140 else if(A == 3 && Z == 1 && S == 0 ) return G4Triton::Triton();
141 else if(A == 3 && Z == 2 && S == 0 ) return G4He3::He3();
142 else if(A == 3 && Z == 1 && S == -1) return G4HyperTriton::Definition();
143 else if(A == 4 && Z == 2 && S == 0 ) return G4Alpha::Alpha();
144 else if(A == 4 && Z == 1 && S == -1) return G4HyperH4::Definition();
145 else if(A == 4 && Z == 2 && S == -1) return G4HyperAlpha::Definition();
146 else if(A == 4 && Z == 1 && S == -2) return G4DoubleHyperH4::Definition();
147 else if(A == 4 && Z == 0 && S == -2) return G4DoubleHyperDoubleNeutron::Definition();
148 else if(A == 5 && Z == 2 && S == -1) return G4HyperHe5::Definition();
149 else if(A > 0 && Z > 0 && A > Z )
150 { // Returns ground state ion definition.
151 auto ionfromtable = G4IonTable::GetIonTable()->GetIon(Z, A, std::abs(S), 0); // S is the number of lambdas
152 if(ionfromtable)
153 return ionfromtable;
154 else
155 {
156 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
157 << " to G4ParticleDefinition, trouble ahead" << G4endl;
158 return 0;
159 }
160 }
161 else
162 { // Error, unrecognized particle
163 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S
164 << " to G4ParticleDefinition, trouble ahead" << G4endl;
165 return 0;
166 }
167}
168
169G4ReactionProduct* G4AblaInterface::toG4Particle(G4int A, G4int Z, G4int S,
170 G4double kinE, G4double px,
171 G4double py, G4double pz) const {
172 G4ParticleDefinition* def = toG4ParticleDefinition(A, Z, S);
173 if(def == 0)
174 { // Check if we have a valid particle definition
175 return 0;
176 }
177
178 const G4double energy = kinE * MeV;
179 const G4ThreeVector momentum(px, py, pz);
180 const G4ThreeVector momentumDirection = momentum.unit();
181 G4DynamicParticle p(def, momentumDirection, energy);
183 (*r) = p;
184 return r;
185}
186
187void G4AblaInterface::ModelDescription(std::ostream& outFile) const
188{
189 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
190}
191
192void G4AblaInterface::DeExciteModelDescription(std::ostream& outFile) const
193{
194 outFile
195 << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
196 << "the gamma emission and the evaporation of neutrons, light charged\n"
197 << "particles and IMFs, as well as fission where applicable. The code\n"
198 << "included in Geant4 is a C++ translation of the original Fortran\n"
199 << "code ABLA07. Although the model has been recently extended to\n"
200 << "hypernuclei by including the evaporation of lambda particles.\n"
201 << "More details about the physics are available in the Geant4\n"
202 << "Physics Reference Manual and in the reference articles.\n\n"
203 << "References:\n"
204 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of "
205 "Joint\n"
206 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
207 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. Leray, "
208 "Y. Yariv,\n"
209 << "A. Mengoni, A. Stanculescu, and G. Mank (IAEA INDC(NDS)-530, Vienna, "
210 "2008), pp. 181–221.\n\n"
211 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, "
212 "021602 (2018)\n\n";
213}
214
G4double S(G4double temp)
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
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 void InitialiseModel() final
virtual void ModelDescription(std::ostream &outFile) const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
G4AblaInterface(G4ExcitationHandler *ptr=nullptr)
virtual ~G4AblaInterface()
Definition: G4Abla.hh:54
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:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4DoubleHyperDoubleNeutron * Definition()
static G4DoubleHyperH4 * Definition()
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:307
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:330
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
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 G4HyperTriton * Definition()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
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 * Proton()
Definition: G4Proton.cc:92
void SetCreatorModelID(const G4int mod)
static G4Triton * Triton()
Definition: G4Triton.cc:93
G4ExcitationHandler * GetExcitationHandler() const
void SetExcitationHandler(G4ExcitationHandler *ptr)
void clear()
G4double enerj[VARNTPSIZE]
G4double pylab[VARNTPSIZE]
G4int svv[VARNTPSIZE]
G4int avv[VARNTPSIZE]
G4double pzlab[VARNTPSIZE]
G4int zvv[VARNTPSIZE]
G4double pxlab[VARNTPSIZE]
void clear()
G4double energy(const ThreeVector &p, const G4double m)