ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4AblaInterface.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4AblaInterface.cc
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"
41 #include "G4ParticleDefinition.hh"
43 #include "G4ReactionProduct.hh"
44 #include "G4DynamicParticle.hh"
45 #include "G4IonTable.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4PhysicalConstants.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 {
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 
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 
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 
129  G4double kinE,
130  G4double px,
131  G4double py, G4double pz) const {
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 
146 void G4AblaInterface::ModelDescription(std::ostream& outFile) const {
147  outFile << "ABLA++ does not provide an implementation of the ApplyYourself method!\n\n";
148 }
149 
150 void 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