Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4FTFModel Class Reference

#include <G4FTFModel.hh>

Inheritance diagram for G4FTFModel:
G4VPartonStringModel G4VHighEnergyGenerator

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 69 of file G4FTFModel.cc.

References python.hepunit::MeV, python.hepunit::perCent, G4VHighEnergyGenerator::SetEnergyMomentumCheckLevels(), and G4VPartonStringModel::SetThisPointer().

69  :
70  G4VPartonStringModel( modelName ),
71  theExcitation( new G4DiffractiveExcitation() ),
72  theElastic( new G4ElasticHNScattering() ),
73  theAnnihilation( new G4FTFAnnihilation() )
74 {
76  theParameters = 0;
77  NumberOfInvolvedNucleonsOfTarget = 0;
78  NumberOfInvolvedNucleonsOfProjectile= 0;
79 
80  LowEnergyLimit = 5000.0*MeV;
81  HighEnergyInter = true;
82 
83  G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
84  ProjectileResidual4Momentum = tmp;
85  ProjectileResidualMassNumber = 0;
86  ProjectileResidualCharge = 0;
87  ProjectileResidualExcitationEnergy = 0.0;
88 
89  TargetResidual4Momentum = tmp;
90  TargetResidualMassNumber = 0;
91  TargetResidualCharge = 0;
92  TargetResidualExcitationEnergy = 0.0;
93 
95 }
G4VPartonStringModel(const G4String &modelName="Parton String Model")
void SetThisPointer(G4VPartonStringModel *aPointer)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double AbsoluteLevel)
float perCent
Definition: hepunit.py:239
G4FTFModel::~G4FTFModel ( )

Definition at line 105 of file G4FTFModel.cc.

References G4Nucleon::GetSplitableHadron().

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

Member Function Documentation

G4V3DNucleus * G4FTFModel::GetProjectileNucleus ( ) const
inlinevirtual

Reimplemented from G4VPartonStringModel.

Definition at line 138 of file G4FTFModel.hh.

References G4VParticipants::GetProjectileNucleus().

Referenced by GetStrings().

138  {
139  return theParticipants.GetProjectileNucleus();
140 }
virtual G4V3DNucleus * GetProjectileNucleus() const
G4ExcitedStringVector * G4FTFModel::GetStrings ( )
virtual

Implements G4VPartonStringModel.

Definition at line 257 of file G4FTFModel.cc.

References G4cin, G4cout, G4endl, G4FTFParticipants::GetInteraction(), G4FTFParticipants::GetList(), G4InteractionContent::GetProjectile(), GetProjectileNucleus(), G4Nucleon::GetSplitableHadron(), G4FTFParticipants::Next(), and G4FTFParticipants::StartLoop().

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

Definition at line 133 of file G4FTFModel.hh.

References G4VParticipants::GetWoundedNucleus().

133  {
134  return theParticipants.GetWoundedNucleus();
135 }
virtual G4V3DNucleus * GetWoundedNucleus() const
G4V3DNucleus * G4FTFModel::GetWoundedNucleus ( ) const
inlinevirtual

Implements G4VPartonStringModel.

Definition at line 128 of file G4FTFModel.hh.

References G4VParticipants::GetWoundedNucleus().

128  {
129  return theParticipants.GetWoundedNucleus();
130 }
virtual G4V3DNucleus * GetWoundedNucleus() const
void G4FTFModel::Init ( const G4Nucleus aNucleus,
const G4DynamicParticle aProjectile 
)
virtual

Implements G4VPartonStringModel.

Definition at line 141 of file G4FTFModel.cc.

