Warning, file /include/Geant4/G4QGSParticipants.hh was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
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 #ifndef G4QGSParticipants_h
0027 #define G4QGSParticipants_h 1
0028
0029 #include "Randomize.hh"
0030 #include "G4VParticipants.hh"
0031 #include "G4Nucleon.hh"
0032 #include "G4InteractionContent.hh"
0033 #include "G4QGSDiffractiveExcitation.hh"
0034 #include "G4SingleDiffractiveExcitation.hh"
0035 #include "G4PartonPair.hh"
0036 #include "G4QGSMSplitableHadron.hh"
0037 #include "G4V3DNucleus.hh"
0038
0039 #include "G4VSplitableHadron.hh"
0040
0041 #include "G4Reggeons.hh"
0042 #include "G4QuarkExchange.hh"
0043
0044 class G4QGSParticipants : public G4VParticipants
0045 {
0046 public:
0047 G4QGSParticipants();
0048 G4QGSParticipants(const G4QGSParticipants &right);
0049 const G4QGSParticipants & operator=(const G4QGSParticipants &right);
0050 virtual ~G4QGSParticipants();
0051
0052 G4bool operator==(const G4QGSParticipants &right) const;
0053 G4bool operator!=(const G4QGSParticipants &right) const;
0054
0055 virtual void DoLorentzBoost(G4ThreeVector aBoost)
0056 {
0057 theCurrentVelocity = -aBoost;
0058 if(theNucleus) theNucleus->DoLorentzBoost(aBoost);
0059 theBoost = aBoost;
0060 }
0061
0062 G4PartonPair* GetNextPartonPair();
0063 void BuildInteractions(const G4ReactionProduct &thePrimary);
0064 void StartPartonPairLoop();
0065
0066 private:
0067 G4V3DNucleus* GetTargetNucleus() const;
0068 G4V3DNucleus* GetProjectileNucleus() const;
0069
0070 void PrepareInitialState( const G4ReactionProduct& thePrimary );
0071 void GetList( const G4ReactionProduct& thePrimary );
0072
0073 void StoreInvolvedNucleon();
0074 void ReggeonCascade();
0075 G4bool PutOnMassShell();
0076 void GetResiduals();
0077
0078 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
0079
0080 G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
0081 G4LorentzVector& residualMomentum, G4double& sumMasses,
0082 G4double& residualExcitationEnergy, G4double& residualMass,
0083 G4int& residualMassNumber, G4int& residualCharge );
0084
0085
0086 G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
0087 G4Nucleon* involvedNucleons[], G4double& sumMasses );
0088
0089 G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
0090 G4double dCor, G4V3DNucleus* nucleus,
0091 const G4LorentzVector& pResidual,
0092 const G4double residualMass, const G4int residualMassNumber,
0093 const G4int numberOfInvolvedNucleons,
0094 G4Nucleon* involvedNucleons[], G4double& mass2 );
0095
0096 G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
0097 const G4double projectileMass2, const G4double targetMass2,
0098 const G4double nucleusY, const G4bool isProjectileNucleus,
0099 const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
0100 G4double& targetWminus, G4double& projectileWplus, G4bool& success );
0101
0102 G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
0103 const G4LorentzRotation& boostFromCmsToLab,
0104 const G4double residualMass, const G4int residualMassNumber,
0105 const G4int numberOfInvolvedNucleons,
0106 G4Nucleon* involvedNucleons[],
0107 G4LorentzVector& residual4Momentum );
0108
0109 void CreateStrings();
0110
0111 private:
0112
0113 void SetCofNuclearDestruction( const G4double aValue );
0114 void SetR2ofNuclearDestruction( const G4double aValue );
0115
0116 void SetExcitationEnergyPerWoundedNucleon( const G4double aValue );
0117
0118 void SetDofNuclearDestruction( const G4double aValue );
0119 void SetPt2ofNuclearDestruction( const G4double aValue );
0120 void SetMaxPt2ofNuclearDestruction( const G4double aValue );
0121
0122
0123 G4double GetCofNuclearDestruction();
0124 G4double GetR2ofNuclearDestruction();
0125
0126 G4double GetExcitationEnergyPerWoundedNucleon();
0127
0128 G4double GetDofNuclearDestruction();
0129 G4double GetPt2ofNuclearDestruction();
0130 G4double GetMaxPt2ofNuclearDestruction();
0131
0132 protected:
0133 virtual G4VSplitableHadron* SelectInteractions(const G4ReactionProduct &thePrimary);
0134 void SplitHadrons();
0135 void PerformSoftCollisions();
0136 void PerformDiffractiveCollisions();
0137 G4bool DeterminePartonMomenta();
0138
0139 protected:
0140 struct DeleteInteractionContent {void operator()(G4InteractionContent*aC){delete aC;}};
0141 std::vector<G4InteractionContent*> theInteractions;
0142 struct DeleteSplitableHadron{void operator()(G4VSplitableHadron*aS){delete aS;}};
0143 std::vector<G4VSplitableHadron*> theTargets;
0144 struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
0145 std::vector<G4PartonPair*> thePartonPairs;
0146
0147 G4QuarkExchange theQuarkExchange;
0148 G4SingleDiffractiveExcitation theSingleDiffExcitation;
0149 G4QGSDiffractiveExcitation theDiffExcitaton;
0150 G4int ModelMode;
0151
0152 G4ThreeVector theBoost;
0153 G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta);
0154
0155 protected:
0156
0157 enum { SOFT, DIFFRACTIVE };
0158 enum { ALL, WITHOUT_R, NON_DIFF };
0159 enum { PrD, TrD, DD, NonD, Qexc };
0160
0161 const G4int nCutMax;
0162 const G4double ThresholdParameter;
0163 const G4double QGSMThreshold;
0164 const G4double theNucleonRadius;
0165
0166 G4ThreeVector theCurrentVelocity;
0167 G4QGSMSplitableHadron* theProjectileSplitable;
0168
0169 private:
0170 G4ReactionProduct theProjectile;
0171
0172 G4Reggeons* Regge;
0173 G4int InteractionMode;
0174
0175 G4double alpha;
0176 G4double beta;
0177
0178 G4double sigmaPt;
0179
0180 G4Nucleon* TheInvolvedNucleonsOfTarget[250];
0181 G4int NumberOfInvolvedNucleonsOfTarget;
0182
0183 G4Nucleon* TheInvolvedNucleonsOfProjectile[250];
0184 G4int NumberOfInvolvedNucleonsOfProjectile;
0185
0186 G4LorentzVector ProjectileResidual4Momentum;
0187 G4int ProjectileResidualMassNumber;
0188 G4int ProjectileResidualCharge;
0189 G4double ProjectileResidualExcitationEnergy;
0190
0191 G4LorentzVector TargetResidual4Momentum;
0192 G4int TargetResidualMassNumber;
0193 G4int TargetResidualCharge;
0194 G4double TargetResidualExcitationEnergy;
0195
0196 private:
0197
0198 G4double CofNuclearDestruction;
0199 G4double R2ofNuclearDestruction;
0200
0201 G4double ExcitationEnergyPerWoundedNucleon;
0202
0203 G4double DofNuclearDestruction;
0204 G4double Pt2ofNuclearDestruction;
0205 G4double MaxPt2ofNuclearDestruction;
0206 };
0207
0208 inline void G4QGSParticipants::StartPartonPairLoop()
0209 {
0210 }
0211
0212 inline G4PartonPair* G4QGSParticipants::GetNextPartonPair()
0213 {
0214 if (thePartonPairs.empty()) return 0;
0215 G4PartonPair * result = thePartonPairs.back();
0216 thePartonPairs.pop_back();
0217 return result;
0218 }
0219
0220 inline void G4QGSParticipants::SplitHadrons()
0221 {
0222 unsigned int i;
0223 for(i = 0; i < theInteractions.size(); i++)
0224 {
0225 theInteractions[i]->SplitHadrons();
0226 }
0227 }
0228
0229 inline G4V3DNucleus* G4QGSParticipants::GetTargetNucleus() const {
0230 return theNucleus;
0231 }
0232
0233 inline G4V3DNucleus* G4QGSParticipants::GetProjectileNucleus() const {
0234 return 0;
0235 }
0236
0237
0238
0239
0240 inline void G4QGSParticipants::SetCofNuclearDestruction( const G4double aValue ) {
0241 CofNuclearDestruction = aValue;
0242 }
0243
0244 inline void G4QGSParticipants::SetR2ofNuclearDestruction( const G4double aValue ) {
0245 R2ofNuclearDestruction = aValue;
0246 }
0247
0248 inline void G4QGSParticipants::SetExcitationEnergyPerWoundedNucleon( const G4double aValue ) {
0249 ExcitationEnergyPerWoundedNucleon = aValue;
0250 }
0251
0252 inline void G4QGSParticipants::SetDofNuclearDestruction( const G4double aValue ) {
0253 DofNuclearDestruction = aValue;
0254 }
0255
0256 inline void G4QGSParticipants::SetPt2ofNuclearDestruction( const G4double aValue ) {
0257 Pt2ofNuclearDestruction = aValue;
0258 }
0259
0260 inline void G4QGSParticipants::SetMaxPt2ofNuclearDestruction( const G4double aValue ) {
0261 MaxPt2ofNuclearDestruction = aValue;
0262 }
0263
0264
0265 inline G4double G4QGSParticipants::GetCofNuclearDestruction() {
0266 return CofNuclearDestruction;
0267 }
0268
0269 inline G4double G4QGSParticipants::GetR2ofNuclearDestruction() {
0270 return R2ofNuclearDestruction;
0271 }
0272
0273 inline G4double G4QGSParticipants::GetExcitationEnergyPerWoundedNucleon() {
0274 return ExcitationEnergyPerWoundedNucleon;
0275 }
0276
0277 inline G4double G4QGSParticipants::GetDofNuclearDestruction() {
0278 return DofNuclearDestruction;
0279 }
0280
0281 inline G4double G4QGSParticipants::GetPt2ofNuclearDestruction() {
0282 return Pt2ofNuclearDestruction;
0283 }
0284
0285 inline G4double G4QGSParticipants::GetMaxPt2ofNuclearDestruction() {
0286 return MaxPt2ofNuclearDestruction;
0287 }
0288
0289 #endif
0290