Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FTFModel Class Reference

#include <G4FTFModel.hh>

+ Inheritance diagram for G4FTFModel:

Public Member Functions

 G4FTFModel (const G4String &modelName="FTF")
 
 ~G4FTFModel () override
 
G4V3DNucleusGetTargetNucleus () const
 
G4V3DNucleusGetWoundedNucleus () const override
 
G4V3DNucleusGetProjectileNucleus () const override
 
void ModelDescription (std::ostream &) const override
 
 G4FTFModel (const G4FTFModel &right)=delete
 
const G4FTFModeloperator= (const G4FTFModel &right)=delete
 
G4bool operator== (const G4FTFModel &right) const =delete
 
G4bool operator!= (const G4FTFModel &right) const =delete
 
void SetImpactParameter (const G4double b_value)
 
G4double GetImpactParameter () const
 
void SetBminBmax (const G4double bmin_value, const G4double bmax_value)
 
G4bool SampleBinInterval () const
 
G4double GetBmin () const
 
G4double GetBmax () const
 
G4int GetNumberOfProjectileSpectatorNucleons () const
 
G4int GetNumberOfTargetSpectatorNucleons () const
 
G4int GetNumberOfNNcollisions () const
 
- Public Member Functions inherited from G4VPartonStringModel
 G4VPartonStringModel (const G4String &modelName="Parton String Model")
 
 ~G4VPartonStringModel () override
 
 G4VPartonStringModel (const G4VPartonStringModel &right)=delete
 
const G4VPartonStringModeloperator= (const G4VPartonStringModel &right)=delete
 
G4bool operator== (const G4VPartonStringModel &right) const =delete
 
G4bool operator!= (const G4VPartonStringModel &right) const =delete
 
void SetFragmentationModel (G4VStringFragmentation *aModel)
 
G4KineticTrackVectorScatter (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary) override
 
void ModelDescription (std::ostream &outFile) const override
 
G4V3DNucleusGetProjectileNucleus () const override
 
- Public Member Functions inherited from G4VHighEnergyGenerator
 G4VHighEnergyGenerator (const G4String &modelName="HighEnergyGenerator")
 
 ~G4VHighEnergyGenerator () override
 
 G4VHighEnergyGenerator (const G4VHighEnergyGenerator &right)=delete
 
const G4VHighEnergyGeneratoroperator= (const G4VHighEnergyGenerator &right)=delete
 
G4bool operator== (const G4VHighEnergyGenerator &right) const =delete
 
G4bool operator!= (const G4VHighEnergyGenerator &right) const =delete
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Protected Member Functions

void Init (const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile) override
 
G4ExcitedStringVectorGetStrings () override
 
- Protected Member Functions inherited from G4VPartonStringModel
G4bool EnergyAndMomentumCorrector (G4KineticTrackVector *Output, G4LorentzVector &TotalCollisionMomentum)
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 

Additional Inherited Members

- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 61 of file G4FTFModel.hh.

Constructor & Destructor Documentation

◆ G4FTFModel() [1/2]

G4FTFModel::G4FTFModel ( const G4String & modelName = "FTF")

Definition at line 70 of file G4FTFModel.cc.

70 :
71 G4VPartonStringModel( modelName ),
72 theExcitation( new G4DiffractiveExcitation() ),
73 theElastic( new G4ElasticHNScattering() ),
74 theAnnihilation( new G4FTFAnnihilation() )
75{
76 // ---> JVY theParameters = 0;
77 theParameters = new G4FTFParameters();
78 //
79 NumberOfInvolvedNucleonsOfTarget = 0;
80 NumberOfInvolvedNucleonsOfProjectile= 0;
81 for ( G4int i = 0; i < 250; ++i ) {
82 TheInvolvedNucleonsOfTarget[i] = 0;
83 TheInvolvedNucleonsOfProjectile[i] = 0;
84 }
85
86 //LowEnergyLimit = 2000.0*MeV;
87 LowEnergyLimit = 1000.0*MeV;
88
89 HighEnergyInter = true;
90
91 G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
92 ProjectileResidual4Momentum = tmp;
93 ProjectileResidualMassNumber = 0;
94 ProjectileResidualCharge = 0;
95 ProjectileResidualLambdaNumber = 0;
96 ProjectileResidualExcitationEnergy = 0.0;
97
98 TargetResidual4Momentum = tmp;
99 TargetResidualMassNumber = 0;
100 TargetResidualCharge = 0;
101 TargetResidualExcitationEnergy = 0.0;
102
103 Bimpact = -1.0;
104 BinInterval = false;
105 Bmin = 0.0;
106 Bmax = 0.0;
107 NumberOfProjectileSpectatorNucleons = 0;
108 NumberOfTargetSpectatorNucleons = 0;
109 NumberOfNNcollisions = 0;
110
111 SetEnergyMomentumCheckLevels( 2.0*perCent, 150.0*MeV );
112}
int G4int
Definition G4Types.hh:85
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
G4VPartonStringModel(const G4String &modelName="Parton String Model")

