Geant4 Cross Reference |
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 47 #include "G4VPartonStringModel.hh" 48 #include "G4FTFParameters.hh" 49 #include "G4FTFParticipants.hh" 50 #include "G4ExcitedStringVector.hh" 51 #include "G4DiffractiveExcitation.hh" 52 #include "G4ElasticHNScattering.hh" 53 #include "G4FTFAnnihilation.hh" 54 #include "G4Proton.hh" 55 #include "G4Neutron.hh" 56 57 class G4VSplitableHadron; 58 class G4ExcitedString; 59 60 61 class G4FTFModel : public G4VPartonStringModel { 62 public: 63 G4FTFModel( const G4String& modelName = "FTF" ); 64 ~G4FTFModel() override; 65 66 G4V3DNucleus* GetTargetNucleus() const; 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 ); 78 G4double GetImpactParameter() const; 79 void SetBminBmax( const G4double bmin_value, const G4double bmax_value ); 80 G4bool SampleBinInterval() const; 81 G4double GetBmin() const; 82 G4double GetBmax() const; 83 G4int GetNumberOfProjectileSpectatorNucleons() const; 84 G4int GetNumberOfTargetSpectatorNucleons() const; 85 G4int GetNumberOfNNcollisions() const; 86 87 protected: 88 void Init( const G4Nucleus& aNucleus, 89 const G4DynamicParticle& aProjectile ) override; 90 G4ExcitedStringVector* GetStrings() override; 91 92 private: 93 void StoreInvolvedNucleon(); 94 void ReggeonCascade(); 95 G4bool PutOnMassShell(); 96 G4bool ExciteParticipants(); 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: 106 struct CommonVariables { 107 G4int TResidualMassNumber = 0, TResidualCharge = 0, PResidualMassNumber = 0, 108 PResidualCharge = 0, PResidualLambdaNumber = 0; 109 G4double SqrtS = 0.0, S = 0.0, SumMasses = 0.0, 110 TResidualExcitationEnergy = 0.0, TResidualMass = 0.0, TNucleonMass = 0.0, 111 PResidualExcitationEnergy = 0.0, PResidualMass = 0.0, PNucleonMass = 0.0, 112 Mprojectile = 0.0, M2projectile = 0.0, Pzprojectile = 0.0, Eprojectile = 0.0, 113 WplusProjectile = 0.0, 114 Mtarget = 0.0, M2target = 0.0, Pztarget = 0.0, Etarget = 0.0, WminusTarget = 0.0, 115 Mt2targetNucleon = 0.0, PztargetNucleon = 0.0, EtargetNucleon = 0.0, 116 Mt2projectileNucleon = 0.0, PzprojectileNucleon = 0.0, EprojectileNucleon = 0.0, 117 YtargetNucleus = 0.0, YprojectileNucleus = 0.0, 118 XminusNucleon = 0.0, XplusNucleon = 0.0, XminusResidual = 0.0, XplusResidual = 0.0; 119 G4ThreeVector PtNucleon, PtResidual, PtNucleonP, PtResidualP, PtNucleonT, PtResidualT; 120 G4LorentzVector Psum, Pprojectile, Ptmp, Ptarget, TResidual4Momentum, PResidual4Momentum; 121 G4LorentzRotation toCms, toLab; 122 }; 123 G4int AdjustNucleonsAlgorithm_beforeSampling( G4int interactionCase, 124 G4VSplitableHadron* SelectedAntiBaryon, 125 G4Nucleon* ProjectileNucleon, 126 G4VSplitableHadron* SelectedTargetNucleon, 127 G4Nucleon* TargetNucleon, 128 G4bool Annihilation, 129 CommonVariables& common ); 130 G4bool AdjustNucleonsAlgorithm_Sampling( G4int interactionCase, 131 CommonVariables& common ); 132 void AdjustNucleonsAlgorithm_afterSampling( G4int interactionCase, 133 G4VSplitableHadron* SelectedAntiBaryon, 134 G4VSplitableHadron* SelectedTargetNucleon, 135 CommonVariables& common ); 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 173 G4ReactionProduct theProjectile; 174 G4FTFParticipants theParticipants; 175 176 G4Nucleon* TheInvolvedNucleonsOfTarget[250]; 177 G4int NumberOfInvolvedNucleonsOfTarget; 178 179 G4Nucleon* TheInvolvedNucleonsOfProjectile[250]; 180 G4int NumberOfInvolvedNucleonsOfProjectile; 181 182 G4FTFParameters* theParameters; 183 G4DiffractiveExcitation* theExcitation; 184 G4ElasticHNScattering* theElastic; 185 G4FTFAnnihilation* theAnnihilation; 186 187 std::vector< G4VSplitableHadron* > theAdditionalString; 188 189 G4double LowEnergyLimit; 190 G4bool HighEnergyInter; 191 192 G4LorentzVector ProjectileResidual4Momentum; 193 G4int ProjectileResidualMassNumber; 194 G4int ProjectileResidualCharge; 195 G4int ProjectileResidualLambdaNumber; // Number of (anti-)lambdas for projectile (anti-)hypernucleus 196 G4double ProjectileResidualExcitationEnergy; 197 198 G4LorentzVector TargetResidual4Momentum; 199 G4int TargetResidualMassNumber; 200 G4int TargetResidualCharge; 201 G4double TargetResidualExcitationEnergy; 202 203 G4double Bimpact; 204 G4bool BinInterval; 205 G4double Bmin; 206 G4double Bmax; 207 G4int NumberOfProjectileSpectatorNucleons; 208 G4int NumberOfTargetSpectatorNucleons; 209 G4int NumberOfNNcollisions; 210 }; 211 212 inline G4V3DNucleus* G4FTFModel::GetWoundedNucleus() const { 213 return theParticipants.GetWoundedNucleus(); 214 } 215 216 inline G4V3DNucleus* G4FTFModel::GetTargetNucleus() const { 217 return theParticipants.GetWoundedNucleus(); 218 } 219 220 inline G4V3DNucleus* G4FTFModel::GetProjectileNucleus() const { 221 return theParticipants.GetProjectileNucleus(); 222 } 223 224 inline void G4FTFModel::SetImpactParameter( const G4double b_value ) { 225 Bimpact = b_value; 226 } 227 228 inline G4double G4FTFModel::GetImpactParameter() const { 229 return Bimpact; 230 } 231 232 inline 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 240 inline G4bool G4FTFModel::SampleBinInterval() const { 241 return BinInterval; 242 } 243 244 inline G4double G4FTFModel::GetBmin() const { 245 return Bmin; 246 } 247 248 inline G4double G4FTFModel::GetBmax() const { 249 return Bmax; 250 } 251 252 inline G4int G4FTFModel::GetNumberOfProjectileSpectatorNucleons() const { 253 return NumberOfProjectileSpectatorNucleons; 254 } 255 256 inline G4int G4FTFModel::GetNumberOfTargetSpectatorNucleons() const { 257 return NumberOfTargetSpectatorNucleons; 258 } 259 260 inline G4int G4FTFModel::GetNumberOfNNcollisions() const { 261 return NumberOfNNcollisions; 262 } 263 264 #endif 265