Geant4 10.7.0
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
 
- 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
 
virtual G4V3DNucleusGetWoundedNucleus () const =0
 
virtual G4V3DNucleusGetProjectileNucleus () const
 
virtual G4KineticTrackVectorScatter (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
 
void ModelDescription (std::ostream &) const override
 
- 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 ModelDescription (std::ostream &outFile) const
 
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
virtual void Init (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
 
virtual G4ExcitedStringVectorGetStrings ()=0
 
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 ProjectileResidualExcitationEnergy = 0.0;
96
97 TargetResidual4Momentum = tmp;
98 TargetResidualMassNumber = 0;
99 TargetResidualCharge = 0;
100 TargetResidualExcitationEnergy = 0.0;
101
102 SetEnergyMomentumCheckLevels( 2.0*perCent, 150.0*MeV );
103}
int G4int
Definition: G4Types.hh:85
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)

◆ ~G4FTFModel()

G4FTFModel::~G4FTFModel ( )
override

Definition at line 113 of file G4FTFModel.cc.

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

◆ G4FTFModel() [2/2]

G4FTFModel::G4FTFModel ( const G4FTFModel right)
delete

Member Function Documentation

◆ GetProjectileNucleus()

G4V3DNucleus * G4FTFModel::GetProjectileNucleus ( ) const
inlineoverridevirtual

Reimplemented from G4VHighEnergyGenerator.

Definition at line 201 of file G4FTFModel.hh.

201 {
202 return theParticipants.GetProjectileNucleus();
203}
virtual G4V3DNucleus * GetProjectileNucleus() const

Referenced by GetStrings().

◆ GetStrings()

G4ExcitedStringVector * G4FTFModel::GetStrings ( )
overrideprotectedvirtual

Implements G4VPartonStringModel.

Definition at line 290 of file G4FTFModel.cc.

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

◆ GetTargetNucleus()

G4V3DNucleus * G4FTFModel::GetTargetNucleus ( ) const
inline

Definition at line 197 of file G4FTFModel.hh.

197 {
198 return theParticipants.GetWoundedNucleus();
199}
virtual G4V3DNucleus * GetWoundedNucleus() const

◆ GetWoundedNucleus()

G4V3DNucleus * G4FTFModel::GetWoundedNucleus ( ) const
inlineoverridevirtual

Implements G4VHighEnergyGenerator.

Definition at line 193 of file G4FTFModel.hh.

193 {
194 return theParticipants.GetWoundedNucleus();
195}

◆ Init()

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

Implements G4VPartonStringModel.

Definition at line 153 of file G4FTFModel.cc.

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

◆ ModelDescription()

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

Reimplemented from G4VHighEnergyGenerator.

Definition at line 3069 of file G4FTFModel.cc.

3069 {
3070 desc << " FTF (Fritiof) Model \n"
3071 << "The FTF model is based on the well-known FRITIOF \n"
3072 << "model (B. Andersson et al., Nucl. Phys. B281, 289 \n"
3073 << "(1987)). Its first program implementation was given\n"
3074 << "by B. Nilsson-Almquist and E. Stenlund (Comp. Phys.\n"
3075 << "Comm. 43, 387 (1987)). The Fritiof model assumes \n"
3076 << "that all hadron-hadron interactions are binary \n"
3077 << "reactions, h_1+h_2->h_1'+h_2' where h_1' and h_2' \n"
3078 << "are excited states of the hadrons with continuous \n"
3079 << "mass spectra. The excited hadrons are considered as\n"
3080 << "QCD-strings, and the corresponding LUND-string \n"
3081 << "fragmentation model is applied for a simulation of \n"
3082 << "their decays. \n"
3083 << " The Fritiof model assumes that in the course of \n"
3084 << "a hadron-nucleus interaction a string originated \n"
3085 << "from the projectile can interact with various intra\n"
3086 << "nuclear nucleons and becomes into highly excited \n"
3087 << "states. The probability of multiple interactions is\n"
3088 << "calculated in the Glauber approximation. A cascading\n"
3089 << "of secondary particles was neglected as a rule. Due\n"
3090 << "to these, the original Fritiof model fails to des- \n"
3091 << "cribe a nuclear destruction and slow particle spectra.\n"
3092 << " In order to overcome the difficulties we enlarge\n"
3093 << "the model by the reggeon theory inspired model of \n"
3094 << "nuclear desctruction (Kh. Abdel-Waged and V.V. Uzhi-\n"
3095 << "nsky, Phys. Atom. Nucl. 60, 828 (1997); Yad. Fiz. 60, 925\n"
3096 << "(1997)). Momenta of the nucleons ejected from a nuc-\n"
3097 << "leus in the reggeon cascading are sampled according\n"
3098 << "to a Fermi motion algorithm presented in (EMU-01 \n"
3099 << "Collaboration (M.I. Adamovich et al.) Zeit. fur Phys.\n"
3100 << "A358, 337 (1997)). \n"
3101 << " New features were also added to the Fritiof model\n"
3102 << "implemented in Geant4: a simulation of elastic had-\n"
3103 << "ron-nucleon scatterings, a simulation of binary \n"
3104 << "reactions like NN>NN* in hadron-nucleon interactions,\n"
3105 << "a separate simulation of single diffractive and non-\n"
3106 << " diffractive events. These allowed to describe after\n"
3107 << "model parameter tuning a wide set of experimental \n"
3108 << "data. \n";
3109}

◆ 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

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