◆ ~G4FTFModel()

G4FTFModel::~G4FTFModel ( )
override

Definition at line 122 of file G4FTFModel.cc.

122 {
123 // Because FTF model can be called for various particles
124 //
125 // ---> NOTE (JVY): This statement below is no longer true !!!
126 // theParameters must be erased at the end of each call.
127 // Thus the delete is also in G4FTFModel::GetStrings() method.
128 // ---> JVY
129 //
130 if ( theParameters != 0 ) delete theParameters;
131 if ( theExcitation != 0 ) delete theExcitation;
132 if ( theElastic != 0 ) delete theElastic;
133 if ( theAnnihilation != 0 ) delete theAnnihilation;
134
135 // Erasing of strings created at annihilation.
136 if ( theAdditionalString.size() != 0 ) {
137 std::for_each( theAdditionalString.begin(), theAdditionalString.end(),
139 }
140 theAdditionalString.clear();
141
142 // Erasing of target involved nucleons.
143 if ( NumberOfInvolvedNucleonsOfTarget != 0 ) {
144 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; ++i ) {
145 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfTarget[i]->GetSplitableHadron();
146 if ( aNucleon ) delete aNucleon;
147 }
148 }
149
150 // Erasing of projectile involved nucleons.
151 if ( NumberOfInvolvedNucleonsOfProjectile != 0 ) {
152 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfProjectile; ++i ) {
153 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfProjectile[i]->GetSplitableHadron();
154 if ( aNucleon ) delete aNucleon;
155 }
156 }
157}
G4VSplitableHadron * GetSplitableHadron() const
Definition G4Nucleon.hh:97

◆ G4FTFModel() [2/2]

G4FTFModel::G4FTFModel ( const G4FTFModel & right)
delete

Member Function Documentation

◆ GetBmax()

G4double G4FTFModel::GetBmax ( ) const
inline

Definition at line 248 of file G4FTFModel.hh.

248 {
249 return Bmax;
250}

Referenced by Init().

◆ GetBmin()

G4double G4FTFModel::GetBmin ( ) const
inline

Definition at line 244 of file G4FTFModel.hh.

244 {
245 return Bmin;
246}

Referenced by Init().

◆ GetImpactParameter()

G4double G4FTFModel::GetImpactParameter ( ) const
inline

Definition at line 228 of file G4FTFModel.hh.

228 {
229 return Bimpact;
230}

◆ GetNumberOfNNcollisions()

G4int G4FTFModel::GetNumberOfNNcollisions ( ) const
inline

Definition at line 260 of file G4FTFModel.hh.

260 {
261 return NumberOfNNcollisions;
262}

◆ GetNumberOfProjectileSpectatorNucleons()

G4int G4FTFModel::GetNumberOfProjectileSpectatorNucleons ( ) const
inline

Definition at line 252 of file G4FTFModel.hh.

252 {
253 return NumberOfProjectileSpectatorNucleons;
254}

◆ GetNumberOfTargetSpectatorNucleons()

G4int G4FTFModel::GetNumberOfTargetSpectatorNucleons ( ) const
inline

Definition at line 256 of file G4FTFModel.hh.

256 {
257 return NumberOfTargetSpectatorNucleons;
258}

◆ GetProjectileNucleus()

G4V3DNucleus * G4FTFModel::GetProjectileNucleus ( ) const
inlineoverridevirtual

Reimplemented from G4VHighEnergyGenerator.

Definition at line 220 of file G4FTFModel.hh.

