Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FTFModel Class Reference

#include <G4FTFModel.hh>

Inheritance diagram for G4FTFModel:
Collaboration diagram for G4FTFModel:

Public Member Functions

 G4FTFModel (const G4String &modelName="FTF")
 
 ~G4FTFModel ()
 
void Init (const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile)
 
G4ExcitedStringVectorGetStrings ()
 
G4V3DNucleusGetWoundedNucleus () const
 
G4V3DNucleusGetTargetNucleus () const
 
G4V3DNucleusGetProjectileNucleus () const
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VPartonStringModel
 G4VPartonStringModel (const G4String &modelName="Parton String Model")
 
virtual ~G4VPartonStringModel ()
 
void SetFragmentationModel (G4VStringFragmentation *aModel)
 
G4KineticTrackVectorScatter (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
- Public Member Functions inherited from G4VHighEnergyGenerator
 G4VHighEnergyGenerator (const G4String &modelName="High Energy Generator")
 
virtual ~G4VHighEnergyGenerator ()
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double AbsoluteLevel)
 
virtual G4String GetModelName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VPartonStringModel
void SetThisPointer (G4VPartonStringModel *aPointer)
 
G4bool EnergyAndMomentumCorrector (G4KineticTrackVector *Output, G4LorentzVector &TotalCollisionMomentum)
 

Detailed Description

Definition at line 63 of file G4FTFModel.hh.

Constructor & Destructor Documentation

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

Definition at line 72 of file G4FTFModel.cc.

72  :
73  G4VPartonStringModel( modelName ),
74  theExcitation( new G4DiffractiveExcitation() ),
75  theElastic( new G4ElasticHNScattering() ),
76  theAnnihilation( new G4FTFAnnihilation() )
77 {
79  theParameters = 0;
80  NumberOfInvolvedNucleonsOfTarget = 0;
81  NumberOfInvolvedNucleonsOfProjectile= 0;
82  for ( G4int i = 0; i < 250; i++ ) {
83  TheInvolvedNucleonsOfTarget[i] = 0;
84  TheInvolvedNucleonsOfProjectile[i] = 0;
85  }
86 
87  //LowEnergyLimit = 2000.0*MeV;
88  LowEnergyLimit = 1000.0*MeV;
89 
90  HighEnergyInter = true;
91 
92  G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
93  ProjectileResidual4Momentum = tmp;
94  ProjectileResidualMassNumber = 0;
95  ProjectileResidualCharge = 0;
96  ProjectileResidualExcitationEnergy = 0.0;
97 
98  TargetResidual4Momentum = tmp;
99  TargetResidualMassNumber = 0;
100  TargetResidualCharge = 0;
101  TargetResidualExcitationEnergy = 0.0;
102 
104 }
G4VPartonStringModel(const G4String &modelName="Parton String Model")
static constexpr double perCent
Definition: G4SIunits.hh:332
int G4int
Definition: G4Types.hh:78
void SetThisPointer(G4VPartonStringModel *aPointer)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double AbsoluteLevel)
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

G4FTFModel::~G4FTFModel ( )

Definition at line 114 of file G4FTFModel.cc.

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

Here is the call graph for this function:

Member Function Documentation

G4V3DNucleus * G4FTFModel::GetProjectileNucleus ( ) const
inlinevirtual

Reimplemented from G4VPartonStringModel.

Definition at line 170 of file G4FTFModel.hh.

170  {
171  return theParticipants.GetProjectileNucleus();
172 }
virtual G4V3DNucleus * GetProjectileNucleus() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4ExcitedStringVector * G4FTFModel::GetStrings ( )
virtual

Implements G4VPartonStringModel.

Definition at line 272 of file G4FTFModel.cc.

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

Here is the call graph for this function:

G4V3DNucleus * G4FTFModel::GetTargetNucleus ( ) const
inline

Definition at line 165 of file G4FTFModel.hh.

165  {
166  return theParticipants.GetWoundedNucleus();
167 }
virtual G4V3DNucleus * GetWoundedNucleus() const

Here is the call graph for this function:

G4V3DNucleus * G4FTFModel::GetWoundedNucleus ( ) const
inlinevirtual

Implements G4VPartonStringModel.

Definition at line 160 of file G4FTFModel.hh.

160  {
161  return theParticipants.GetWoundedNucleus();
162 }
virtual G4V3DNucleus * GetWoundedNucleus() const

Here is the call graph for this function:

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

Implements G4VPartonStringModel.

Definition at line 150 of file G4FTFModel.cc.

