File indexing completed on 2025-01-18 09:58:18
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036 #ifndef G4FTFModel_h
0037 #define G4FTFModel_h 1
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047 #include "G4VPartonStringModel.hh"
0048 #include "G4FTFParameters.hh"
0049 #include "G4FTFParticipants.hh"
0050 #include "G4ExcitedStringVector.hh"
0051 #include "G4DiffractiveExcitation.hh"
0052 #include "G4ElasticHNScattering.hh"
0053 #include "G4FTFAnnihilation.hh"
0054 #include "G4Proton.hh"
0055 #include "G4Neutron.hh"
0056
0057 class G4VSplitableHadron;
0058 class G4ExcitedString;
0059
0060
0061 class G4FTFModel : public G4VPartonStringModel {
0062 public:
0063 G4FTFModel( const G4String& modelName = "FTF" );
0064 ~G4FTFModel() override;
0065
0066 G4V3DNucleus* GetTargetNucleus() const;
0067 G4V3DNucleus* GetWoundedNucleus() const override;
0068 G4V3DNucleus* GetProjectileNucleus() const override;
0069
0070 void ModelDescription( std::ostream& ) const override;
0071
0072 G4FTFModel( const G4FTFModel& right ) = delete;
0073 const G4FTFModel& operator=( const G4FTFModel& right ) = delete;
0074 G4bool operator==( const G4FTFModel& right ) const = delete;
0075 G4bool operator!=( const G4FTFModel& right ) const = delete;
0076
0077 void SetImpactParameter( const G4double b_value );
0078 G4double GetImpactParameter() const;
0079 void SetBminBmax( const G4double bmin_value, const G4double bmax_value );
0080 G4bool SampleBinInterval() const;
0081 G4double GetBmin() const;
0082 G4double GetBmax() const;
0083 G4int GetNumberOfProjectileSpectatorNucleons() const;
0084 G4int GetNumberOfTargetSpectatorNucleons() const;
0085 G4int GetNumberOfNNcollisions() const;
0086
0087 protected:
0088 void Init( const G4Nucleus& aNucleus,
0089 const G4DynamicParticle& aProjectile ) override;
0090 G4ExcitedStringVector* GetStrings() override;
0091
0092 private:
0093 void StoreInvolvedNucleon();
0094 void ReggeonCascade();
0095 G4bool PutOnMassShell();
0096 G4bool ExciteParticipants();
0097 void BuildStrings( G4ExcitedStringVector* strings );
0098 void GetResiduals();
0099
0100 G4bool AdjustNucleons( G4VSplitableHadron* SelectedAntiBaryon,
0101 G4Nucleon* ProjectileNucleon,
0102 G4VSplitableHadron* SelectedTargetNucleon,
0103 G4Nucleon* TargetNucleon,
0104 G4bool Annihilation );
0105
0106 struct CommonVariables {
0107 G4int TResidualMassNumber = 0, TResidualCharge = 0, PResidualMassNumber = 0,
0108 PResidualCharge = 0, PResidualLambdaNumber = 0;
0109 G4double SqrtS = 0.0, S = 0.0, SumMasses = 0.0,
0110 TResidualExcitationEnergy = 0.0, TResidualMass = 0.0, TNucleonMass = 0.0,
0111 PResidualExcitationEnergy = 0.0, PResidualMass = 0.0, PNucleonMass = 0.0,
0112 Mprojectile = 0.0, M2projectile = 0.0, Pzprojectile = 0.0, Eprojectile = 0.0,
0113 WplusProjectile = 0.0,
0114 Mtarget = 0.0, M2target = 0.0, Pztarget = 0.0, Etarget = 0.0, WminusTarget = 0.0,
0115 Mt2targetNucleon = 0.0, PztargetNucleon = 0.0, EtargetNucleon = 0.0,
0116 Mt2projectileNucleon = 0.0, PzprojectileNucleon = 0.0, EprojectileNucleon = 0.0,
0117 YtargetNucleus = 0.0, YprojectileNucleus = 0.0,
0118 XminusNucleon = 0.0, XplusNucleon = 0.0, XminusResidual = 0.0, XplusResidual = 0.0;
0119 G4ThreeVector PtNucleon, PtResidual, PtNucleonP, PtResidualP, PtNucleonT, PtResidualT;
0120 G4LorentzVector Psum, Pprojectile, Ptmp, Ptarget, TResidual4Momentum, PResidual4Momentum;
0121 G4LorentzRotation toCms, toLab;
0122 };
0123 G4int AdjustNucleonsAlgorithm_beforeSampling( G4int interactionCase,
0124 G4VSplitableHadron* SelectedAntiBaryon,
0125 G4Nucleon* ProjectileNucleon,
0126 G4VSplitableHadron* SelectedTargetNucleon,
0127 G4Nucleon* TargetNucleon,
0128 G4bool Annihilation,
0129 CommonVariables& common );
0130 G4bool AdjustNucleonsAlgorithm_Sampling( G4int interactionCase,
0131 CommonVariables& common );
0132 void AdjustNucleonsAlgorithm_afterSampling( G4int interactionCase,
0133 G4VSplitableHadron* SelectedAntiBaryon,
0134 G4VSplitableHadron* SelectedTargetNucleon,
0135 CommonVariables& common );
0136
0137 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
0138
0139 G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
0140 G4LorentzVector& residualMomentum, G4double& sumMasses,
0141 G4double& residualExcitationEnergy, G4double& residualMass,
0142 G4int& residualMassNumber, G4int& residualCharge );
0143
0144
0145 G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
0146 G4Nucleon* involvedNucleons[], G4double& sumMasses );
0147
0148
0149 G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
0150 G4double dCor, G4V3DNucleus* nucleus,
0151 const G4LorentzVector& pResidual,
0152 const G4double residualMass, const G4int residualMassNumber,
0153 const G4int numberOfInvolvedNucleons,
0154 G4Nucleon* involvedNucleons[], G4double& mass2 );
0155
0156
0157
0158 G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
0159 const G4double projectileMass2, const G4double targetMass2,
0160 const G4double nucleusY, const G4bool isProjectileNucleus,
0161 const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
0162 G4double& targetWminus, G4double& projectileWplus, G4bool& success );
0163
0164
0165 G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
0166 const G4LorentzRotation& boostFromCmsToLab,
0167 const G4double residualMass, const G4int residualMassNumber,
0168 const G4int numberOfInvolvedNucleons,
0169 G4Nucleon* involvedNucleons[],
0170 G4LorentzVector& residual4Momentum );
0171
0172
0173 G4ReactionProduct theProjectile;
0174 G4FTFParticipants theParticipants;
0175
0176 G4Nucleon* TheInvolvedNucleonsOfTarget[250];
0177 G4int NumberOfInvolvedNucleonsOfTarget;
0178
0179 G4Nucleon* TheInvolvedNucleonsOfProjectile[250];
0180 G4int NumberOfInvolvedNucleonsOfProjectile;
0181
0182 G4FTFParameters* theParameters;
0183 G4DiffractiveExcitation* theExcitation;
0184 G4ElasticHNScattering* theElastic;
0185 G4FTFAnnihilation* theAnnihilation;
0186
0187 std::vector< G4VSplitableHadron* > theAdditionalString;
0188
0189 G4double LowEnergyLimit;
0190 G4bool HighEnergyInter;
0191
0192 G4LorentzVector ProjectileResidual4Momentum;
0193 G4int ProjectileResidualMassNumber;
0194 G4int ProjectileResidualCharge;
0195 G4int ProjectileResidualLambdaNumber;
0196 G4double ProjectileResidualExcitationEnergy;
0197
0198 G4LorentzVector TargetResidual4Momentum;
0199 G4int TargetResidualMassNumber;
0200 G4int TargetResidualCharge;
0201 G4double TargetResidualExcitationEnergy;
0202
0203 G4double Bimpact;
0204 G4bool BinInterval;
0205 G4double Bmin;
0206 G4double Bmax;
0207 G4int NumberOfProjectileSpectatorNucleons;
0208 G4int NumberOfTargetSpectatorNucleons;
0209 G4int NumberOfNNcollisions;
0210 };
0211
0212 inline G4V3DNucleus* G4FTFModel::GetWoundedNucleus() const {
0213 return theParticipants.GetWoundedNucleus();
0214 }
0215
0216 inline G4V3DNucleus* G4FTFModel::GetTargetNucleus() const {
0217 return theParticipants.GetWoundedNucleus();
0218 }
0219
0220 inline G4V3DNucleus* G4FTFModel::GetProjectileNucleus() const {
0221 return theParticipants.GetProjectileNucleus();
0222 }
0223
0224 inline void G4FTFModel::SetImpactParameter( const G4double b_value ) {
0225 Bimpact = b_value;
0226 }
0227
0228 inline G4double G4FTFModel::GetImpactParameter() const {
0229 return Bimpact;
0230 }
0231
0232 inline void G4FTFModel::SetBminBmax( const G4double bmin_value, const G4double bmax_value ) {
0233 BinInterval = false;
0234 if ( bmin_value < 0.0 || bmax_value < 0.0 || bmax_value < bmin_value ) return;
0235 BinInterval = true;
0236 Bmin = bmin_value;
0237 Bmax = bmax_value;
0238 }
0239
0240 inline G4bool G4FTFModel::SampleBinInterval() const {
0241 return BinInterval;
0242 }
0243
0244 inline G4double G4FTFModel::GetBmin() const {
0245 return Bmin;
0246 }
0247
0248 inline G4double G4FTFModel::GetBmax() const {
0249 return Bmax;
0250 }
0251
0252 inline G4int G4FTFModel::GetNumberOfProjectileSpectatorNucleons() const {
0253 return NumberOfProjectileSpectatorNucleons;
0254 }
0255
0256 inline G4int G4FTFModel::GetNumberOfTargetSpectatorNucleons() const {
0257 return NumberOfTargetSpectatorNucleons;
0258 }
0259
0260 inline G4int G4FTFModel::GetNumberOfNNcollisions() const {
0261 return NumberOfNNcollisions;
0262 }
0263
0264 #endif