220 {
221 return theParticipants.GetProjectileNucleus();
222}
virtual G4V3DNucleus * GetProjectileNucleus() const

Referenced by GetStrings().

◆ GetStrings()

G4ExcitedStringVector * G4FTFModel::GetStrings ( )
overrideprotectedvirtual

Implements G4VPartonStringModel.

Definition at line 298 of file G4FTFModel.cc.

298 {
299
300 #ifdef debugFTFmodel
301 G4cout << "G4FTFModel::GetStrings() " << G4endl;
302 #endif
303
305 theParticipants.GetList( theProjectile, theParameters );
306
307 SetImpactParameter( theParticipants.GetImpactParameter() );
308
309 StoreInvolvedNucleon();
310
311 G4bool Success( true );
312
313 if ( HighEnergyInter ) {
314 ReggeonCascade();
315
316 #ifdef debugFTFmodel
317 G4cout << "FTF PutOnMassShell " << G4endl;
318 #endif
319
320 Success = PutOnMassShell();
321
322 #ifdef debugFTFmodel
323 G4cout << "FTF PutOnMassShell Success? " << Success << G4endl;
324 #endif
325
326 }
327
328 #ifdef debugFTFmodel
329 G4cout << "FTF ExciteParticipants " << G4endl;
330 #endif
331
332 if ( Success ) Success = ExciteParticipants();
333
334 #ifdef debugFTFmodel
335 G4cout << "FTF ExciteParticipants Success? " << Success << G4endl;
336 #endif
337
338 if ( Success ) {
339
340 #ifdef debugFTFmodel
341 G4cout << "FTF BuildStrings ";
342 #endif
343
344 BuildStrings( theStrings );
345
346 #ifdef debugFTFmodel
347 G4cout << "FTF BuildStrings " << theStrings << " OK" << G4endl
348 << "FTF GetResiduals of Nuclei " << G4endl;
349 #endif
350
351 GetResiduals();
352
353 /*
354 if ( theParameters != 0 ) {
355 delete theParameters;
356 theParameters = 0;
357 }
358 */
359 } else if ( ! GetProjectileNucleus() ) {
360 // Erase the hadron projectile
361 std::vector< G4VSplitableHadron* > primaries;
362 theParticipants.StartLoop();
363 while ( theParticipants.Next() ) { /* Loop checking, 10.08.2015, A.Ribon */
364 const G4InteractionContent& interaction = theParticipants.GetInteraction();
365 // Do not allow for duplicates
366 if ( primaries.end() ==
367 std::find( primaries.begin(), primaries.end(), interaction.GetProjectile() ) ) {
368 primaries.push_back( interaction.GetProjectile() );
369 }
370 }
371 std::for_each( primaries.begin(), primaries.end(), DeleteVSplitableHadron() );
372 primaries.clear();
373 }
374
375 // Cleaning of the memory
376 G4VSplitableHadron* aNucleon = 0;
377
378 // Erase the projectile nucleons
379 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfProjectile; ++i ) {
380 aNucleon = TheInvolvedNucleonsOfProjectile[i]->GetSplitableHadron();
381 if ( aNucleon ) delete aNucleon;
382 }
383 NumberOfInvolvedNucleonsOfProjectile = 0;
384
385 // Erase the target nucleons
386 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; ++i ) {
387 aNucleon = TheInvolvedNucleonsOfTarget[i]->GetSplitableHadron();
388 if ( aNucleon ) delete aNucleon;
389 }
390 NumberOfInvolvedNucleonsOfTarget = 0;
391
392 #ifdef debugFTFmodel
393 G4cout << "End of FTF. Go to fragmentation" << G4endl
394 << "To continue - enter 1, to stop - ^C" << G4endl;
395 #endif
396
397 theParticipants.Clean();
398
399 return theStrings;
400}
std::vector< G4ExcitedString * > G4ExcitedStringVector
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetImpactParameter(const G4double b_value)
G4V3DNucleus * GetProjectileNucleus() const override
void GetList(const G4ReactionProduct &thePrimary, G4FTFParameters *theParameters)
G4double GetImpactParameter() const
G4InteractionContent & GetInteraction()
G4VSplitableHadron * GetProjectile() const

◆ GetTargetNucleus()

G4V3DNucleus * G4FTFModel::GetTargetNucleus ( ) const
inline

Definition at line 216 of file G4FTFModel.hh.