150  {
151 
152  theProjectile = aProjectile;
153 
154  G4double PlabPerParticle( 0.0 ); // Laboratory momentum Pz per particle/nucleon
155 
156  #ifdef debugFTFmodel
157  G4cout << "FTF init Proj Name " << theProjectile.GetDefinition()->GetParticleName() << G4endl
158  << "FTF init Proj Mass " << theProjectile.GetMass()
159  << " " << theProjectile.GetMomentum() << G4endl
160  << "FTF init Proj B Q " << theProjectile.GetDefinition()->GetBaryonNumber()
161  << " " << (G4int) theProjectile.GetDefinition()->GetPDGCharge() << G4endl
162  << "FTF init Target A Z " << aNucleus.GetA_asInt()
163  << " " << aNucleus.GetZ_asInt() << G4endl;
164  #endif
165 
166  theParticipants.Clean();
167 
168  theParticipants.SetProjectileNucleus( 0 );
169 
170  G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
171  ProjectileResidualMassNumber = 0;
172  ProjectileResidualCharge = 0;
173  ProjectileResidualExcitationEnergy = 0.0;
174  ProjectileResidual4Momentum = tmp;
175 
176  TargetResidualMassNumber = aNucleus.GetA_asInt();
177  TargetResidualCharge = aNucleus.GetZ_asInt();
178  TargetResidualExcitationEnergy = 0.0;
179  TargetResidual4Momentum = tmp;
180  G4double TargetResidualMass = G4ParticleTable::GetParticleTable()->GetIonTable()
181  ->GetIonMass( TargetResidualCharge, TargetResidualMassNumber );
182 
183  TargetResidual4Momentum.setE( TargetResidualMass );
184 
185  if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
186  // Projectile is a hadron : meson or baryon
187  PlabPerParticle = theProjectile.GetMomentum().z();
188  ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
189  ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
190  ProjectileResidualExcitationEnergy = 0.0;
191  //G4double ProjectileResidualMass = theProjectile.GetMass();
192  ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
193  ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
194  if ( PlabPerParticle < LowEnergyLimit ) {
195  HighEnergyInter = false;
196  } else {
197  HighEnergyInter = true;
198  }
199  } else {
200  if ( theProjectile.GetDefinition()->GetBaryonNumber() > 1 ) {
201  // Projectile is a nucleus
202  theParticipants.InitProjectileNucleus(theProjectile.GetDefinition()->GetBaryonNumber(),
203  G4int(theProjectile.GetDefinition()->GetPDGCharge()));
204  ProjectileResidualMassNumber = theProjectile.GetDefinition()->GetBaryonNumber();
205  ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
206  PlabPerParticle = theProjectile.GetMomentum().z() /
207  theProjectile.GetDefinition()->GetBaryonNumber();
208  if ( PlabPerParticle < LowEnergyLimit ) {
209  HighEnergyInter = false;
210  } else {
211  HighEnergyInter = true;
212  }
213  } else if ( theProjectile.GetDefinition()->GetBaryonNumber() < -1 ) {
214  // Projectile is an anti-nucleus
215  theParticipants.InitProjectileNucleus(
216  std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ),
217  std::abs( G4int( theProjectile.GetDefinition()->GetPDGCharge() ) ) );
218  theParticipants.theProjectileNucleus->StartLoop();
219  G4Nucleon* aNucleon;
220  while ( ( aNucleon = theParticipants.theProjectileNucleus->GetNextNucleon() ) ) { /* Loop checking, 10.08.2015, A.Ribon */
221  if ( aNucleon->GetDefinition() == G4Proton::Proton() ) {
223  } else if ( aNucleon->GetDefinition() == G4Neutron::Neutron() ) {
225  }
226  }
227  ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
228  ProjectileResidualCharge = std::abs( G4int(theProjectile.GetDefinition()->GetPDGCharge()) );
229  PlabPerParticle = theProjectile.GetMomentum().z() /
230  std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
231  if ( PlabPerParticle < LowEnergyLimit ) {
232  HighEnergyInter = false;
233  } else {
234  HighEnergyInter = true;
235  }
236  }
237 
238  G4ThreeVector BoostVector = theProjectile.GetMomentum() / theProjectile.GetTotalEnergy();
239  theParticipants.theProjectileNucleus->DoLorentzBoost( BoostVector );
240  theParticipants.theProjectileNucleus->DoLorentzContraction( BoostVector );
241  ProjectileResidualExcitationEnergy = 0.0;
242  //G4double ProjectileResidualMass = theProjectile.GetMass();
243  ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
244  ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
245  }
246 
247  // Init target nucleus
248  theParticipants.Init( aNucleus.GetA_asInt(), aNucleus.GetZ_asInt() );
249  //theParticipants.Init( aNucleus.GetA_asInt(), 0 ); // For h+neutron
250 
251  if ( theParameters != 0 ) delete theParameters;
252  theParameters = new G4FTFParameters( theProjectile.GetDefinition(), aNucleus.GetA_asInt(),
253  aNucleus.GetZ_asInt(), PlabPerParticle );
254 
255  if ( theAdditionalString.size() != 0 ) {
256  std::for_each( theAdditionalString.begin(), theAdditionalString.end(),
258  }
259  theAdditionalString.clear();
260 
261  #ifdef debugFTFmodel
262  G4cout << "FTF end of Init" << G4endl << G4endl;
263  #endif
264 
265  //if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 &&
266  // aNucleus.GetA_asInt() < 2 ) theParameters->SetProbabilityOfElasticScatt( 0.0 );
267 }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
void SetParticleType(G4Proton *aProton)
Definition: G4Nucleon.hh:77
virtual G4bool StartLoop()=0
virtual void DoLorentzContraction(const G4LorentzVector &theBoost)=0
virtual void SetProjectileNucleus(G4V3DNucleus *aNucleus)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
const G4ParticleDefinition * GetDefinition() const
G4IonTable * GetIonTable() const
G4GLOB_DLL std::ostream G4cout
virtual void Init(G4int theZ, G4int theA)
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1335
G4V3DNucleus * theProjectileNucleus
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetTotalEnergy() const
virtual void InitProjectileNucleus(G4int theZ, G4int theA)
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
static G4ParticleTable * GetParticleTable()
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
virtual G4Nucleon * GetNextNucleon()=0
void setVect(const Hep3Vector &)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double GetMass() const
static G4AntiNeutron * AntiNeutron()

