Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
42
44 : G4HadronicInteraction( name )
45 , theTransport( nullptr ), theHighEnergyGenerator( nullptr )
46 , theQuasielastic( nullptr ), theCosmicCoalescence( nullptr )
47 , theStringModelID( -1 ), theFTFQuasiElasticID( -1 ), theFTFDiffractiveID( -1 )
48{
49 theParticleChange = new G4HadFinalState;
50 theStringModelID = G4PhysicsModelCatalog::GetModelID( "model_" + name );
51 theFTFQuasiElasticID = G4PhysicsModelCatalog::GetModelID( "model_QuasiElastic_FTF" );
52 theFTFDiffractiveID = G4PhysicsModelCatalog::GetModelID( "model_Diffractive_FTF" );
53}
54
55
57{
58 delete theParticleChange;
59}
60
61
62void G4TheoFSGenerator::ModelDescription( std::ostream& outFile ) const
63{
64 outFile << GetModelName() <<" consists of a " << theHighEnergyGenerator->GetModelName()
65 << " string model and a stage to de-excite the excited nuclear fragment.\n<p>"
66 << "The string model simulates the interaction of\n"
67 << "an incident hadron with a nucleus, forming \n"
68 << "excited strings, decays these strings into hadrons,\n"
69 << "and leaves an excited nucleus. \n"
70 << "<p>The string model:\n";
71 theHighEnergyGenerator->ModelDescription( outFile );
72 outFile << "\n<p>";
73 theTransport->PropagateModelDescription( outFile );
74}
75
76
78{
79 // init particle change
80 theParticleChange->Clear();
81 theParticleChange->SetStatusChange( stopAndKill );
82 G4double timePrimary = thePrimary.GetGlobalTime();
83
84 // Temporarily dummy treatment of heavy (charm and bottom) hadron projectiles at low energies.
85 // Cascade models are currently not applicable for heavy hadrons and string models cannot
86 // handle them properly at low energies - let's say safely below ~100 MeV.
87 // In these cases, we return as final state the initial state unchanged.
88 // For most applications, this is a safe simplification, giving that the nearly all
89 // slowly moving charm and bottom hadrons decay before any hadronic interaction can occur.
90 // Note that we prefer not to use G4HadronicParameters::GetMinEnergyTransitionFTF_Cascade()
91 // (typically ~3 GeV) because FTFP works reasonably well below such a value.
92 const G4double energyThresholdForCharmAndBottomHadrons = 100.0*CLHEP::MeV;
93 if ( thePrimary.GetKineticEnergy() < energyThresholdForCharmAndBottomHadrons &&
94 ( thePrimary.GetDefinition()->GetQuarkContent( 4 ) != 0 || // Has charm constituent quark
95 thePrimary.GetDefinition()->GetAntiQuarkContent( 4 ) != 0 || // Has anti-charm constituent anti-quark
96 thePrimary.GetDefinition()->GetQuarkContent( 5 ) != 0 || // Has bottom constituent quark
97 thePrimary.GetDefinition()->GetAntiQuarkContent( 5 ) != 0 ) ) { // Has anti-bottom constituent anti-quark
98 theParticleChange->SetStatusChange( isAlive );
99 theParticleChange->SetEnergyChange( thePrimary.GetKineticEnergy() );
100 theParticleChange->SetMomentumChange( thePrimary.Get4Momentum().vect().unit() );
101 return theParticleChange;
102 }
103
104 // Temporarily dummy treatment of light hypernuclei projectiles at low energies.
105 // Bertini and Binary intra-nuclear cascade models are currently not applicable
106 // for hypernuclei and string models cannot handle them properly at low energies -
107 // let's say safely below ~100 MeV.
108 // In these cases, we return as final state the initial state unchanged.
109 // For most applications, this is a safe simplification, giving that the nearly all
110 // slowly moving hypernuclei decay before any hadronic interaction can occur.
111 // Note that we prefer not to use G4HadronicParameters::GetMinEnergyTransitionFTF_Cascade()
112 // (typically ~3 GeV) because FTFP works reasonably well below such a value.
113 // Note that for anti-hypernuclei, FTF model works fine down to zero kinetic energy,
114 // so there is no need of any extra, dummy treatment.
115 const G4double energyThresholdForHyperNuclei = 100.0*CLHEP::MeV;
116 if ( thePrimary.GetKineticEnergy() < energyThresholdForHyperNuclei &&
117 thePrimary.GetDefinition()->IsHypernucleus() ) {
118 theParticleChange->SetStatusChange( isAlive );
119 theParticleChange->SetEnergyChange( thePrimary.GetKineticEnergy() );
120 theParticleChange->SetMomentumChange( thePrimary.Get4Momentum().vect().unit() );
121 return theParticleChange;
122 }
123
124 // check if models have been registered, and use default, in case this is not true @@
125
126 const G4DynamicParticle aPart( thePrimary.GetDefinition(), thePrimary.Get4Momentum().vect() );
127
128 if ( theQuasielastic )
129 {
130 if ( theQuasielastic->GetFraction( theNucleus, aPart ) > G4UniformRand() )
131 {
132 //G4cout<<"___G4TheoFSGenerator: before Scatter (1) QE=" << theQuasielastic<<G4endl;
133 G4KineticTrackVector* result= theQuasielastic->Scatter( theNucleus, aPart );
134 //G4cout << "^^G4TheoFSGenerator: after Scatter (1) " << G4endl;
135 if ( result )
136 {
137 for ( auto & ptr : *result )
138 {
139 G4DynamicParticle* aNew = new G4DynamicParticle( ptr->GetDefinition(),
140 ptr->Get4Momentum().e(),
141 ptr->Get4Momentum().vect() );
142 theParticleChange->AddSecondary( aNew, ptr->GetCreatorModelID() );
143 delete ptr;
144 }
145 delete result;
146 }
147 else
148 {
149 theParticleChange->SetStatusChange( isAlive );
150 theParticleChange->SetEnergyChange( thePrimary.GetKineticEnergy() );
151 theParticleChange->SetMomentumChange( thePrimary.Get4Momentum().vect().unit() );
152 }
153 return theParticleChange;
154 }
155 }
156
157 // get result from high energy model
158 G4KineticTrackVector* theInitialResult = theHighEnergyGenerator->Scatter( theNucleus, aPart );
159
160 // Assign the creator model ID.
161 // Only for the FTF string model, an inelastic interaction can be classified as quasi-elastic
162 // or diffractive. For the QGS model, the two boolean methods below return always "false".
163 G4int theModelID = theStringModelID;
164 if ( theHighEnergyGenerator->IsQuasiElasticInteraction() ) {
165 theModelID = theFTFQuasiElasticID;
166 } else if ( theHighEnergyGenerator->IsDiffractiveInteraction() ) {
167 theModelID = theFTFDiffractiveID;
168 }
169 for ( auto & ptr : *theInitialResult ) {
170 ptr->SetCreatorModelID( theModelID );
171 }
172
173 //#define DEBUG_initial_result
174 #ifdef DEBUG_initial_result
175 G4double E_out( 0 );
177 for ( auto & ptr : *theInitialResult )
178 {
179 //G4cout << "TheoFS secondary, mom " << ptr->GetDefinition()->GetParticleName()
180 // << " " << ptr->Get4Momentum() << G4endl;
181 E_out += ptr->Get4Momentum().e();
182 }
183 G4double init_mass = ionTable->GetIonMass( theNucleus.GetZ_asInt(), theNucleus.GetA_asInt() );
184 G4double init_E = aPart.Get4Momentum().e();
185 // residual nucleus
186
187 const std::vector< G4Nucleon >& thy = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
188
189 G4int resZ( 0 ), resA( 0 );
190 G4double delta_m( 0.0 );
191 for ( auto & nuc : thy )
192 {
193 if ( nuc.AreYouHit() ) {
194 ++resA;
195 if ( nuc.GetDefinition() == G4Proton::Proton() ) {
196 ++resZ;
197 delta_m += CLHEP::proton_mass_c2;
198 } else {
199 delta_m += CLHEP::neutron_mass_c2;
200 }
201 }
202 }
203
204 G4double final_mass( 0.0 );
205 if ( theNucleus.GetA_asInt() ) {
206 final_mass = ionTable->GetIonMass( theNucleus.GetZ_asInt() - resZ,theNucleus.GetA_asInt() - resA );
207 }
208 G4double E_excit = init_mass + init_E - final_mass - E_out;
209 G4cout << " Corrected delta mass " << init_mass - final_mass - delta_m << G4endl;
210 G4cout << "initial E, mass = " << init_E << ", " << init_mass << G4endl;
211 G4cout << " final E, mass = " << E_out <<", " << final_mass << " excitation_E " << E_excit << G4endl;
212 #endif
213
214 G4ReactionProductVector* theTransportResult = nullptr;
215
216 G4V3DNucleus* theProjectileNucleus = theHighEnergyGenerator->GetProjectileNucleus();
217 if ( theProjectileNucleus == nullptr )
218 {
219 G4int hitCount = 0;
220 const std::vector< G4Nucleon >& they = theHighEnergyGenerator->GetWoundedNucleus()->GetNucleons();
221 for ( auto & nuc : they )
222 {
223 if ( nuc.AreYouHit() ) ++hitCount;
224 }
225 if ( hitCount != theHighEnergyGenerator->GetWoundedNucleus()->GetMassNumber() )
226 {
227 theTransport->SetPrimaryProjectile( thePrimary );
228 theTransportResult = theTransport->Propagate( theInitialResult,
229 theHighEnergyGenerator->GetWoundedNucleus() );
230 if ( theTransportResult == nullptr ) {
231 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
232 throw G4HadronicException( __FILE__, __LINE__, "Null ptr from transport propagate" );
233 }
234 }
235 else
236 {
237 theTransportResult = theDecay.Propagate( theInitialResult,
238 theHighEnergyGenerator->GetWoundedNucleus() );
239 if ( theTransportResult == nullptr ) {
240 G4cout << "G4TheoFSGenerator: null ptr from decay propagate " << G4endl;
241 throw G4HadronicException( __FILE__, __LINE__, "Null ptr from decay propagate" );
242 }
243 }
244 }
245 else
246 {
247 theTransport->SetPrimaryProjectile( thePrimary );
248 theTransportResult = theTransport->PropagateNuclNucl( theInitialResult,
249 theHighEnergyGenerator->GetWoundedNucleus(),
250 theProjectileNucleus );
251 if ( theTransportResult == nullptr ) {
252 G4cout << "G4TheoFSGenerator: null ptr from transport propagate " << G4endl;
253 throw G4HadronicException( __FILE__, __LINE__, "Null ptr from transport propagate" );
254 }
255 }
256
257 // If enabled, apply the Cosmic Rays (CR) coalescence to the list of secondaries produced so far.
258 // This algorithm can form deuterons and antideuterons by coalescence of, respectively,
259 // proton-neutron and antiproton-antineutron pairs close in momentum space.
260 // This can be useful in particular for Cosmic Ray applications.
261 if ( G4HadronicParameters::Instance()->EnableCRCoalescence() ) {
262 if ( theCosmicCoalescence == nullptr ) {
263 theCosmicCoalescence = (G4CRCoalescence*)
264 G4HadronicInteractionRegistry::Instance()->FindModel( "G4CRCoalescence" );
265 if ( theCosmicCoalescence == nullptr ) {
266 theCosmicCoalescence = new G4CRCoalescence;
267 }
268 }
269 theCosmicCoalescence->SetP0Coalescence( thePrimary, theHighEnergyGenerator->GetModelName() );
270 theCosmicCoalescence->GenerateDeuterons( theTransportResult );
271 }
272
273 // Fill particle change
274 for ( auto & ptr : *theTransportResult )
275 {
276 G4DynamicParticle* aNewDP = new G4DynamicParticle( ptr->GetDefinition(),
277 ptr->GetTotalEnergy(),
278 ptr->GetMomentum() );
279 G4HadSecondary aNew = G4HadSecondary( aNewDP );
280 G4double time = std::max( ptr->GetFormationTime(), 0.0 );
281 aNew.SetTime( timePrimary + time );
282 aNew.SetCreatorModelID( ptr->GetCreatorModelID() );
283 aNew.SetParentResonanceDef( ptr->GetParentResonanceDef() );
284 aNew.SetParentResonanceID( ptr->GetParentResonanceID() );
285 theParticleChange->AddSecondary( aNew );
286 delete ptr;
287 }
288
289 // some garbage collection
290 delete theTransportResult;
291
292 // Done
293 return theParticleChange;
294}
295
296
297std::pair<G4double, G4double> G4TheoFSGenerator::GetEnergyMomentumCheckLevels() const
298{
299 if ( theHighEnergyGenerator ) {
300 return theHighEnergyGenerator->GetEnergyMomentumCheckLevels();
301 } else {
302 return std::pair<G4double, G4double>( DBL_MAX, DBL_MAX );
303 }
304}
@ isAlive
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
G4LorentzVector Get4Momentum() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetParentResonanceDef(const G4ParticleDefinition *parentDef)
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
void SetParentResonanceID(const G4int parentID)
static G4HadronicInteractionRegistry * Instance()
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
static G4HadronicParameters * Instance()
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
G4int GetA_asInt() const
Definition G4Nucleus.hh:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
G4int GetQuarkContent(G4int flavor) const
G4bool IsHypernucleus() const
G4int GetAntiQuarkContent(G4int flavor) const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition G4Proton.cc:90
G4TheoFSGenerator(const G4String &name="TheoFSGenerator")
std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const override
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override
void ModelDescription(std::ostream &outFile) const override
~G4TheoFSGenerator() override
#define DBL_MAX
Definition templates.hh:62