References G4AntiNeutron::AntiNeutron(), G4AntiProton::AntiProton(), G4V3DNucleus::DoLorentzBoost(), G4V3DNucleus::DoLorentzContraction(), G4cout, G4endl, G4Nucleus::GetA_asInt(), G4ParticleDefinition::GetBaryonNumber(), G4Nucleon::GetDefinition(), G4ReactionProduct::GetDefinition(), G4IonTable::GetIonMass(), G4ParticleTable::GetIonTable(), G4ReactionProduct::GetMass(), G4ReactionProduct::GetMomentum(), G4V3DNucleus::GetNextNucleon(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ReactionProduct::GetTotalEnergy(), G4Nucleus::GetZ_asInt(), G4VParticipants::Init(), G4VParticipants::InitProjectileNucleus(), G4Neutron::Neutron(), G4Proton::Proton(), CLHEP::HepLorentzVector::setE(), G4Nucleon::SetParticleType(), G4VParticipants::SetProjectileNucleus(), CLHEP::HepLorentzVector::setVect(), G4V3DNucleus::StartLoop(), G4VParticipants::theProjectileNucleus, and CLHEP::Hep3Vector::z().

141  {
142 
143  theProjectile = aProjectile;
144 
145  G4double PlabPerParticle( 0.0 ); // Laboratory momentum Pz per particle/nucleon
146 
147  #ifdef debugFTFmodel
148  G4cout << "FTF init Proj Name " << theProjectile.GetDefinition()->GetParticleName() << G4endl
149  << "FTF init Proj Mass " << theProjectile.GetMass()
150  << " " << theProjectile.GetMomentum() << G4endl
151  << "FTF init Proj B Q " << theProjectile.GetDefinition()->GetBaryonNumber()
152  << " " << (G4int) theProjectile.GetDefinition()->GetPDGCharge() << G4endl
153  << "FTF init Target A Z " << aNucleus.GetA_asInt()
154  << " " << aNucleus.GetZ_asInt() << G4endl;
155  #endif
156 
157  theParticipants.SetProjectileNucleus( 0 );
158 
159  G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
160  ProjectileResidualMassNumber = 0;
161  ProjectileResidualCharge = 0;
162  ProjectileResidualExcitationEnergy = 0.0;
163  ProjectileResidual4Momentum = tmp;
164 
165  TargetResidualMassNumber = aNucleus.GetA_asInt();
166  TargetResidualCharge = aNucleus.GetZ_asInt();
167  TargetResidualExcitationEnergy = 0.0;
168  TargetResidual4Momentum = tmp;
169  G4double TargetResidualMass = G4ParticleTable::GetParticleTable()->GetIonTable()
170  ->GetIonMass( TargetResidualCharge, TargetResidualMassNumber );
171 
172  TargetResidual4Momentum.setE( TargetResidualMass );
173 
174  if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
175  // Projectile is a hadron : meson or baryon
176  PlabPerParticle = theProjectile.GetMomentum().z();
177  ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
178  ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
179  ProjectileResidualExcitationEnergy = 0.0;
180  // G4double ProjectileResidualMass = theProjectile.GetMass();
181  ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
182  ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
183  if ( PlabPerParticle < LowEnergyLimit ) {
184  HighEnergyInter = false;
185  } else {
186  HighEnergyInter = true;
187  }
188  } else {
189  if ( theProjectile.GetDefinition()->GetBaryonNumber() > 1 ) {
190  // Projectile is a nucleus
191  theParticipants.InitProjectileNucleus(theProjectile.GetDefinition()->GetBaryonNumber(),
192  G4int(theProjectile.GetDefinition()->GetPDGCharge()));
193  ProjectileResidualMassNumber = theProjectile.GetDefinition()->GetBaryonNumber();
194  ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
195  PlabPerParticle = theProjectile.GetMomentum().z() /
196  theProjectile.GetDefinition()->GetBaryonNumber();
197  if ( PlabPerParticle < LowEnergyLimit ) {
198  HighEnergyInter = false;
199  } else {
200  HighEnergyInter = true;
201  }
202  } else if ( theProjectile.GetDefinition()->GetBaryonNumber() < -1 ) {
203  // Projectile is an anti-nucleus
204  theParticipants.InitProjectileNucleus(
205  std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ),
206  std::abs( G4int( theProjectile.GetDefinition()->GetPDGCharge() ) ) );
207  theParticipants.theProjectileNucleus->StartLoop();
208  G4Nucleon* aNucleon;
209  while ( ( aNucleon = theParticipants.theProjectileNucleus->GetNextNucleon() ) ) {
210  if ( aNucleon->GetDefinition() == G4Proton::Proton() ) {
212  } else if ( aNucleon->GetDefinition() == G4Neutron::Neutron() ) {
214  }
215  }
216  ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
217  ProjectileResidualCharge = std::abs( G4int(theProjectile.GetDefinition()->GetPDGCharge()) );
218  PlabPerParticle = theProjectile.GetMomentum().z() /
219  std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
220  if ( PlabPerParticle < LowEnergyLimit ) {
221  HighEnergyInter = false;
222  } else {
223  HighEnergyInter = true;
224  }
225  }
226  G4ThreeVector BoostVector = theProjectile.GetMomentum() / theProjectile.GetTotalEnergy();
227  theParticipants.theProjectileNucleus->DoLorentzBoost( BoostVector );
228  theParticipants.theProjectileNucleus->DoLorentzContraction( BoostVector );
229  ProjectileResidualExcitationEnergy = 0.0;
230  //G4double ProjectileResidualMass = theProjectile.GetMass();
231  ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
232  ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
233  }
234 
235  // Init target nucleus
236  theParticipants.Init( aNucleus.GetA_asInt(), aNucleus.GetZ_asInt() );
237 
238  if ( theParameters != 0 ) delete theParameters;
239  theParameters = new G4FTFParameters( theProjectile.GetDefinition(), aNucleus.GetA_asInt(),
240  aNucleus.GetZ_asInt(), PlabPerParticle );
241 
242  if ( theAdditionalString.size() != 0 ) {
243  std::for_each( theAdditionalString.begin(), theAdditionalString.end(),
245  }
246  theAdditionalString.clear();
247 
248  #ifdef debugFTFmodel
249  G4cout << "FTF end of Init" << G4endl << G4endl;
250  #endif
251 
252 }
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
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:1232
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 G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
virtual G4Nucleon * GetNextNucleon()=0
void setVect(const Hep3Vector &)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double GetMass() const
static G4AntiNeutron * AntiNeutron()
void G4FTFModel::ModelDescription ( std::ostream &  desc) const
virtual

Reimplemented from G4VPartonStringModel.

Definition at line 3224 of file G4FTFModel.cc.

References G4endl.

3224  {
3225  desc << "please add description here" << G4endl;
3226 }
#define G4endl
Definition: G4ios.hh:61

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