216 {
217 return theParticipants.GetWoundedNucleus();
218}
virtual G4V3DNucleus * GetWoundedNucleus() const

◆ GetWoundedNucleus()

G4V3DNucleus * G4FTFModel::GetWoundedNucleus ( ) const
inlineoverridevirtual

Implements G4VHighEnergyGenerator.

Definition at line 212 of file G4FTFModel.hh.

212 {
213 return theParticipants.GetWoundedNucleus();
214}

◆ Init()

void G4FTFModel::Init ( const G4Nucleus & aNucleus,
const G4DynamicParticle & aProjectile )
overrideprotectedvirtual

Implements G4VPartonStringModel.

Definition at line 162 of file G4FTFModel.cc.

162 {
163
164 theProjectile = aProjectile;
165
166 G4double PlabPerParticle( 0.0 ); // Laboratory momentum Pz per particle/nucleon
167
168 #ifdef debugFTFmodel
169 G4cout << "FTF init Proj Name " << theProjectile.GetDefinition()->GetParticleName() << G4endl
170 << "FTF init Proj Mass " << theProjectile.GetMass()
171 << " " << theProjectile.GetMomentum() << G4endl
172 << "FTF init Proj B Q " << theProjectile.GetDefinition()->GetBaryonNumber()
173 << " " << (G4int) theProjectile.GetDefinition()->GetPDGCharge() << G4endl
174 << "FTF init Target A Z " << aNucleus.GetA_asInt()
175 << " " << aNucleus.GetZ_asInt() << G4endl;
176 #endif
177
178 theParticipants.Clean();
179
180 theParticipants.SetProjectileNucleus( 0 );
181
182 G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
183 ProjectileResidualMassNumber = 0;
184 ProjectileResidualCharge = 0;
185 ProjectileResidualLambdaNumber = 0;
186 ProjectileResidualExcitationEnergy = 0.0;
187 ProjectileResidual4Momentum = tmp;
188
189 TargetResidualMassNumber = aNucleus.GetA_asInt();
190 TargetResidualCharge = aNucleus.GetZ_asInt();
191 TargetResidualExcitationEnergy = 0.0;
192 TargetResidual4Momentum = tmp;
194 ->GetIonMass( TargetResidualCharge, TargetResidualMassNumber );
195
196 TargetResidual4Momentum.setE( TargetResidualMass );
197
198 if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
199 // Projectile is a hadron : meson or baryon
200 ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
201 ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
202 PlabPerParticle = theProjectile.GetMomentum().z();
203 ProjectileResidualExcitationEnergy = 0.0;
204 //G4double ProjectileResidualMass = theProjectile.GetMass();
205 ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
206 ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
207 if ( PlabPerParticle < LowEnergyLimit ) {
208 HighEnergyInter = false;
209 } else {
210 HighEnergyInter = true;
211 }
212 } else {
213 if ( theProjectile.GetDefinition()->GetBaryonNumber() > 1 ) {
214 // Projectile is a nucleus
215 ProjectileResidualMassNumber = theProjectile.GetDefinition()->GetBaryonNumber();
216 ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
217 ProjectileResidualLambdaNumber = theProjectile.GetDefinition()->GetNumberOfLambdasInHypernucleus();
218 PlabPerParticle = theProjectile.GetMomentum().z() / ProjectileResidualMassNumber;
219 if ( PlabPerParticle < LowEnergyLimit ) {
220 HighEnergyInter = false;
221 } else {
222 HighEnergyInter = true;
223 }
224 theParticipants.InitProjectileNucleus( ProjectileResidualMassNumber, ProjectileResidualCharge,
225 ProjectileResidualLambdaNumber );
226 } else if ( theProjectile.GetDefinition()->GetBaryonNumber() < -1 ) {
227 // Projectile is an anti-nucleus
228 ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
229 ProjectileResidualCharge = std::abs( G4int( theProjectile.GetDefinition()->GetPDGCharge() ) );
230 ProjectileResidualLambdaNumber = theProjectile.GetDefinition()->GetNumberOfAntiLambdasInAntiHypernucleus();
231 PlabPerParticle = theProjectile.GetMomentum().z() / ProjectileResidualMassNumber;
232 if ( PlabPerParticle < LowEnergyLimit ) {
233 HighEnergyInter = false;
234 } else {
235 HighEnergyInter = true;
236 }
237 theParticipants.InitProjectileNucleus( ProjectileResidualMassNumber, ProjectileResidualCharge,
238 ProjectileResidualLambdaNumber );
239 theParticipants.GetProjectileNucleus()->StartLoop();
240 G4Nucleon* aNucleon;
241 while ( ( aNucleon = theParticipants.GetProjectileNucleus()->GetNextNucleon() ) ) { /* Loop checking, 10.08.2015, A.Ribon */
242 if ( aNucleon->GetDefinition() == G4Proton::Definition() ) {
244 } else if ( aNucleon->GetDefinition() == G4Neutron::Definition() ) {
246 } else if ( aNucleon->GetDefinition() == G4Lambda::Definition() ) {
248 }
249 }
250 }
251
252 G4ThreeVector BoostVector = theProjectile.GetMomentum() / theProjectile.GetTotalEnergy();
253 theParticipants.GetProjectileNucleus()->DoLorentzBoost( BoostVector );
254 theParticipants.GetProjectileNucleus()->DoLorentzContraction( BoostVector );
255 ProjectileResidualExcitationEnergy = 0.0;
256 //G4double ProjectileResidualMass = theProjectile.GetMass();
257 ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
258 ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
259 }
260
261 // Init target nucleus (assumed to be never a hypernucleus)
262 theParticipants.Init( aNucleus.GetA_asInt(), aNucleus.GetZ_asInt() );
263
264 NumberOfProjectileSpectatorNucleons = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
265 NumberOfTargetSpectatorNucleons = aNucleus.GetA_asInt();
266 NumberOfNNcollisions = 0;
267
268 // reset/recalculate everything for the new interaction
269 theParameters->InitForInteraction( theProjectile.GetDefinition(), aNucleus.GetA_asInt(),
270 aNucleus.GetZ_asInt(), PlabPerParticle );
271
272 if ( theAdditionalString.size() != 0 ) {
273 std::for_each( theAdditionalString.begin(), theAdditionalString.end(),
275 }
276 theAdditionalString.clear();
277
278 #ifdef debugFTFmodel
279 G4cout << "FTF end of Init" << G4endl << G4endl;
280 #endif
281
282 // In the case of Hydrogen target, for non-ion hadron projectiles,
283 // do NOT simulate quasi-elastic (by forcing to 0 the probability of
284 // elastic scatering in theParameters - which is used only by FTF).
285 // This is necessary because in this case quasi-elastic on a target nucleus
286 // with only one nucleon would be identical to the hadron elastic scattering,
287 // and the latter is already included in the elastic process
288 // (i.e. G4HadronElasticProcess).
289 if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 &&
290 aNucleus.GetA_asInt() < 2 ) theParameters->SetProbabilityOfElasticScatt( 0.0 );
291
292 if ( SampleBinInterval() ) theParticipants.SetBminBmax( GetBmin(), GetBmax() );
293}
double G4double
Definition G4Types.hh:83
double z() const
void setVect(const Hep3Vector &)
static G4AntiLambda * Definition()
static G4AntiNeutron * Definition()
static G4AntiProton * Definition()
G4bool SampleBinInterval() const
G4double GetBmin() const
G4double GetBmax() const
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void InitForInteraction(const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
void SetBminBmax(const G4double bmin_value, const G4double bmax_value)
G4double GetIonMass(G4int Z, G4int A, G4int nL=0, G4int lvl=0) const
static G4Lambda * Definition()
Definition G4Lambda.cc:48
static G4Neutron * Definition()
Definition G4Neutron.cc:50
void SetParticleType(G4Proton *aProton)
Definition G4Nucleon.hh:77
virtual const G4ParticleDefinition * GetDefinition() const
Definition G4Nucleon.hh:86
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
G4int GetNumberOfLambdasInHypernucleus() const
G4int GetNumberOfAntiLambdasInAntiHypernucleus() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Definition()
Definition G4Proton.cc:45
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetMass() const
virtual G4Nucleon * GetNextNucleon()=0
virtual G4bool StartLoop()=0
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
virtual void DoLorentzContraction(const G4LorentzVector &theBoost)=0
virtual void InitProjectileNucleus(G4int theZ, G4int theA, G4int numberOfLambdasOrAntiLambdas=0)
virtual void Init(G4int theZ, G4int theA)
virtual void SetProjectileNucleus(G4V3DNucleus *aNucleus)

◆ ModelDescription()

void G4FTFModel::ModelDescription ( std::ostream & desc) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 3151 of file G4FTFModel.cc.

3151 {
3152 desc << " FTF (Fritiof) Model \n"
3153 << "The FTF model is based on the well-known FRITIOF \n"
3154 << "model (B. Andersson et al., Nucl. Phys. B281, 289 \n"
3155 << "(1987)). Its first program implementation was given\n"
3156 << "by B. Nilsson-Almquist and E. Stenlund (Comp. Phys.\n"
3157 << "Comm. 43, 387 (1987)). The Fritiof model assumes \n"
3158 << "that all hadron-hadron interactions are binary \n"
3159 << "reactions, h_1+h_2->h_1'+h_2' where h_1' and h_2' \n"
3160 << "are excited states of the hadrons with continuous \n"
3161 << "mass spectra. The excited hadrons are considered as\n"
3162 << "QCD-strings, and the corresponding LUND-string \n"
3163 << "fragmentation model is applied for a simulation of \n"
3164 << "their decays. \n"
3165 << " The Fritiof model assumes that in the course of \n"
3166 << "a hadron-nucleus interaction a string originated \n"
3167 << "from the projectile can interact with various intra\n"
3168 << "nuclear nucleons and becomes into highly excited \n"
3169 << "states. The probability of multiple interactions is\n"
3170 << "calculated in the Glauber approximation. A cascading\n"
3171 << "of secondary particles was neglected as a rule. Due\n"
3172 << "to these, the original Fritiof model fails to des- \n"
3173 << "cribe a nuclear destruction and slow particle spectra.\n"
3174 << " In order to overcome the difficulties we enlarge\n"
3175 << "the model by the reggeon theory inspired model of \n"
3176 << "nuclear desctruction (Kh. Abdel-Waged and V.V. Uzhi-\n"
3177 << "nsky, Phys. Atom. Nucl. 60, 828 (1997); Yad. Fiz. 60, 925\n"
3178 << "(1997)). Momenta of the nucleons ejected from a nuc-\n"
3179 << "leus in the reggeon cascading are sampled according\n"
3180 << "to a Fermi motion algorithm presented in (EMU-01 \n"
3181 << "Collaboration (M.I. Adamovich et al.) Zeit. fur Phys.\n"
3182 << "A358, 337 (1997)). \n"
3183 << " New features were also added to the Fritiof model\n"
3184 << "implemented in Geant4: a simulation of elastic had-\n"
3185 << "ron-nucleon scatterings, a simulation of binary \n"
3186 << "reactions like NN>NN* in hadron-nucleon interactions,\n"
3187 << "a separate simulation of single diffractive and non-\n"
3188 << " diffractive events. These allowed to describe after\n"
3189 << "model parameter tuning a wide set of experimental \n"
3190 << "data. \n";
3191}

◆ operator!=()

G4bool G4FTFModel::operator!= ( const G4FTFModel & right) const
delete

◆ operator=()

const G4FTFModel & G4FTFModel::operator= ( const G4FTFModel & right)
delete

◆ operator==()

G4bool G4FTFModel::operator== ( const G4FTFModel & right) const
delete

◆ SampleBinInterval()

G4bool G4FTFModel::SampleBinInterval ( ) const
inline

Definition at line 240 of file G4FTFModel.hh.

240 {
241 return BinInterval;
242}

Referenced by Init().

◆ SetBminBmax()

void G4FTFModel::SetBminBmax ( const G4double bmin_value,
const G4double bmax_value )
inline

Definition at line 232 of file G4FTFModel.hh.

232 {
233 BinInterval = false;
234 if ( bmin_value < 0.0 || bmax_value < 0.0 || bmax_value < bmin_value ) return;
235 BinInterval = true;
236 Bmin = bmin_value;
237 Bmax = bmax_value;
238}

◆ SetImpactParameter()

void G4FTFModel::SetImpactParameter ( const G4double b_value)
inline

Definition at line 224 of file G4FTFModel.hh.

224 {
225 Bimpact = b_value;
226}

Referenced by GetStrings().


The documentation for this class was generated from the following files: