Geant4-11
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 void SetImpactParameter( const G4double b_value );
79 void SetBminBmax( const G4double bmin_value, const G4double bmax_value );
81 G4double GetBmin() const;
82 G4double GetBmax() const;
86
87 protected:
88 void Init( const G4Nucleus& aNucleus,
89 const G4DynamicParticle& aProjectile ) override;
91
92 private:
94 void ReggeonCascade();
97 void BuildStrings( G4ExcitedStringVector* strings );
98 void GetResiduals();
99
100 G4bool AdjustNucleons( G4VSplitableHadron* SelectedAntiBaryon,
101 G4Nucleon* ProjectileNucleon,
102 G4VSplitableHadron* SelectedTargetNucleon,
103 G4Nucleon* TargetNucleon,
104 G4bool Annihilation );
105 // The "AdjustNucleons" method uses the following struct and 3 new utility methods:
109 G4double SqrtS = 0.0, S = 0.0, SumMasses = 0.0,
114 Mtarget = 0.0, M2target = 0.0, Pztarget = 0.0, Etarget = 0.0, WminusTarget = 0.0,
122 };
124 G4VSplitableHadron* SelectedAntiBaryon,
125 G4Nucleon* ProjectileNucleon,
126 G4VSplitableHadron* SelectedTargetNucleon,
127 G4Nucleon* TargetNucleon,
128 G4bool Annihilation,
132 void AdjustNucleonsAlgorithm_afterSampling( G4int interactionCase,
133 G4VSplitableHadron* SelectedAntiBaryon,
134 G4VSplitableHadron* SelectedTargetNucleon,
136
137 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
138
139 G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
140 G4LorentzVector& residualMomentum, G4double& sumMasses,
141 G4double& residualExcitationEnergy, G4double& residualMass,
142 G4int& residualMassNumber, G4int& residualCharge );
143 // Utility method used by PutOnMassShell.
144
145 G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
146 G4Nucleon* involvedNucleons[], G4double& sumMasses );
147 // Utility method used by PutOnMassShell.
148
149 G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
150 G4double dCor, G4V3DNucleus* nucleus,
151 const G4LorentzVector& pResidual,
152 const G4double residualMass, const G4int residualMassNumber,
153 const G4int numberOfInvolvedNucleons,
154 G4Nucleon* involvedNucleons[], G4double& mass2 );
155
156 // Utility method used by PutOnMassShell.
157
158 G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
159 const G4double projectileMass2, const G4double targetMass2,
160 const G4double nucleusY, const G4bool isProjectileNucleus,
161 const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
162 G4double& targetWminus, G4double& projectileWplus, G4bool& success );
163 // Utility method used by PutOnMassShell.
164
165 G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
166 const G4LorentzRotation& boostFromCmsToLab,
167 const G4double residualMass, const G4int residualMassNumber,
168 const G4int numberOfInvolvedNucleons,
169 G4Nucleon* involvedNucleons[],
170 G4LorentzVector& residual4Momentum );
171 // Utility method used by PutOnMassShell.
172
175
178
181
186
187 std::vector< G4VSplitableHadron* > theAdditionalString;
188
191
195 G4int ProjectileResidualLambdaNumber; // Number of (anti-)lambdas for projectile (anti-)hypernucleus
197
202
210};
211
214}
215
218}
219
222}
223
224inline void G4FTFModel::SetImpactParameter( const G4double b_value ) {
225 Bimpact = b_value;
226}
227
229 return Bimpact;
230}
231
232inline void G4FTFModel::SetBminBmax( const G4double bmin_value, const G4double bmax_value ) {
233 BinInterval = false;
234 if ( bmin_value < 0.0 || bmax_value < 0.0 || bmax_value < bmin_value ) return;
235 BinInterval = true;
236 Bmin = bmin_value;
237 Bmax = bmax_value;
238}
239
241 return BinInterval;
242}
243
245 return Bmin;
246}
247
249 return Bmax;
250}
251
254}
255
258}
259
262}
263
264#endif
std::vector< G4ExcitedString * > G4ExcitedStringVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
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:3020
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:2790
G4bool SampleBinInterval() const
Definition: G4FTFModel.hh:240
G4DiffractiveExcitation * theExcitation
Definition: G4FTFModel.hh:183
G4Nucleon * TheInvolvedNucleonsOfProjectile[250]
Definition: G4FTFModel.hh:179
G4ElasticHNScattering * theElastic
Definition: G4FTFModel.hh:184
const G4FTFModel & operator=(const G4FTFModel &right)=delete
G4bool ComputeNucleusProperties(G4V3DNucleus *nucleus, G4LorentzVector &nucleusMomentum, G4LorentzVector &residualMomentum, G4double &sumMasses, G4double &residualExcitationEnergy, G4double &residualMass, G4int &residualMassNumber, G4int &residualCharge)
Definition: G4FTFModel.cc:2631
G4double LowEnergyLimit
Definition: G4FTFModel.hh:189
void SetImpactParameter(const G4double b_value)
Definition: G4FTFModel.hh:224
G4bool GenerateDeltaIsobar(const G4double sqrtS, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4double &sumMasses)
Definition: G4FTFModel.cc:2725
G4V3DNucleus * GetTargetNucleus() const
Definition: G4FTFModel.hh:216
std::vector< G4VSplitableHadron * > theAdditionalString
Definition: G4FTFModel.hh:187
G4int NumberOfTargetSpectatorNucleons
Definition: G4FTFModel.hh:208
G4int NumberOfInvolvedNucleonsOfProjectile
Definition: G4FTFModel.hh:180
G4double ProjectileResidualExcitationEnergy
Definition: G4FTFModel.hh:196
G4FTFModel(const G4String &modelName="FTF")
Definition: G4FTFModel.cc:70
G4LorentzVector TargetResidual4Momentum
Definition: G4FTFModel.hh:198
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
Definition: G4FTFModel.cc:2609
G4int NumberOfNNcollisions
Definition: G4FTFModel.hh:209
G4int ProjectileResidualCharge
Definition: G4FTFModel.hh:194
G4Nucleon * TheInvolvedNucleonsOfTarget[250]
Definition: G4FTFModel.hh:176
G4bool HighEnergyInter
Definition: G4FTFModel.hh:190
G4int GetNumberOfTargetSpectatorNucleons() const
Definition: G4FTFModel.hh:256
G4ExcitedStringVector * GetStrings() override
Definition: G4FTFModel.cc:298
G4ReactionProduct theProjectile
Definition: G4FTFModel.hh:173
G4bool operator!=(const G4FTFModel &right) const =delete
G4int ProjectileResidualMassNumber
Definition: G4FTFModel.hh:193
void BuildStrings(G4ExcitedStringVector *strings)
Definition: G4FTFModel.cc:1975
G4bool PutOnMassShell()
Definition: G4FTFModel.cc:559
G4int AdjustNucleonsAlgorithm_beforeSampling(G4int interactionCase, G4VSplitableHadron *SelectedAntiBaryon, G4Nucleon *ProjectileNucleon, G4VSplitableHadron *SelectedTargetNucleon, G4Nucleon *TargetNucleon, G4bool Annihilation, CommonVariables &common)
Definition: G4FTFModel.cc:1148
void ReggeonCascade()
Definition: G4FTFModel.cc:456
void SetBminBmax(const G4double bmin_value, const G4double bmax_value)
Definition: G4FTFModel.hh:232
G4FTFAnnihilation * theAnnihilation
Definition: G4FTFModel.hh:185
G4int NumberOfInvolvedNucleonsOfTarget
Definition: G4FTFModel.hh:177
G4bool ExciteParticipants()
Definition: G4FTFModel.cc:847
G4V3DNucleus * GetWoundedNucleus() const override
Definition: G4FTFModel.hh:212
G4double Bimpact
Definition: G4FTFModel.hh:203
G4double Bmin
Definition: G4FTFModel.hh:205
G4FTFModel(const G4FTFModel &right)=delete
G4FTFParameters * theParameters
Definition: G4FTFModel.hh:182
G4int TargetResidualMassNumber
Definition: G4FTFModel.hh:199
void GetResiduals()
Definition: G4FTFModel.cc:2290
G4double GetImpactParameter() const
Definition: G4FTFModel.hh:228
void StoreInvolvedNucleon()
Definition: G4FTFModel.cc:405
G4FTFParticipants theParticipants
Definition: G4FTFModel.hh:174
~G4FTFModel() override
Definition: G4FTFModel.cc:122
G4bool BinInterval
Definition: G4FTFModel.hh:204
G4V3DNucleus * GetProjectileNucleus() const override
Definition: G4FTFModel.hh:220
void Init(const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile) override
Definition: G4FTFModel.cc:162
G4int GetNumberOfNNcollisions() const
Definition: G4FTFModel.hh:260
G4int ProjectileResidualLambdaNumber
Definition: G4FTFModel.hh:195
G4double GetBmin() const
Definition: G4FTFModel.hh:244
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:2946
void ModelDescription(std::ostream &) const override
Definition: G4FTFModel.cc:3092
G4bool AdjustNucleonsAlgorithm_Sampling(G4int interactionCase, CommonVariables &common)
Definition: G4FTFModel.cc:1501
G4LorentzVector ProjectileResidual4Momentum
Definition: G4FTFModel.hh:192
G4double TargetResidualExcitationEnergy
Definition: G4FTFModel.hh:201
void AdjustNucleonsAlgorithm_afterSampling(G4int interactionCase, G4VSplitableHadron *SelectedAntiBaryon, G4VSplitableHadron *SelectedTargetNucleon, CommonVariables &common)
Definition: G4FTFModel.cc:1836
G4int GetNumberOfProjectileSpectatorNucleons() const
Definition: G4FTFModel.hh:252
G4bool operator==(const G4FTFModel &right) const =delete
G4int TargetResidualCharge
Definition: G4FTFModel.hh:200
G4double Bmax
Definition: G4FTFModel.hh:206
G4int NumberOfProjectileSpectatorNucleons
Definition: G4FTFModel.hh:207
G4double GetBmax() const
Definition: G4FTFModel.hh:248
G4bool AdjustNucleons(G4VSplitableHadron *SelectedAntiBaryon, G4Nucleon *ProjectileNucleon, G4VSplitableHadron *SelectedTargetNucleon, G4Nucleon *TargetNucleon, G4bool Annihilation)
Definition: G4FTFModel.cc:1032
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual G4V3DNucleus * GetWoundedNucleus() const
G4LorentzVector Pprojectile
Definition: G4FTFModel.hh:120
G4LorentzRotation toLab
Definition: G4FTFModel.hh:121
G4LorentzVector TResidual4Momentum
Definition: G4FTFModel.hh:120
G4LorentzVector PResidual4Momentum
Definition: G4FTFModel.hh:120
G4LorentzRotation toCms
Definition: G4FTFModel.hh:121
static int FASTCALL common(PROLOG_STATE *state, int tok)
Definition: xmlrole.cc:1309