Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QGSParticipants.cc
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#include <utility>
27
28#include "G4QGSParticipants.hh"
29#include "globals.hh"
30#include "G4SystemOfUnits.hh"
31#include "G4LorentzVector.hh"
32
33// Class G4QGSParticipants
34
35// HPW Feb 1999
36// Promoting model parameters from local variables class properties
38
39G4QGSParticipants::G4QGSParticipants() : theDiffExcitaton(), //0.7*GeV, 250*MeV, 250*MeV),
40 ModelMode(SOFT),
41 //nCutMax(7),ThresholdParameter(0.45*GeV),
42 nCutMax(7),ThresholdParameter(0.000*GeV),
43 QGSMThreshold(3*GeV),theNucleonRadius(1.5*fermi)
44
45{
46}
47
49: G4VParticipants(),ModelMode(right.ModelMode), nCutMax(right.nCutMax),
50 ThresholdParameter(right.ThresholdParameter), QGSMThreshold(right.QGSMThreshold),
51 theNucleonRadius(right.theNucleonRadius)
52{
53}
54
55
57{
58}
59
61{
62
63 // Find the collisions and collition conditions
64 G4VSplitableHadron* aProjectile = SelectInteractions(thePrimary);
65
66 // now build the parton pairs. HPW
68
69 // soft collisions first HPW, ordering is vital
71
72 // the rest is diffractive HPW
74
75 // clean-up, if necessary
76 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
77 theInteractions.clear();
78 std::for_each(theTargets.begin(), theTargets.end(), DeleteSplitableHadron());
79 theTargets.clear();
80 delete aProjectile;
81}
82
84{
85 G4VSplitableHadron* aProjectile = new G4QGSMSplitableHadron(thePrimary, TRUE); // @@@ check the TRUE
86 G4PomeronCrossSection theProbability(thePrimary.GetDefinition()); // @@@ should be data member
87 G4double outerRadius = theNucleus->GetOuterRadius();
88 // Check reaction threshold
89
91 G4Nucleon * pNucleon = theNucleus->GetNextNucleon();
92 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
93 //--DEBUG--G4cout << " qgspart- " << aPrimaryMomentum << " # " << aPrimaryMomentum.mag()
94 //--DEBUG-- << pNucleon->Get4Momentum() << " # " << (pNucleon->Get4Momentum()).mag()<< G4endl;
95 G4double s_nucleus = (aPrimaryMomentum + pNucleon->Get4Momentum()).mag2();
96 G4double ThresholdMass = thePrimary.GetMass() + pNucleon->GetDefinition()->GetPDGMass();
98 if (sqr(ThresholdMass + ThresholdParameter) > s_nucleus)
99 {
101 //throw G4HadronicException(__FILE__, __LINE__, "Initial energy is too low. The 4-vectors of the input are inconsistant with the particle masses.");
102 }
103 if (sqr(ThresholdMass + QGSMThreshold) > s_nucleus) // thus only diffractive in cascade!
104 {
106 }
107
108 // first find the collisions HPW
109 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
110 theInteractions.clear();
111 G4int totalCuts = 0;
112
113 #ifdef debug_QGS
114 G4double eK = thePrimary.GetKineticEnergy()/GeV;
115 #endif
116 #ifdef debug_G4QGSParticipants
117 G4double impactUsed = 0;
118 G4LorentzVector intNuclMom;
119 #endif
120
121 while(theInteractions.size() == 0)
122 {
123 // choose random impact parameter HPW
124 std::pair<G4double, G4double> theImpactParameter;
125 theImpactParameter = theNucleus->ChooseImpactXandY(outerRadius+theNucleonRadius);
126 G4double impactX = theImpactParameter.first;
127 G4double impactY = theImpactParameter.second;
128
129 // loop over nucleons to find collisions
131 G4int nucleonCount = 0; // debug
133 #ifdef debug_G4QGSParticipants
134 intNuclMom=aPrimaryMomentum;
135 #endif
136 while( (pNucleon = theNucleus->GetNextNucleon()) )
137 {
138 if(totalCuts>1.5*thePrimary.GetKineticEnergy()/GeV)
139 {
140 break;
141 }
142 nucleonCount++; // debug
143 // Needs to be moved to Probability class @@@
144 G4LorentzVector nucleonMomentum=pNucleon->Get4Momentum();
145 nucleonMomentum.setE(nucleonMomentum.e()-pNucleon->GetBindingEnergy());
146 G4double s_nucleon = (aPrimaryMomentum + nucleonMomentum).mag2();
147 G4double Distance2 = sqr(impactX - pNucleon->GetPosition().x()) +
148 sqr(impactY - pNucleon->GetPosition().y());
149 G4double Probability = theProbability.GetInelasticProbability(s_nucleon, Distance2);
150 // test for inelastic collision
151 G4double rndNumber = G4UniformRand();
152 // ModelMode = DIFFRACTIVE;
153 if (Probability > rndNumber)
154 {
155 #ifdef debug_G4QGSParticipants
156 G4cout << "DEBUG p="<< Probability<<" r="<<rndNumber<<" d="<<std::sqrt(Distance2)<<G4endl;
157 G4cout << " qgspart+ " << aPrimaryMomentum << " # " << aPrimaryMomentum.mag()
158 << pNucleon->Get4Momentum() << " # " << (pNucleon->Get4Momentum()).mag()<< G4endl;
159 intNuclMom += nucleonMomentum;
160 #endif
161 pNucleon->SetMomentum(nucleonMomentum);
162 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(*pNucleon);
164 theTargets.push_back(aTarget);
165 pNucleon->Hit(aTarget);
166 if ((theProbability.GetDiffractiveProbability(s_nucleon, Distance2)/Probability > G4UniformRand()
167 &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ))
168 {
169 // diffractive interaction occurs
171 {
172 theSingleDiffExcitation.ExciteParticipants(aProjectile, aTarget);
173 }
174 else
175 {
176 theDiffExcitaton.ExciteParticipants(aProjectile, aTarget);
177 }
178 G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
179 aInteraction->SetTarget(aTarget);
180 theInteractions.push_back(aInteraction);
181 aInteraction->SetNumberOfDiffractiveCollisions(1);
182 totalCuts += 1;
183 }
184 else
185 {
186 // nondiffractive soft interaction occurs
187 // sample nCut+1 (cut Pomerons) pairs of strings can be produced
188 G4int nCut;
189 G4double * running = new G4double[nCutMax];
190 running[0] = 0;
191 for(nCut = 0; nCut < nCutMax; nCut++)
192 {
193 running[nCut] = theProbability.GetCutPomeronProbability(s_nucleon, Distance2, nCut + 1);
194 if(nCut!=0) running[nCut] += running[nCut-1];
195 }
196 G4double random = running[nCutMax-1]*G4UniformRand();
197 for(nCut = 0; nCut < nCutMax; nCut++)
198 {
199 if(running[nCut] > random) break;
200 }
201 delete [] running;
202 nCut = 0;
203 aTarget->IncrementCollisionCount(nCut+1);
204 aProjectile->IncrementCollisionCount(nCut+1);
205 G4InteractionContent * aInteraction = new G4InteractionContent(aProjectile);
206 aInteraction->SetTarget(aTarget);
207 aInteraction->SetNumberOfSoftCollisions(nCut+1);
208 theInteractions.push_back(aInteraction);
209 totalCuts += nCut+1;
210 #ifdef debug_G4QGSParticipants
211 impactUsed=Distance2;
212 #endif
213 }
214 }
215 }
216 #ifdef debug_G4QGSParticipants
217 G4cout << G4endl<<"NUCLEONCOUNT "<<nucleonCount<<G4endl;
218 G4cout << " Interact 4-Vect " << intNuclMom << G4endl;
219 #endif
220 }
221 #ifdef debug_G4QGSParticipants
222 G4cout << G4endl<<"CUTDEBUG "<< totalCuts <<G4endl;
223 G4cout << "Impact Parameter used = "<<impactUsed<<G4endl;
224 #endif
225 return aProjectile;
226}
227
229{
230 // remove the "G4PartonPair::PROJECTILE", etc., which are not necessary. @@@
231 unsigned int i;
232 for(i = 0; i < theInteractions.size(); i++)
233 {
234 G4InteractionContent* anIniteraction = theInteractions[i];
235 G4VSplitableHadron* aProjectile = anIniteraction->GetProjectile();
236 G4Parton* aParton = aProjectile->GetNextParton();
237 G4PartonPair * aPartonPair;
238 // projectile first HPW
239 if (aParton)
240 {
241 aPartonPair = new G4PartonPair(aParton, aProjectile->GetNextAntiParton(),
244 #ifdef debug_G4QGSPart_PDiffColl
245 G4cout << "DiffPair Pro " << aPartonPair->GetParton1()->GetPDGcode() << " "
246 << aPartonPair->GetParton1()->Get4Momentum() << " "
247 << aPartonPair->GetParton1()->GetX() << " " << G4endl;
248 G4cout << " " << aPartonPair->GetParton2()->GetPDGcode() << " "
249 << aPartonPair->GetParton2()->Get4Momentum() << " "
250 << aPartonPair->GetParton2()->GetX() << " " << G4endl;
251 #endif
252 thePartonPairs.push_back(aPartonPair);
253 }
254 // then target HPW
255 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
256 aParton = aTarget->GetNextParton();
257 if (aParton)
258 {
259 aPartonPair = new G4PartonPair(aParton, aTarget->GetNextAntiParton(),
262 #ifdef debug_G4QGSPart_PDiffColl
263 G4cout << "DiffPair Tgt " << aPartonPair->GetParton1()->GetPDGcode() << " "
264 << aPartonPair->GetParton1()->Get4Momentum() << " "
265 << aPartonPair->GetParton1()->GetX() << " " << G4endl;
266 G4cout << " " << aPartonPair->GetParton2()->GetPDGcode() << " "
267 << aPartonPair->GetParton2()->Get4Momentum() << " "
268 << aPartonPair->GetParton2()->GetX() << " " << G4endl;
269 #endif
270 thePartonPairs.push_back(aPartonPair);
271 }
272 }
273}
274
276{
277 std::vector<G4InteractionContent*>::iterator i;
278 G4LorentzVector str4Mom;
279 i = theInteractions.begin();
280 while ( i != theInteractions.end() )
281 {
282 G4InteractionContent* anIniteraction = *i;
283 G4PartonPair * aPair = NULL;
284 if (anIniteraction->GetNumberOfSoftCollisions())
285 {
286 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
287 G4VSplitableHadron* pTarget = anIniteraction->GetTarget();
288 for (G4int j = 0; j < anIniteraction->GetNumberOfSoftCollisions(); j++)
289 {
290 aPair = new G4PartonPair(pTarget->GetNextParton(), pProjectile->GetNextAntiParton(),
292 #ifdef debug_G4QGSPart_PSoftColl
293 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
294 << aPair->GetParton1()->Get4Momentum() << " "
295 << aPair->GetParton1()->GetX() << " " << G4endl;
296 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
297 << aPair->GetParton2()->Get4Momentum() << " "
298 << aPair->GetParton2()->GetX() << " " << G4endl;
299 #endif
300 #ifdef debug_G4QGSParticipants
301 str4Mom += aPair->GetParton1()->Get4Momentum();
302 str4Mom += aPair->GetParton2()->Get4Momentum();
303 #endif
304 thePartonPairs.push_back(aPair);
305 aPair = new G4PartonPair(pProjectile->GetNextParton(), pTarget->GetNextAntiParton(),
307 #ifdef debug_G4QGSPart_PSoftColl
308 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
309 << aPair->GetParton1()->Get4Momentum() << " "
310 << aPair->GetParton1()->GetX() << " " << G4endl;
311 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
312 << aPair->GetParton2()->Get4Momentum() << " "
313 << aPair->GetParton2()->GetX() << " " << G4endl;
314 #endif
315 #ifdef debug_G4QGSParticipants
316 str4Mom += aPair->GetParton1()->Get4Momentum();
317 str4Mom += aPair->GetParton2()->Get4Momentum();
318 #endif
319 thePartonPairs.push_back(aPair);
320 }
321 delete *i;
322 i=theInteractions.erase(i); // i now points to the next interaction
323 } else {
324 i++;
325 }
326 }
327 #ifdef debug_G4QGSPart_PSoftColl
328 G4cout << " string 4 mom " << str4Mom << G4endl;
329 #endif
330}
G4int G4QGSParticipants_NPart
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double x() const
double y() const
void SetNumberOfDiffractiveCollisions(int)
G4VSplitableHadron * GetProjectile() const
void SetTarget(G4VSplitableHadron *aTarget)
G4VSplitableHadron * GetTarget() const
void SetMomentum(G4LorentzVector &aMomentum)
Definition: G4Nucleon.hh:70
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:90
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
virtual G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
G4Parton * GetParton2(void)
Definition: G4PartonPair.hh:79
G4Parton * GetParton1(void)
Definition: G4PartonPair.hh:74
G4double GetX()
Definition: G4Parton.hh:86
G4int GetPDGcode() const
Definition: G4Parton.hh:124
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:140
G4double GetInelasticProbability(const G4double s, const G4double impactsquare)
G4double GetDiffractiveProbability(const G4double s, const G4double impactsquare)
G4double GetCutPomeronProbability(const G4double s, const G4double impactsquare, const G4int nPomerons)
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
const G4double ThresholdParameter
G4SingleDiffractiveExcitation theSingleDiffExcitation
std::vector< G4InteractionContent * > theInteractions
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
G4QGSDiffractiveExcitation theDiffExcitaton
virtual ~G4QGSParticipants()
const G4double QGSMThreshold
const G4double theNucleonRadius
void BuildInteractions(const G4ReactionProduct &thePrimary)
std::vector< G4PartonPair * > thePartonPairs
std::vector< G4VSplitableHadron * > theTargets
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetMass() const
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
virtual G4double GetOuterRadius()=0
virtual G4Nucleon * GetNextNucleon()=0
virtual G4bool StartLoop()=0
std::pair< G4double, G4double > ChooseImpactXandY(G4double maxImpact)
Definition: G4V3DNucleus.hh:87
G4V3DNucleus * theNucleus
virtual G4Parton * GetNextParton()=0
virtual G4Parton * GetNextAntiParton()=0
void IncrementCollisionCount(G4int aCount)
#define TRUE
Definition: globals.hh:55
T sqr(const T &x)
Definition: templates.hh:145