Geant4  10.00.p01
G4FTFModel.hh
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 // $Id: G4FTFModel.hh 74627 2013-10-17 07:04:38Z gcosmo $
28 // GEANT4 tag $Name: $
29 //
30 // Class Description
31 // Final state production code for hadron inelastic scattering above 3 GeV
32 // based on the modeling ansatz used in FRITIOF.
33 // To be used in your physics list in case you need this physics.
34 // In this case you want to register an object of this class with an object
35 // of G4TheoFSGenerator.
36 // Class Description - End
37 
38 #ifndef G4FTFModel_h
39 #define G4FTFModel_h 1
40 
41 // ------------------------------------------------------------
42 // GEANT 4 class header file
43 //
44 // ---------------- G4FTFModel ----------------
45 // by Gunter Folger, May 1998.
46 // class implementing the excitation in the FTF Parton String Model
47 // ------------------------------------------------------------
48 
49 #include "G4VPartonStringModel.hh"
50 #include "G4FTFParameters.hh"
51 #include "G4FTFParticipants.hh"
52 #include "G4ExcitedStringVector.hh"
54 #include "G4ElasticHNScattering.hh"
55 #include "G4FTFAnnihilation.hh"
56 #include "G4Proton.hh"
57 #include "G4Neutron.hh"
58 
59 class G4VSplitableHadron;
60 class G4ExcitedString;
61 
62 
64 
65  public:
66  G4FTFModel( const G4String& modelName = "FTF" );
67  ~G4FTFModel();
68 
69  void Init( const G4Nucleus& aNucleus, const G4DynamicParticle& aProjectile );
74 
75  virtual void ModelDescription( std::ostream& ) const;
76 
77  private:
78  G4FTFModel( const G4FTFModel& right );
79  const G4FTFModel& operator=( const G4FTFModel& right );
80  int operator==( const G4FTFModel& right ) const;
81  int operator!=( const G4FTFModel& right ) const;
82 
83  void StoreInvolvedNucleon();
84  void ReggeonCascade();
88  void GetResiduals();
89  G4bool AdjustNucleons( G4VSplitableHadron* SelectedAntiBaryon,
90  G4Nucleon* ProjectileNucleon,
91  G4VSplitableHadron* SelectedTargetNucleon,
92  G4Nucleon* TargetNucleon,
93  G4bool Annihilation );
94  G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
95 
98 
101 
104 
109 
110  std::vector< G4VSplitableHadron* > theAdditionalString;
111 
114 
119 
124 
125 };
126 
127 
130 }
131 
132 
135 }
136 
137 
140 }
141 
142 #endif
G4ElasticHNScattering * theElastic
Definition: G4FTFModel.hh:107
std::vector< G4ExcitedString * > G4ExcitedStringVector
G4V3DNucleus * GetWoundedNucleus() const
Definition: G4FTFModel.hh:128
G4bool PutOnMassShell()
Definition: G4FTFModel.cc:513
G4FTFModel(const G4String &modelName="FTF")
Definition: G4FTFModel.cc:69
CLHEP::Hep3Vector G4ThreeVector
G4bool AdjustNucleons(G4VSplitableHadron *SelectedAntiBaryon, G4Nucleon *ProjectileNucleon, G4VSplitableHadron *SelectedTargetNucleon, G4Nucleon *TargetNucleon, G4bool Annihilation)
Definition: G4FTFModel.cc:1641
void ReggeonCascade()
Definition: G4FTFModel.cc:409
G4ReactionProduct theProjectile
Definition: G4FTFModel.hh:96
int operator!=(const G4FTFModel &right) const
G4V3DNucleus * GetTargetNucleus() const
Definition: G4FTFModel.hh:133
G4int TargetResidualMassNumber
Definition: G4FTFModel.hh:121
int G4int
Definition: G4Types.hh:78
G4FTFParameters * theParameters
Definition: G4FTFModel.hh:105
G4FTFAnnihilation * theAnnihilation
Definition: G4FTFModel.hh:108
G4bool ExciteParticipants()
Definition: G4FTFModel.cc:1468
G4Nucleon * TheInvolvedNucleonsOfTarget[250]
Definition: G4FTFModel.hh:99
const G4FTFModel & operator=(const G4FTFModel &right)
std::vector< G4VSplitableHadron * > theAdditionalString
Definition: G4FTFModel.hh:110
void GetResiduals()
Definition: G4FTFModel.cc:3034
G4FTFParticipants theParticipants
Definition: G4FTFModel.hh:97
G4V3DNucleus * GetProjectileNucleus() const
Definition: G4FTFModel.hh:138
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
Definition: G4FTFModel.cc:3205
G4int ProjectileResidualMassNumber
Definition: G4FTFModel.hh:116
G4DiffractiveExcitation * theExcitation
Definition: G4FTFModel.hh:106
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual G4V3DNucleus * GetWoundedNucleus() const
G4double LowEnergyLimit
Definition: G4FTFModel.hh:112
G4ExcitedStringVector * BuildStrings()
Definition: G4FTFModel.cc:2762
G4Nucleon * TheInvolvedNucleonsOfProjectile[250]
Definition: G4FTFModel.hh:102
void Init(const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile)
Definition: G4FTFModel.cc:141
virtual void ModelDescription(std::ostream &) const
Definition: G4FTFModel.cc:3224
G4int ProjectileResidualCharge
Definition: G4FTFModel.hh:117
G4int NumberOfInvolvedNucleonsOfProjectile
Definition: G4FTFModel.hh:103
void StoreInvolvedNucleon()
Definition: G4FTFModel.cc:358
G4int TargetResidualCharge
Definition: G4FTFModel.hh:122
G4LorentzVector TargetResidual4Momentum
Definition: G4FTFModel.hh:120
G4bool HighEnergyInter
Definition: G4FTFModel.hh:113
double G4double
Definition: G4Types.hh:76
G4LorentzVector ProjectileResidual4Momentum
Definition: G4FTFModel.hh:115
G4double ProjectileResidualExcitationEnergy
Definition: G4FTFModel.hh:118
G4ExcitedStringVector * GetStrings()
Definition: G4FTFModel.cc:257
G4int NumberOfInvolvedNucleonsOfTarget
Definition: G4FTFModel.hh:100
G4double TargetResidualExcitationEnergy
Definition: G4FTFModel.hh:123
CLHEP::HepLorentzVector G4LorentzVector
int operator==(const G4FTFModel &right) const