Here is the call graph for this function:

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

Reimplemented from G4VPartonStringModel.

Definition at line 3366 of file G4FTFModel.cc.

3366  {
3367  desc << " FTF (Fritiof) Model \n"
3368  << "The FTF model is based on the well-known FRITIOF \n"
3369  << "model (B. Andersson et al., Nucl. Phys. B281, 289 \n"
3370  << "(1987)). Its first program implementation was given\n"
3371  << "by B. Nilsson-Almquist and E. Stenlund (Comp. Phys.\n"
3372  << "Comm. 43, 387 (1987)). The Fritiof model assumes \n"
3373  << "that all hadron-hadron interactions are binary \n"
3374  << "reactions, h_1+h_2->h_1'+h_2' where h_1' and h_2' \n"
3375  << "are excited states of the hadrons with continuous \n"
3376  << "mass spectra. The excited hadrons are considered as\n"
3377  << "QCD-strings, and the corresponding LUND-string \n"
3378  << "fragmentation model is applied for a simulation of \n"
3379  << "their decays. \n"
3380  << " The Fritiof model assumes that in the course of \n"
3381  << "a hadron-nucleus interaction a string originated \n"
3382  << "from the projectile can interact with various intra\n"
3383  << "nuclear nucleons and becomes into highly excited \n"
3384  << "states. The probability of multiple interactions is\n"
3385  << "calculated in the Glauber approximation. A cascading\n"
3386  << "of secondary particles was neglected as a rule. Due\n"
3387  << "to these, the original Fritiof model fails to des- \n"
3388  << "cribe a nuclear destruction and slow particle spectra.\n"
3389  << " In order to overcome the difficulties we enlarge\n"
3390  << "the model by the reggeon theory inspired model of \n"
3391  << "nuclear desctruction (Kh. Abdel-Waged and V.V. Uzhi-\n"
3392  << "nsky, Phys. Atom. Nucl. 60, 828 (1997); Yad. Fiz. 60, 925\n"
3393  << "(1997)). Momenta of the nucleons ejected from a nuc-\n"
3394  << "leus in the reggeon cascading are sampled according\n"
3395  << "to a Fermi motion algorithm presented in (EMU-01 \n"
3396  << "Collaboration (M.I. Adamovich et al.) Zeit. fur Phys.\n"
3397  << "A358, 337 (1997)). \n"
3398  << " New features were also added to the Fritiof model\n"
3399  << "implemented in Geant4: a simulation of elastic had-\n"
3400  << "ron-nucleon scatterings, a simulation of binary \n"
3401  << "reactions like NN>NN* in hadron-nucleon interactions,\n"
3402  << "a separate simulation of single diffractive and non-\n"
3403  << " diffractive events. These allowed to describe after\n"
3404  << "model parameter tuning a wide set of experimental \n"
3405  << "data. \n";
3406 }

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