Geant4 10.7.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, 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#define ABLAXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
38#ifdef ABLAXX_IN_GEANT4_MODE
39
40#include "G4AblaInterface.hh"
43#include "G4ReactionProduct.hh"
44#include "G4DynamicParticle.hh"
45#include "G4IonTable.hh"
46#include "G4SystemOfUnits.hh"
48#include <iostream>
49#include <cmath>
50
52 G4VPreCompoundModel(NULL, "ABLA"),
53 ablaResult(new G4VarNtp),
54 volant(new G4Volant),
55 theABLAModel(new G4Abla(volant, ablaResult)),
56 eventNumber(0)
57{
58 theABLAModel->initEvapora();
59 theABLAModel->SetParameters();
60}
61
63 delete volant;
64 delete ablaResult;
65 delete theABLAModel;
66}
67
69 volant->clear();
70 ablaResult->clear();
71
72 const G4int ARem = aFragment.GetA_asInt();
73 const G4int ZRem = aFragment.GetZ_asInt();
74 const G4double eStarRem = aFragment.GetExcitationEnergy() / MeV;
75 const G4double jRem = aFragment.GetAngularMomentum().mag() / hbar_Planck;
76 const G4LorentzVector &pRem = aFragment.GetMomentum();
77 const G4double pxRem = pRem.x() / MeV;
78 const G4double pyRem = pRem.y() / MeV;
79 const G4double pzRem = pRem.z() / MeV;
80
81 eventNumber++;
82
83 theABLAModel->DeexcitationAblaxx(ARem, ZRem,
84 eStarRem,
85 jRem,
86 pxRem,
87 pyRem,
88 pzRem,
89 eventNumber);
90
92
93 for(int j = 0; j < ablaResult->ntrack; ++j) { // Copy ABLA result to the EventInfo
94
95 G4ReactionProduct *product = toG4Particle(ablaResult->avv[j],
96 ablaResult->zvv[j],
97 ablaResult->svv[j],
98 ablaResult->enerj[j],
99 ablaResult->pxlab[j],
100 ablaResult->pylab[j],
101 ablaResult->pzlab[j]);
102 if(product)
103 result->push_back(product);
104 }
105 return result;
106}
107
108G4ParticleDefinition *G4AblaInterface::toG4ParticleDefinition(G4int A, G4int Z, G4int S) const {
109 if (A == 1 && Z == 1 && S == 0) return G4Proton::Proton();
110 else if(A == 1 && Z == 0 && S == 0) return G4Neutron::Neutron();
111 else if(A == 1 && Z == 0 && S == -1) return G4Lambda::Lambda();
112 else if(A == -1 && Z == 1 && S == 0) return G4PionPlus::PionPlus();
113 else if(A == -1 && Z == -1 && S == 0) return G4PionMinus::PionMinus();
114 else if(A == -1 && Z == 0 && S == 0) return G4PionZero::PionZero();
115 else if(A == 0 && Z == 0 && S == 0) return G4Gamma::Gamma();
116 else if(A == 2 && Z == 1 && S == 0) return G4Deuteron::Deuteron();
117 else if(A == 3 && Z == 1 && S == 0) return G4Triton::Triton();
118 else if(A == 3 && Z == 2 && S == 0) return G4He3::He3();
119 else if(A == 4 && Z == 2 && S == 0) return G4Alpha::Alpha();
120 else if(A > 0 && Z > 0 && A > Z) { // Returns ground state ion definition.
121 return G4IonTable::GetIonTable()->GetIon(Z, A, std::abs(S));//S is the number of lambdas
122 } else { // Error, unrecognized particle
123 G4cout << "Can't convert particle with A=" << A << ", Z=" << Z << ", S=" << S << " to G4ParticleDefinition, trouble ahead" << G4endl;
124 return 0;
125 }
126}
127
128G4ReactionProduct *G4AblaInterface::toG4Particle(G4int A, G4int Z, G4int S,
129 G4double kinE,
130 G4double px,
131 G4double py, G4double pz) const {
132 G4ParticleDefinition *def = toG4ParticleDefinition(A, Z, S);
133 if(def == 0) { // Check if we have a valid particle definition
134 return 0;
135 }
136
137 const G4double energy = kinE * MeV;
138 const G4ThreeVector momentum(px, py, pz);
139 const G4ThreeVector momentumDirection = momentum.unit();
140 G4DynamicParticle p(def, momentumDirection, energy);
142 (*r) = p;
143 return r;
144}
145
146void G4AblaInterface::ModelDescription(std::ostream& outFile) const {
147 outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
148}
149
150void G4AblaInterface::DeExciteModelDescription(std::ostream& outFile) const {
151 outFile
152 << "ABLA++ is a statistical model for nuclear de-excitation. It simulates\n"
153 << "the gamma emission and the evaporation of neutrons, light charged particles\n"
154 << "and IMFs, as well as fission where applicable. The code included in Geant4\n"
155 << "is a C++ translation of the original Fortran code ABLA07. Although the model\n"
156 << "has been recently extended to hypernuclei by including the evaporation of lambda\n"
157 << "particles. More details about the physics are available in the\n"
158 << "Geant4 Physics Reference Manual and in the reference articles.\n\n"
159 << "References:\n"
160 << "(1) A. Kelic, M. V. Ricciardi, and K. H. Schmidt, in Proceedings of Joint\n"
161 << "ICTP-IAEA Advanced Workshop on Model Codes for Spallation Reactions,\n"
162 << "ICTP Trieste, Italy, 4–8 February 2008, edited by D. Filges, S. Leray, Y. Yariv,\n"
163 << "A. Mengoni, A. Stanculescu, and G. Mank (IAEA INDC(NDS)-530, Vienna, 2008), pp. 181–221.\n\n"
164 << "(2) J.L. Rodriguez-Sanchez, J.-C. David et al., Phys. Rev. C 98, 021602 (2018)\n\n";
165}
166
167#endif // ABLAXX_IN_GEANT4_MODE
double S(double temp)
double A(double temperature)
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
double mag() const
virtual void DeExciteModelDescription(std::ostream &outFile) const
virtual void ModelDescription(std::ostream &outFile) const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
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
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:275
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
G4int GetZ_asInt() const
Definition: G4Fragment.hh:263
G4ThreeVector GetAngularMomentum() const
Definition: G4Fragment.cc:254
G4int GetA_asInt() const
Definition: G4Fragment.hh:258
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4He3 * He3()
Definition: G4He3.cc:93
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 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
static G4Triton * Triton()
Definition: G4Triton.cc:94
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)