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