Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4QGSParticipants.hh
이 파일의 문서화 페이지로 가기
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 #ifndef G4QGSParticipants_h
27 #define G4QGSParticipants_h 1
28 
29 #include "Randomize.hh"
30 #include "G4VParticipants.hh"
31 #include "G4Nucleon.hh"
32 #include "G4InteractionContent.hh"
35 #include "G4PartonPair.hh"
36 #include "G4QGSMSplitableHadron.hh"
37 #include "G4V3DNucleus.hh"
38 
39 #include "G4VSplitableHadron.hh" // Uzhi
40 
41 #include "G4Reggeons.hh"
42 #include "G4QuarkExchange.hh"
43 
45 {
46 public:
49  const G4QGSParticipants & operator=(const G4QGSParticipants &right);
50  virtual ~G4QGSParticipants();
51 
52  int operator==(const G4QGSParticipants &right) const;
53  int operator!=(const G4QGSParticipants &right) const;
54 
55  virtual void DoLorentzBoost(G4ThreeVector aBoost)
56  {
57  theCurrentVelocity = -aBoost;
59  theBoost = aBoost;
60  }
61 
63  void BuildInteractions(const G4ReactionProduct &thePrimary);
64  void StartPartonPairLoop();
65 
66 private:
69 
70  void PrepareInitialState( const G4ReactionProduct& thePrimary );
71  void GetList( const G4ReactionProduct& thePrimary );
72 
73  void StoreInvolvedNucleon();
74  void ReggeonCascade();
76  void GetResiduals();
77 
78  G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
79 
80  G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum,
81  G4LorentzVector& residualMomentum, G4double& sumMasses,
82  G4double& residualExcitationEnergy, G4double& residualMass,
83  G4int& residualMassNumber, G4int& residualCharge );
84  // Utility methods used by PutOnMassShell.
85 
86  G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
87  G4Nucleon* involvedNucleons[], G4double& sumMasses );
88 
89  G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
90  G4double dCor, G4V3DNucleus* nucleus,
91  const G4LorentzVector& pResidual,
92  const G4double residualMass, const G4int residualMassNumber,
93  const G4int numberOfInvolvedNucleons,
94  G4Nucleon* involvedNucleons[], G4double& mass2 );
95 
96  G4bool CheckKinematics( const G4double sValue, const G4double sqrtS,
97  const G4double projectileMass2, const G4double targetMass2,
98  const G4double nucleusY, const G4bool isProjectileNucleus,
99  const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
100  G4double& targetWminus, G4double& projectileWplus, G4bool& success );
101 
102  G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus,
103  const G4LorentzRotation& boostFromCmsToLab,
104  const G4double residualMass, const G4int residualMassNumber,
105  const G4int numberOfInvolvedNucleons,
106  G4Nucleon* involvedNucleons[],
107  G4LorentzVector& residual4Momentum );
108 
109  void CreateStrings();
110 
111 private:
112  // Set parameters of nuclear destruction
113  void SetCofNuclearDestruction( const G4double aValue );
114  void SetR2ofNuclearDestruction( const G4double aValue );
115 
116  void SetExcitationEnergyPerWoundedNucleon( const G4double aValue );
117 
118  void SetDofNuclearDestruction( const G4double aValue );
119  void SetPt2ofNuclearDestruction( const G4double aValue );
120  void SetMaxPt2ofNuclearDestruction( const G4double aValue );
121 
122  // Get parameters of nuclear destruction
125 
127 
131 
132 protected:
133  virtual G4VSplitableHadron* SelectInteractions(const G4ReactionProduct &thePrimary);
134  void SplitHadrons();
135  void PerformSoftCollisions();
138 
139 protected:
141  std::vector<G4InteractionContent*> theInteractions;
143  std::vector<G4VSplitableHadron*> theTargets;
144  struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
145  std::vector<G4PartonPair*> thePartonPairs;
146 
147  G4QuarkExchange theQuarkExchange; // Uzhi 20 Oct. 2016
151 
153  G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta);
154 protected:
155  // model parameters HPW
156  enum { SOFT, DIFFRACTIVE };
157  enum { ALL, WITHOUT_R, NON_DIFF }; // Interaction modes
158  enum { PrD, TrD, DD, NonD, Qexc }; // Interaction types
159 
160  const G4int nCutMax;
164 
165  G4ThreeVector theCurrentVelocity; // Uzhi 17 Apr. 2015 Is it needed?
167 private:
169 
170  G4Reggeons* Regge; // Uzhi 18 Oct. 2016
172 
175 
177 
180 
183 
188 
193 
194 private:
195  // Parameters of nuclear destruction
196  G4double CofNuclearDestruction; // Cnd of nuclear destruction
198 
200 
201  G4double DofNuclearDestruction; // D for momentum sampling
204 
205 };
206 
208 {
209 }
210 
212 {
213  if (thePartonPairs.empty()) return 0;
215  thePartonPairs.pop_back();
216  return result;
217 }
218 
219 
221 {
222  unsigned int i;
223  for(i = 0; i < theInteractions.size(); i++)
224  {
225  theInteractions[i]->SplitHadrons();
226  }
227 }
228 //--------------------------------------
229 //Uzhi Copy from FTF Model.hh
230 /*
231 inline G4V3DNucleus* G4QGSParticipants::GetWoundedNucleus() const {
232  return theNucleus;
233 }
234 */
235 
237  return theNucleus;
238 }
239 
241  return 0;
242 }
243 
244 // Uzhi Start copy from FTFparameters
245 // Set parameters of nuclear destruction
246 
248  CofNuclearDestruction = aValue;
249 }
250 
252  R2ofNuclearDestruction = aValue;
253 }
254 
257 }
258 
260  DofNuclearDestruction = aValue;
261 }
262 
264  Pt2ofNuclearDestruction = aValue;
265 }
266 
269 }
270 
271 // Get parameters of nuclear destruction
273  return CofNuclearDestruction;
274 }
275 
277  return R2ofNuclearDestruction;
278 }
279 
282 }
283 
285  return DofNuclearDestruction;
286 }
287 
290 }
291 
294 }
295 //Uzhi End copy from FTFparameters
296 #endif
297 
298 
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)
G4double GetMaxPt2ofNuclearDestruction()
G4double GetExcitationEnergyPerWoundedNucleon()
G4QGSMSplitableHadron * theProjectileSplitable
void SetDofNuclearDestruction(const G4double aValue)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
G4bool FinalizeKinematics(const G4double w, const G4bool isProjectileNucleus, const G4LorentzRotation &boostFromCmsToLab, const G4double residualMass, const G4int residualMassNumber, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4LorentzVector &residual4Momentum)
G4int NumberOfInvolvedNucleonsOfProjectile
const G4double QGSMThreshold
std::vector< G4InteractionContent * > theInteractions
G4double GetDofNuclearDestruction()
G4double MaxPt2ofNuclearDestruction
int operator!=(const G4QGSParticipants &right) const
void SetPt2ofNuclearDestruction(const G4double aValue)
G4PartonPair * GetNextPartonPair()
G4ThreeVector theBoost
G4LorentzVector ProjectileResidual4Momentum
virtual ~G4QGSParticipants()
void BuildInteractions(const G4ReactionProduct &thePrimary)
G4double Pt2ofNuclearDestruction
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
G4LorentzVector TargetResidual4Momentum
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4V3DNucleus * GetProjectileNucleus() const
G4double ExcitationEnergyPerWoundedNucleon
G4Nucleon * TheInvolvedNucleonsOfTarget[250]
G4V3DNucleus * GetTargetNucleus() const
G4bool GenerateDeltaIsobar(const G4double sqrtS, const G4int numberOfInvolvedNucleons, G4Nucleon *involvedNucleons[], G4double &sumMasses)
void GetList(const G4ReactionProduct &thePrimary)
G4bool ComputeNucleusProperties(G4V3DNucleus *nucleus, G4LorentzVector &nucleusMomentum, G4LorentzVector &residualMomentum, G4double &sumMasses, G4double &residualExcitationEnergy, G4double &residualMass, G4int &residualMassNumber, G4int &residualCharge)
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
virtual void DoLorentzBoost(G4ThreeVector aBoost)
std::vector< G4VSplitableHadron * > theTargets
G4double G4ParticleHPJENDLHEData::G4double result
G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta)
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)
G4double TargetResidualExcitationEnergy
G4ThreeVector theCurrentVelocity
int G4int
Definition: G4Types.hh:78
G4double GetCofNuclearDestruction()
G4SingleDiffractiveExcitation theSingleDiffExcitation
G4int NumberOfInvolvedNucleonsOfTarget
G4Nucleon * TheInvolvedNucleonsOfProjectile[250]
int operator==(const G4QGSParticipants &right) const
G4double ProjectileResidualExcitationEnergy
void SetCofNuclearDestruction(const G4double aValue)
const G4double ThresholdParameter
G4double GetR2ofNuclearDestruction()
G4ReactionProduct theProjectile
void PrepareInitialState(const G4ReactionProduct &thePrimary)
G4QGSDiffractiveExcitation theDiffExcitaton
void operator()(G4VSplitableHadron *aS)
std::vector< G4PartonPair * > thePartonPairs
G4double GetPt2ofNuclearDestruction()
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
const G4double theNucleonRadius
void operator()(G4InteractionContent *aC)
G4QuarkExchange theQuarkExchange
const G4QGSParticipants & operator=(const G4QGSParticipants &right)
void SetR2ofNuclearDestruction(const G4double aValue)
G4V3DNucleus * theNucleus