Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FTFModel.hh
이 파일의 문서화 페이지로 가기
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 110870 2018-06-22 12:14:16Z 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  public:
65  G4FTFModel( const G4String& modelName = "FTF" );
66  ~G4FTFModel();
67 
68  void Init( const G4Nucleus& aNucleus, const G4DynamicParticle& aProjectile );
73 
74  virtual void ModelDescription( std::ostream& ) const;
75 
76  private:
77  G4FTFModel( const G4FTFModel& right );
78  const G4FTFModel& operator=( const G4FTFModel& right );
79  int operator==( const G4FTFModel& right ) const;
80  int operator!=( const G4FTFModel& right ) const;
81 
82  void StoreInvolvedNucleon();
83  void ReggeonCascade();
87  void GetResiduals();
88 
89  G4bool AdjustNucleons( G4VSplitableHadron* SelectedAntiBaryon,
90  G4Nucleon* ProjectileNucleon,
91  G4VSplitableHadron* SelectedTargetNucleon,
92  G4Nucleon* TargetNucleon,
93  G4bool Annihilation );
94  // The "AdjustNucleons" method uses the following struct and 3 new utility methods:
95  struct CommonVariables {
98  G4double SqrtS = 0.0, S = 0.0, SumMasses = 0.0,
101  Mprojectile = 0.0, M2projectile = 0.0, Pzprojectile = 0.0, Eprojectile = 0.0,
103  Mtarget = 0.0, M2target = 0.0, Pztarget = 0.0, Etarget = 0.0, WminusTarget = 0.0,
111  };
113  G4VSplitableHadron* SelectedAntiBaryon,
114  G4Nucleon* ProjectileNucleon,
115  G4VSplitableHadron* SelectedTargetNucleon,
116  G4Nucleon* TargetNucleon,
117  G4bool Annihilation,
121  void AdjustNucleonsAlgorithm_afterSampling( G4int interactionCase,
122  G4VSplitableHadron* SelectedAntiBaryon,
123  G4VSplitableHadron* SelectedTargetNucleon,
125 
126  G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
127 
128  G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
129  G4LorentzVector& residualMomentum, G4double& sumMasses,
130  G4double& residualExcitationEnergy, G4double& residualMass,
131  G4int& residualMassNumber, G4int& residualCharge );
132  // Utility method used by PutOnMassShell.
133 
134  G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
135  G4Nucleon* involvedNucleons[], G4double& sumMasses );
136  // Utility method used by PutOnMassShell.
137 
138  G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
139  G4double dCor, G4V3DNucleus* nucleus,
140  const G4LorentzVector& pResidual,
141  const G4double residualMass, const G4int residualMassNumber,
142  const G4int numberOfInvolvedNucleons,
143  G4Nucleon* involvedNucleons[], G4double& mass2 );
144 
145  // Utility method used by PutOnMassShell.
146 
147  G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
148  const G4double projectileMass2, const G4double targetMass2,
149  const G4double nucleusY, const G4bool isProjectileNucleus,
150  const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
151  G4double& targetWminus, G4double& projectileWplus, G4bool& success );
152  // Utility method used by PutOnMassShell.
153 
154  G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
155  const G4LorentzRotation& boostFromCmsToLab,
156  const G4double residualMass, const G4int residualMassNumber,
157  const G4int numberOfInvolvedNucleons,
158  G4Nucleon* involvedNucleons[],
159  G4LorentzVector& residual4Momentum );
160  // Utility method used by PutOnMassShell.
161 
164 
167 
170 
175 
176  std::vector< G4VSplitableHadron* > theAdditionalString;
177 
180 
185 
190 };
191 
192 
195 }
196 
197 
200 }
201 
202 
205 }
206 
207 #endif
208 
G4LorentzVector TResidual4Momentum
Definition: G4FTFModel.hh:109
G4DiffractiveExcitation * theExcitation
Definition: G4FTFModel.hh:172
const G4FTFModel & operator=(const G4FTFModel &right)
G4double TargetResidualExcitationEnergy
Definition: G4FTFModel.hh:189
G4bool FinalizeKinematics(const G4double w, const G4bool isProjectileNucleus, const G4LorentzRotation &boostFromCmsToLab, const G4double residualMass, const G4int residualMassNumber, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4LorentzVector &residual4Momentum)
Definition: G4FTFModel.cc:2985
void GetResiduals()
Definition: G4FTFModel.cc:2255
G4FTFAnnihilation * theAnnihilation
Definition: G4FTFModel.hh:174
G4bool AdjustNucleonsAlgorithm_Sampling(G4int interactionCase, CommonVariables &common)
Definition: G4FTFModel.cc:1476
G4V3DNucleus * GetWoundedNucleus() const
Definition: G4FTFModel.hh:193
G4int NumberOfInvolvedNucleonsOfTarget
Definition: G4FTFModel.hh:166
G4LorentzVector TargetResidual4Momentum
Definition: G4FTFModel.hh:186
std::vector< G4ExcitedString * > G4ExcitedStringVector
void ReggeonCascade()
Definition: G4FTFModel.cc:441
G4LorentzRotation toCms
Definition: G4FTFModel.hh:110
G4bool ExciteParticipants()
Definition: G4FTFModel.cc:851
std::vector< G4VSplitableHadron * > theAdditionalString
Definition: G4FTFModel.hh:176
static int FASTCALL common(PROLOG_STATE *state, int tok)
Definition: xmlrole.cc:1309
virtual void ModelDescription(std::ostream &) const
Definition: G4FTFModel.cc:3057
G4FTFModel(const G4String &modelName="FTF")
Definition: G4FTFModel.cc:72
G4Nucleon * TheInvolvedNucleonsOfProjectile[250]
Definition: G4FTFModel.hh:168
G4V3DNucleus * GetProjectileNucleus() const
Definition: G4FTFModel.hh:203
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
Definition: G4FTFModel.cc:2577
void Init(const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile)
Definition: G4FTFModel.cc:156
G4Nucleon * TheInvolvedNucleonsOfTarget[250]
Definition: G4FTFModel.hh:165
G4int ProjectileResidualMassNumber
Definition: G4FTFModel.hh:182
G4LorentzVector ProjectileResidual4Momentum
Definition: G4FTFModel.hh:181
G4int NumberOfInvolvedNucleonsOfProjectile
Definition: G4FTFModel.hh:169
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ReactionProduct theProjectile
Definition: G4FTFModel.hh:162
void AdjustNucleonsAlgorithm_afterSampling(G4int interactionCase, G4VSplitableHadron *SelectedAntiBaryon, G4VSplitableHadron *SelectedTargetNucleon, CommonVariables &common)
Definition: G4FTFModel.cc:1811
G4ElasticHNScattering * theElastic
Definition: G4FTFModel.hh:173
G4LorentzVector Pprojectile
Definition: G4FTFModel.hh:109
G4int AdjustNucleonsAlgorithm_beforeSampling(G4int interactionCase, G4VSplitableHadron *SelectedAntiBaryon, G4Nucleon *ProjectileNucleon, G4VSplitableHadron *SelectedTargetNucleon, G4Nucleon *TargetNucleon, G4bool Annihilation, CommonVariables &common)
Definition: G4FTFModel.cc:1155
G4int ProjectileResidualCharge
Definition: G4FTFModel.hh:183
G4double LowEnergyLimit
Definition: G4FTFModel.hh:178
G4ExcitedStringVector * GetStrings()
Definition: G4FTFModel.cc:286
virtual G4V3DNucleus * GetWoundedNucleus() const
G4bool SamplingNucleonKinematics(G4double averagePt2, const G4double maxPt2, G4double dCor, G4V3DNucleus *nucleus, const G4LorentzVector &pResidual, const G4double residualMass, const G4int residualMassNumber, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4double &mass2)
Definition: G4FTFModel.cc:2747
G4int TargetResidualMassNumber
Definition: G4FTFModel.hh:187
int G4int
Definition: G4Types.hh:78
G4int TargetResidualCharge
Definition: G4FTFModel.hh:188
virtual G4V3DNucleus * GetProjectileNucleus() const
G4bool HighEnergyInter
Definition: G4FTFModel.hh:179
G4FTFParticipants theParticipants
Definition: G4FTFModel.hh:163
G4LorentzRotation toLab
Definition: G4FTFModel.hh:110
void StoreInvolvedNucleon()
Definition: G4FTFModel.cc:391
G4bool PutOnMassShell()
Definition: G4FTFModel.cc:544
G4V3DNucleus * GetTargetNucleus() const
Definition: G4FTFModel.hh:198
G4double ProjectileResidualExcitationEnergy
Definition: G4FTFModel.hh:184
int operator==(const G4FTFModel &right) const
G4bool AdjustNucleons(G4VSplitableHadron *SelectedAntiBaryon, G4Nucleon *ProjectileNucleon, G4VSplitableHadron *SelectedTargetNucleon, G4Nucleon *TargetNucleon, G4bool Annihilation)
Definition: G4FTFModel.cc:1040
G4bool CheckKinematics(const G4double sValue, const G4double sqrtS, const G4double projectileMass2, const G4double targetMass2, const G4double nucleusY, const G4bool isProjectileNucleus, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4double &targetWminus, G4double &projectileWplus, G4bool &success)
Definition: G4FTFModel.cc:2911
G4bool GenerateDeltaIsobar(const G4double sqrtS, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4double &sumMasses)
Definition: G4FTFModel.cc:2683
G4ExcitedStringVector * BuildStrings()
Definition: G4FTFModel.cc:1950
int operator!=(const G4FTFModel &right) const
G4LorentzVector PResidual4Momentum
Definition: G4FTFModel.hh:109
G4bool ComputeNucleusProperties(G4V3DNucleus *nucleus, G4LorentzVector &nucleusMomentum, G4LorentzVector &residualMomentum, G4double &sumMasses, G4double &residualExcitationEnergy, G4double &residualMass, G4int &residualMassNumber, G4int &residualCharge)
Definition: G4FTFModel.cc:2599
G4FTFParameters * theParameters
Definition: G4FTFModel.hh:171