Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28// Class Description
29// Final state production code for hadron inelastic scattering above 3 GeV
30// based on the modeling ansatz used in FRITIOF.
31// To be used in your physics list in case you need this physics.
32// In this case you want to register an object of this class with an object
33// of G4TheoFSGenerator.
34// Class Description - End
35
36#ifndef G4FTFModel_h
37#define G4FTFModel_h 1
38
39// ------------------------------------------------------------
40// GEANT 4 class header file
41//
42// ---------------- G4FTFModel ----------------
43// by Gunter Folger, May 1998.
44// class implementing the excitation in the FTF Parton String Model
45// ------------------------------------------------------------
46
48#include "G4FTFParameters.hh"
49#include "G4FTFParticipants.hh"
53#include "G4FTFAnnihilation.hh"
54#include "G4Proton.hh"
55#include "G4Neutron.hh"
56
58class G4ExcitedString;
59
60
62 public:
63 G4FTFModel( const G4String& modelName = "FTF" );
64 ~G4FTFModel() override;
65
67 G4V3DNucleus* GetWoundedNucleus() const override;
68 G4V3DNucleus* GetProjectileNucleus() const override;
69
70 void ModelDescription( std::ostream& ) const override;
71
72 G4FTFModel( const G4FTFModel& right ) = delete;
73 const G4FTFModel& operator=( const G4FTFModel& right ) = delete;
74 G4bool operator==( const G4FTFModel& right ) const = delete;
75 G4bool operator!=( const G4FTFModel& right ) const = delete;
76
77 protected:
78 void Init( const G4Nucleus& aNucleus,
79 const G4DynamicParticle& aProjectile ) override;
81
82 private:
83 void StoreInvolvedNucleon();
84 void ReggeonCascade();
85 G4bool PutOnMassShell();
86 G4bool ExciteParticipants();
87 void BuildStrings( G4ExcitedStringVector* strings );
88 void GetResiduals();
89
90 G4bool AdjustNucleons( G4VSplitableHadron* SelectedAntiBaryon,
91 G4Nucleon* ProjectileNucleon,
92 G4VSplitableHadron* SelectedTargetNucleon,
93 G4Nucleon* TargetNucleon,
94 G4bool Annihilation );
95 // The "AdjustNucleons" method uses the following struct and 3 new utility methods:
96 struct CommonVariables {
97 G4int TResidualMassNumber = 0, TResidualCharge = 0, PResidualMassNumber = 0,
98 PResidualCharge = 0;
99 G4double SqrtS = 0.0, S = 0.0, SumMasses = 0.0,
100 TResidualExcitationEnergy = 0.0, TResidualMass = 0.0, TNucleonMass = 0.0,
101 PResidualExcitationEnergy = 0.0, PResidualMass = 0.0, PNucleonMass = 0.0,
102 Mprojectile = 0.0, M2projectile = 0.0, Pzprojectile = 0.0, Eprojectile = 0.0,
103 WplusProjectile = 0.0,
104 Mtarget = 0.0, M2target = 0.0, Pztarget = 0.0, Etarget = 0.0, WminusTarget = 0.0,
105 Mt2targetNucleon = 0.0, PztargetNucleon = 0.0, EtargetNucleon = 0.0,
106 Mt2projectileNucleon = 0.0, PzprojectileNucleon = 0.0, EprojectileNucleon = 0.0,
107 YtargetNucleus = 0.0, YprojectileNucleus = 0.0,
108 XminusNucleon = 0.0, XplusNucleon = 0.0, XminusResidual = 0.0, XplusResidual = 0.0;
109 G4ThreeVector PtNucleon, PtResidual, PtNucleonP, PtResidualP, PtNucleonT, PtResidualT;
110 G4LorentzVector Psum, Pprojectile, Ptmp, Ptarget, TResidual4Momentum, PResidual4Momentum;
111 G4LorentzRotation toCms, toLab;
112 };
113 G4int AdjustNucleonsAlgorithm_beforeSampling( G4int interactionCase,
114 G4VSplitableHadron* SelectedAntiBaryon,
115 G4Nucleon* ProjectileNucleon,
116 G4VSplitableHadron* SelectedTargetNucleon,
117 G4Nucleon* TargetNucleon,
118 G4bool Annihilation,
119 CommonVariables& common );
120 G4bool AdjustNucleonsAlgorithm_Sampling( G4int interactionCase,
121 CommonVariables& common );
122 void AdjustNucleonsAlgorithm_afterSampling( G4int interactionCase,
123 G4VSplitableHadron* SelectedAntiBaryon,
124 G4VSplitableHadron* SelectedTargetNucleon,
125 CommonVariables& common );
126
127 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
128
129 G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
130 G4LorentzVector& residualMomentum, G4double& sumMasses,
131 G4double& residualExcitationEnergy, G4double& residualMass,
132 G4int& residualMassNumber, G4int& residualCharge );
133 // Utility method used by PutOnMassShell.
134
135 G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
136 G4Nucleon* involvedNucleons[], G4double& sumMasses );
137 // Utility method used by PutOnMassShell.
138
139 G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
140 G4double dCor, G4V3DNucleus* nucleus,
141 const G4LorentzVector& pResidual,
142 const G4double residualMass, const G4int residualMassNumber,
143 const G4int numberOfInvolvedNucleons,
144 G4Nucleon* involvedNucleons[], G4double& mass2 );
145
146 // Utility method used by PutOnMassShell.
147
148 G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
149 const G4double projectileMass2, const G4double targetMass2,
150 const G4double nucleusY, const G4bool isProjectileNucleus,
151 const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
152 G4double& targetWminus, G4double& projectileWplus, G4bool& success );
153 // Utility method used by PutOnMassShell.
154
155 G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
156 const G4LorentzRotation& boostFromCmsToLab,
157 const G4double residualMass, const G4int residualMassNumber,
158 const G4int numberOfInvolvedNucleons,
159 G4Nucleon* involvedNucleons[],
160 G4LorentzVector& residual4Momentum );
161 // Utility method used by PutOnMassShell.
162
163 G4ReactionProduct theProjectile;
164 G4FTFParticipants theParticipants;
165
166 G4Nucleon* TheInvolvedNucleonsOfTarget[250];
167 G4int NumberOfInvolvedNucleonsOfTarget;
168
169 G4Nucleon* TheInvolvedNucleonsOfProjectile[250];
170 G4int NumberOfInvolvedNucleonsOfProjectile;
171
172 G4FTFParameters* theParameters;
173 G4DiffractiveExcitation* theExcitation;
174 G4ElasticHNScattering* theElastic;
175 G4FTFAnnihilation* theAnnihilation;
176
177 std::vector< G4VSplitableHadron* > theAdditionalString;
178
179 G4double LowEnergyLimit;
180 G4bool HighEnergyInter;
181
182 G4LorentzVector ProjectileResidual4Momentum;
183 G4int ProjectileResidualMassNumber;
184 G4int ProjectileResidualCharge;
185 G4double ProjectileResidualExcitationEnergy;
186
187 G4LorentzVector TargetResidual4Momentum;
188 G4int TargetResidualMassNumber;
189 G4int TargetResidualCharge;
190 G4double TargetResidualExcitationEnergy;
191};
192
194 return theParticipants.GetWoundedNucleus();
195}
196
198 return theParticipants.GetWoundedNucleus();
199}
200
202 return theParticipants.GetProjectileNucleus();
203}
204
205#endif
206
double S(double temp)
std::vector< G4ExcitedString * > G4ExcitedStringVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4FTFModel & operator=(const G4FTFModel &right)=delete
G4V3DNucleus * GetTargetNucleus() const
Definition: G4FTFModel.hh:197
G4ExcitedStringVector * GetStrings() override
Definition: G4FTFModel.cc:290
G4bool operator!=(const G4FTFModel &right) const =delete
G4V3DNucleus * GetWoundedNucleus() const override
Definition: G4FTFModel.hh:193
G4FTFModel(const G4FTFModel &right)=delete
~G4FTFModel() override
Definition: G4FTFModel.cc:113
G4V3DNucleus * GetProjectileNucleus() const override
Definition: G4FTFModel.hh:201
void Init(const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile) override
Definition: G4FTFModel.cc:153
void ModelDescription(std::ostream &) const override
Definition: G4FTFModel.cc:3069
G4bool operator==(const G4FTFModel &right) const =delete
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual G4V3DNucleus * GetWoundedNucleus() const