Geant4-11
G4TheoFSGenerator.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//
27// G4TheoFSGenerator
28//
29// 20110307 M. Kelsey -- Add call to new theTransport->SetPrimaryProjectile()
30// to provide access to full initial state (for Bertini)
31// 20110805 M. Kelsey -- Follow change to G4V3DNucleus::GetNucleons()
32
33#include "G4DynamicParticle.hh"
34#include "G4TheoFSGenerator.hh"
36#include "G4ReactionProduct.hh"
37#include "G4IonTable.hh"
39#include "G4CRCoalescence.hh"
41
44 , theTransport(nullptr), theHighEnergyGenerator(nullptr)
45 , theQuasielastic(nullptr)
46 , theCosmicCoalescence(nullptr)
47 , theStringModelID(-1)
48{
51}
52
54{
55 delete theParticleChange;
56}
57
58void G4TheoFSGenerator::ModelDescription(std::ostream& outFile) const
59{
60 outFile << GetModelName() <<" consists of a " << theHighEnergyGenerator->GetModelName()
61 << " string model and a stage to de-excite the excited nuclear fragment.\n<p>"
62 << "The string model simulates the interaction of\n"
63 << "an incident hadron with a nucleus, forming \n"
64 << "excited strings, decays these strings into hadrons,\n"
65 << "and leaves an excited nucleus. \n"
66 << "<p>The string model:\n";
68 outFile <<"\n<p>";
70}
71
73{
74 // init particle change
77 G4double timePrimary=thePrimary.GetGlobalTime();
78
79 // Temporarily dummy treatment of heavy (charm and bottom) hadron projectiles at low energies.
80 // Cascade models are currently not applicable for heavy hadrons and string models cannot
81 // handle them properly at low energies - let's say safely below ~100 MeV.
82 // In these cases, we return as final state the initial state unchanged.
83 // For most applications, this is a safe simplification, giving that the nearly all
84 // slowly moving charm and bottom hadrons decay before any hadronic interaction can occur.
85 // Note that we prefer not to use G4HadronicParameters::GetMinEnergyTransitionFTF_Cascade()
86 // (typicall ~3 GeV) because FTFP works reasonably well below such a value.
87 const G4double energyThresholdForCharmAndBottomHadrons = 100.0*CLHEP::MeV;
88 if ( thePrimary.GetKineticEnergy() < energyThresholdForCharmAndBottomHadrons &&
89 ( thePrimary.GetDefinition()->GetQuarkContent( 4 ) != 0 || // Has charm constituent quark
90 thePrimary.GetDefinition()->GetAntiQuarkContent( 4 ) != 0 || // Has anti-charm constituent anti-quark
91 thePrimary.GetDefinition()->GetQuarkContent( 5 ) != 0 || // Has bottom constituent quark
92 thePrimary.GetDefinition()->GetAntiQuarkContent( 5 ) != 0 ) ) { // Has anti-bottom constituent anti-quark
96 return theParticleChange;
97 }
98
99 // check if models have been registered, and use default, in case this is not true @@
100
101 const G4DynamicParticle aPart(thePrimary.GetDefinition(),thePrimary.Get4Momentum().vect());
102
103 if ( theQuasielastic )
104 {
105 if ( theQuasielastic->GetFraction(theNucleus, aPart) > G4UniformRand() )
106 {
107 //G4cout<<"___G4TheoFSGenerator: before Scatter (1) QE=" << theQuasielastic<<G4endl;
108 G4KineticTrackVector *result= theQuasielastic->Scatter(theNucleus, aPart);
109 //G4cout << "^^G4TheoFSGenerator: after Scatter (1) " << G4endl;
110 if (result)
111 {
112 for(auto & ptr : *result)
113 {
114 G4DynamicParticle * aNew =
115 new G4DynamicParticle(ptr->GetDefinition(),
116 ptr->Get4Momentum().e(),
117 ptr->Get4Momentum().vect());
118 theParticleChange->AddSecondary(aNew, ptr->GetCreatorModelID());
119 delete ptr;
120 }
121 delete result;
122 }
123 else
124 {
128 }
129 return theParticleChange;
130 }
131 }
132
133 // get result from high energy model
134 G4KineticTrackVector * theInitialResult =
135 theHighEnergyGenerator->Scatter(theNucleus, aPart);
136
137 // Assign the creator model ID
138 for ( auto & ptr : *theInitialResult ) {
139 ptr->SetCreatorModelID( theStringModelID );
140 }
141
142 //#define DEBUG_initial_result
143 #ifdef DEBUG_initial_result
144 G4double E_out(0);
146 for(auto & ptr : *theInitialResult)
147 {
148 //G4cout << "TheoFS secondary, mom " << ptr->GetDefinition()->GetParticleName()
149 // << " " << ptr->Get4Momentum() << G4endl;
150 E_out += ptr->Get4Momentum().e();
151 }
152 G4double init_mass= ionTable->GetIonMass(theNucleus.GetZ_asInt(),theNucleus.GetA_asInt());
153 G4double init_E=aPart.Get4Momentum().e();
154 // residual nucleus
155
156 const std::vector<G4Nucleon> & thy = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
157
158 G4int resZ(0),resA(0);
159 G4double delta_m(0);
160 for(auto & nuc : thy)
161 {
162 if(nuc.AreYouHit()) {
163 ++resA;
164 if ( nuc.GetDefinition() == G4Proton::Proton() ) {
165 ++resZ;
166 delta_m += CLHEP::proton_mass_c2;
167 } else {
168 delta_m += CLHEP::neutron_mass_c2;
169 }
170 }
171 }
172
173 G4double final_mass(0);
174 if ( theNucleus.GetA_asInt() ) {
175 final_mass=ionTable->GetIonMass(theNucleus.GetZ_asInt()-resZ,theNucleus.GetA_asInt()- resA);
176 }
177 G4double E_excit=init_mass + init_E - final_mass - E_out;
178 G4cout << " Corrected delta mass " << init_mass - final_mass - delta_m << G4endl;
179 G4cout << "initial E, mass = " << init_E << ", " << init_mass << G4endl;
180 G4cout << " final E, mass = " << E_out <<", " << final_mass << " excitation_E " << E_excit << G4endl;
181 #endif
182
183 G4ReactionProductVector * theTransportResult = nullptr;
184
185 // Uzhi Nov. 2012 for nucleus-nucleus collision
186 G4V3DNucleus* theProjectileNucleus = theHighEnergyGenerator->GetProjectileNucleus();
187 if(theProjectileNucleus == nullptr)
188 {
189 G4int hitCount = 0;
190 const std::vector<G4Nucleon>& they = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
191 for(auto & nuc : they)
192 {
193 if(nuc.AreYouHit()) ++hitCount;
194 }
196 {
197 theTransport->SetPrimaryProjectile(thePrimary); // For Bertini Cascade
198 theTransportResult =
199 theTransport->Propagate(theInitialResult,
201 if ( !theTransportResult ) {
202 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
203 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
204 }
205 }
206 else
207 {
208 theTransportResult = theDecay.Propagate(theInitialResult,
210 if ( theTransportResult == nullptr ) {
211 G4cout << "G4TheoFSGenerator: null ptr from decay propagate " << G4endl;
212 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from decay propagate");
213 }
214 }
215 }
216 else
217 {
219 theTransportResult = theTransport->PropagateNuclNucl(theInitialResult,
221 theProjectileNucleus);
222 if ( !theTransportResult ) {
223 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
224 throw G4HadronicException(__FILE__, __LINE__, "Null ptr from transport propagate");
225 }
226 }
227
228 // If enabled, apply the Cosmic Rays (CR) coalescence to the list of secondaries produced so far.
229 // This algorithm can form deuterons and antideuterons by coalescence of, respectively,
230 // proton-neutron and antiproton-antineutron pairs close in momentum space.
231 // This can be useful in particular for Cosmic Ray applications.
232 if ( G4HadronicParameters::Instance()->EnableCRCoalescence() ) {
233 if(nullptr == theCosmicCoalescence) {
236 if(nullptr == theCosmicCoalescence) {
238 }
239 }
241 theCosmicCoalescence->GenerateDeuterons( theTransportResult );
242 }
243
244 // Fill particle change
245 for(auto & ptr : *theTransportResult)
246 {
247 G4DynamicParticle * aNewDP =
248 new G4DynamicParticle(ptr->GetDefinition(),
249 ptr->GetTotalEnergy(),
250 ptr->GetMomentum());
251 G4HadSecondary aNew = G4HadSecondary(aNewDP);
252 G4double time = std::max(ptr->GetFormationTime(), 0.0);
253 aNew.SetTime(timePrimary + time);
254 aNew.SetCreatorModelID(ptr->GetCreatorModelID());
256 delete ptr;
257 }
258
259 // some garbage collection
260 delete theTransportResult;
261
262 // Done
263 return theParticleChange;
264}
265
266std::pair<G4double, G4double> G4TheoFSGenerator::GetEnergyMomentumCheckLevels() const
267{
270 } else {
271 return std::pair<G4double, G4double>(DBL_MAX, DBL_MAX);
272 }
273}
@ isAlive
@ stopAndKill
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
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
void GenerateDeuterons(G4ReactionProductVector *result)
void SetP0Coalescence(const G4HadProjectile &thePrimary, G4String)
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *)
G4LorentzVector Get4Momentum() const
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
const G4String & GetModelName() const
static G4HadronicParameters * Instance()
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
Definition: G4IonTable.cc:1517
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4int GetQuarkContent(G4int flavor) const
G4int GetAntiQuarkContent(G4int flavor) const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4KineticTrackVector * Scatter(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
G4double GetFraction(G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
G4CRCoalescence * theCosmicCoalescence
G4VHighEnergyGenerator * theHighEnergyGenerator
G4TheoFSGenerator(const G4String &name="TheoFSGenerator")
G4QuasiElasticChannel * theQuasielastic
G4VIntraNuclearTransportModel * theTransport
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const override
G4HadFinalState * theParticleChange
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override
void ModelDescription(std::ostream &outFile) const override
~G4TheoFSGenerator() override
G4DecayStrongResonances theDecay
virtual const std::vector< G4Nucleon > & GetNucleons()=0
virtual G4int GetMassNumber()=0
void ModelDescription(std::ostream &) const override
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
virtual G4ReactionProductVector * PropagateNuclNucl(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
void SetPrimaryProjectile(const G4HadProjectile &aPrimary)
virtual void PropagateModelDescription(std::ostream &outFile) const
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0
static constexpr double proton_mass_c2
static constexpr double neutron_mass_c2
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
#define DBL_MAX
Definition: templates.hh:62