Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QGSParticipants Class Reference

#include <G4QGSParticipants.hh>

+ Inheritance diagram for G4QGSParticipants:

Classes

struct  DeleteInteractionContent
 
struct  DeletePartonPair
 
struct  DeleteSplitableHadron
 

Public Member Functions

 G4QGSParticipants ()
 
 G4QGSParticipants (const G4QGSParticipants &right)
 
const G4QGSParticipantsoperator= (const G4QGSParticipants &right)
 
virtual ~G4QGSParticipants ()
 
G4bool operator== (const G4QGSParticipants &right) const
 
G4bool operator!= (const G4QGSParticipants &right) const
 
virtual void DoLorentzBoost (G4ThreeVector aBoost)
 
G4PartonPairGetNextPartonPair ()
 
void BuildInteractions (const G4ReactionProduct &thePrimary)
 
void StartPartonPairLoop ()
 
- Public Member Functions inherited from G4VParticipants
 G4VParticipants ()
 
virtual ~G4VParticipants ()
 
 G4VParticipants (const G4VParticipants &right)=delete
 
const G4VParticipantsoperator= (const G4VParticipants &right)=delete
 
G4bool operator== (const G4VParticipants &right) const =delete
 
G4bool operator!= (const G4VParticipants &right) const =delete
 
virtual void Init (G4int theZ, G4int theA)
 
virtual void SetNucleus (G4V3DNucleus *aNucleus)
 
virtual G4V3DNucleusGetWoundedNucleus () const
 
virtual void InitProjectileNucleus (G4int theZ, G4int theA, G4int numberOfLambdasOrAntiLambdas=0)
 
virtual void SetProjectileNucleus (G4V3DNucleus *aNucleus)
 

Protected Types

enum  { SOFT , DIFFRACTIVE }
 
enum  { ALL , WITHOUT_R , NON_DIFF }
 
enum  {
  PrD , TrD , DD , NonD ,
  Qexc
}
 

Protected Member Functions

virtual G4VSplitableHadronSelectInteractions (const G4ReactionProduct &thePrimary)
 
void SplitHadrons ()
 
void PerformSoftCollisions ()
 
void PerformDiffractiveCollisions ()
 
G4bool DeterminePartonMomenta ()
 
G4double SampleX (G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta)
 

Protected Attributes

std::vector< G4InteractionContent * > theInteractions
 
std::vector< G4VSplitableHadron * > theTargets
 
std::vector< G4PartonPair * > thePartonPairs
 
G4QuarkExchange theQuarkExchange
 
G4SingleDiffractiveExcitation theSingleDiffExcitation
 
G4QGSDiffractiveExcitation theDiffExcitaton
 
G4int ModelMode
 
G4ThreeVector theBoost
 
const G4int nCutMax
 
const G4double ThresholdParameter
 
const G4double QGSMThreshold
 
const G4double theNucleonRadius
 
G4ThreeVector theCurrentVelocity
 
G4QGSMSplitableHadrontheProjectileSplitable
 
- Protected Attributes inherited from G4VParticipants
G4V3DNucleustheNucleus
 
G4V3DNucleustheProjectileNucleus
 

Detailed Description

Definition at line 44 of file G4QGSParticipants.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
protected
Enumerator
SOFT 
DIFFRACTIVE 

Definition at line 157 of file G4QGSParticipants.hh.

◆ anonymous enum

anonymous enum
protected
Enumerator
ALL 
WITHOUT_R 
NON_DIFF 

Definition at line 158 of file G4QGSParticipants.hh.

◆ anonymous enum

anonymous enum
protected
Enumerator
PrD 
TrD 
DD 
NonD 
Qexc 

Definition at line 159 of file G4QGSParticipants.hh.

Constructor & Destructor Documentation

◆ G4QGSParticipants() [1/2]

G4QGSParticipants::G4QGSParticipants ( )

Definition at line 49 of file G4QGSParticipants.cc.

51 ThresholdParameter(0.0*GeV), QGSMThreshold(3.0*GeV),
53 theProjectileSplitable(nullptr), Regge(nullptr),
54 InteractionMode(ALL), alpha(-0.5), beta(2.5), sigmaPt(0.0),
55 NumberOfInvolvedNucleonsOfTarget(0), NumberOfInvolvedNucleonsOfProjectile(0),
56 ProjectileResidual4Momentum(G4LorentzVector()), ProjectileResidualMassNumber(0),
57 ProjectileResidualCharge(0), ProjectileResidualExcitationEnergy(0.0),
58 TargetResidual4Momentum(G4LorentzVector()), TargetResidualMassNumber(0),
59 TargetResidualCharge(0), TargetResidualExcitationEnergy(0.0),
60 CofNuclearDestruction(0.0), R2ofNuclearDestruction(0.0),
61 ExcitationEnergyPerWoundedNucleon(0.0), DofNuclearDestruction(0.0),
62 Pt2ofNuclearDestruction(0.0), MaxPt2ofNuclearDestruction(0.0)
63{
64 for (size_t i=0; i < 250; ++i) {
65 TheInvolvedNucleonsOfTarget[i] = nullptr;
66 TheInvolvedNucleonsOfProjectile[i] = nullptr;
67 }
68 // Parameters setting
69 SetCofNuclearDestruction( 1.0 );
70 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
71 SetDofNuclearDestruction( 0.3 );
72 SetPt2ofNuclearDestruction( 0.075*GeV*GeV );
73 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
74 SetExcitationEnergyPerWoundedNucleon( 40.0*MeV );
75
76 sigmaPt=0.25*sqr(GeV);
77}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
const G4double ThresholdParameter
G4ThreeVector theCurrentVelocity
G4QGSDiffractiveExcitation theDiffExcitaton
G4QGSMSplitableHadron * theProjectileSplitable
const G4double QGSMThreshold
const G4double theNucleonRadius
T sqr(const T &x)
Definition templates.hh:128

◆ G4QGSParticipants() [2/2]

G4QGSParticipants::G4QGSParticipants ( const G4QGSParticipants & right)

Definition at line 79 of file G4QGSParticipants.cc.

86 Regge(right.Regge), InteractionMode(right.InteractionMode),
87 alpha(right.alpha), beta(right.beta), sigmaPt(right.sigmaPt),
88 NumberOfInvolvedNucleonsOfTarget(right.NumberOfInvolvedNucleonsOfTarget),
89 NumberOfInvolvedNucleonsOfProjectile(right.NumberOfInvolvedNucleonsOfProjectile),
90 ProjectileResidual4Momentum(right.ProjectileResidual4Momentum),
91 ProjectileResidualMassNumber(right.ProjectileResidualMassNumber),
92 ProjectileResidualCharge(right.ProjectileResidualCharge),
93 ProjectileResidualExcitationEnergy(right.ProjectileResidualExcitationEnergy),
94 TargetResidual4Momentum(right.TargetResidual4Momentum),
95 TargetResidualMassNumber(right.TargetResidualMassNumber),
96 TargetResidualCharge(right.TargetResidualCharge),
97 TargetResidualExcitationEnergy(right.TargetResidualExcitationEnergy),
98 CofNuclearDestruction(right.CofNuclearDestruction),
99 R2ofNuclearDestruction(right.R2ofNuclearDestruction),
100 ExcitationEnergyPerWoundedNucleon(right.ExcitationEnergyPerWoundedNucleon),
101 DofNuclearDestruction(right.DofNuclearDestruction),
102 Pt2ofNuclearDestruction(right.Pt2ofNuclearDestruction),
103 MaxPt2ofNuclearDestruction(right.MaxPt2ofNuclearDestruction)
104{
105 for (size_t i=0; i < 250; ++i) {
106 TheInvolvedNucleonsOfTarget[i] = right.TheInvolvedNucleonsOfTarget[i];
107 TheInvolvedNucleonsOfProjectile[i] = right.TheInvolvedNucleonsOfProjectile[i];
108 }
109}

◆ ~G4QGSParticipants()

G4QGSParticipants::~G4QGSParticipants ( )
virtual

Definition at line 111 of file G4QGSParticipants.cc.

111{}

Member Function Documentation

◆ BuildInteractions()

void G4QGSParticipants::BuildInteractions ( const G4ReactionProduct & thePrimary)

Definition at line 113 of file G4QGSParticipants.cc.

114{
115 theProjectile = thePrimary;
116
117 Regge = new G4Reggeons(theProjectile.GetDefinition());
118
120
121 NumberOfInvolvedNucleonsOfProjectile= 0;
122 G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
123 ProjectileResidualMassNumber = 0;
124 ProjectileResidualCharge = 0;
125 ProjectileResidualExcitationEnergy = 0.0;
126 ProjectileResidual4Momentum = tmp;
127
128 NumberOfInvolvedNucleonsOfTarget= 0;
129 TargetResidualMassNumber = theNucleus->GetMassNumber();
130 TargetResidualCharge = theNucleus->GetCharge();
131 TargetResidualExcitationEnergy = 0.0;
132
134 G4Nucleon* NuclearNucleon;
135 while ( ( NuclearNucleon = theNucleus->GetNextNucleon() ) ) {
136 tmp+=NuclearNucleon->Get4Momentum();
137 }
138 TargetResidual4Momentum = tmp;
139
140 if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
141 // Projectile is a hadron : meson or baryon
142 ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
143 ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
144 ProjectileResidualExcitationEnergy = 0.0;
145 ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
146 ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
147 }
148
149 #ifdef debugQGSParticipants
150 G4cout <<G4endl<< "G4QGSParticipants::BuildInteractions---------" << G4endl
151 << "thePrimary " << thePrimary.GetDefinition()->GetParticleName()<<" "
152 <<ProjectileResidual4Momentum<<G4endl;
153 G4cout << "Target A and Z " << theNucleus->GetMassNumber() <<" "<<theNucleus->GetCharge()<<" "
154 << TargetResidual4Momentum<<G4endl;
155 #endif
156
157 G4bool Success( true );
158
159 const G4int maxNumberOfLoops = 1000;
160 G4int loopCounter = 0;
161 do
162 {
163 const G4int maxNumberOfInternalLoops = 1000;
164 G4int internalLoopCounter = 0;
165 do
166 {
167 if(std::abs(theProjectile.GetDefinition()->GetPDGEncoding()) < 100.0)
168 {
169 SelectInteractions(theProjectile); // for lepton projectile
170 }
171 else
172 {
173 GetList( theProjectile ); // Get list of participating nucleons for hadron projectile
174 }
175
176 if ( theInteractions.size() == 0 ) return;
177
178 StoreInvolvedNucleon(); // Store participating nucleon
179
180 ReggeonCascade(); // Make reggeon cascading. Involve nucleons in the cascading.
181
182 Success = PutOnMassShell(); // Ascribe momenta to the involved and participating nucleons
183
184 if(!Success) PrepareInitialState( thePrimary );
185
186 } while( (!Success) && ++internalLoopCounter < maxNumberOfInternalLoops );
187
188 if ( internalLoopCounter >= maxNumberOfInternalLoops ) {
189 Success = false;
190 }
191
192 if ( Success ) {
193 #ifdef debugQGSParticipants
194 G4cout<<G4endl<<"PerformDiffractiveCollisions(); if they happend." <<G4endl;
195 #endif
196
198
199 #ifdef debugQGSParticipants
200 G4cout<<G4endl<<"SplitHadrons();" <<G4endl;
201 #endif
202
203 SplitHadrons();
204
206 #ifdef debugQGSParticipants
207 G4cout<<"Perform non-Diffractive Collisions if they happend. Determine Parton Momenta and so on." <<G4endl;
208 #endif
209 Success = DeterminePartonMomenta();
210 }
211
212 if(!Success) PrepareInitialState( thePrimary );
213 }
214 } while( (!Success) && ++loopCounter < maxNumberOfLoops );
215
216 if ( loopCounter >= maxNumberOfLoops ) {
217 Success = false;
218 #ifdef debugQGSParticipants
219 G4cout<<"NOT Successful ======" <<G4endl;
220 #endif
221 }
222
223 if ( Success ) {
224 CreateStrings(); // To create strings
225
226 GetResiduals(); // To calculate residual nucleus properties
227
228 #ifdef debugQGSParticipants
229 G4cout<<"All O.K. ======" <<G4endl;
230 #endif
231 }
232
233 // clean-up, if necessary
234 #ifdef debugQGSParticipants
235 G4cout<<"Clearing "<<G4endl;
236 G4cout<<"theInteractions.size() "<<theInteractions.size()<<G4endl;
237 #endif
238
239 if( Regge ) delete Regge;
240
241 std::for_each( theInteractions.begin(), theInteractions.end(), DeleteInteractionContent() );
242 theInteractions.clear();
243
244 // Erasing of target involved nucleons.
245 #ifdef debugQGSParticipants
246 G4cout<<"Erasing of involved target nucleons "<<NumberOfInvolvedNucleonsOfTarget<<G4endl;
247 #endif
248
249 if ( NumberOfInvolvedNucleonsOfTarget != 0 ) {
250 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; i++ ) {
251 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfTarget[i]->GetSplitableHadron();
252 if ( (aNucleon != 0 ) && (aNucleon->GetStatus() >= 1) ) delete aNucleon;
253 }
254 }
255
256 // Erasing of projectile involved nucleons.
257 if ( NumberOfInvolvedNucleonsOfProjectile != 0 ) {
258 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfProjectile; i++ ) {
259 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfProjectile[i]->GetSplitableHadron();
260 if ( aNucleon ) delete aNucleon;
261 }
262 }
263
264 #ifdef debugQGSParticipants
265 G4cout<<"Delition of target nucleons from soft interactions "<<theTargets.size()
266 <<G4endl<<G4endl;
267 #endif
268 std::for_each(theTargets.begin(), theTargets.end(), DeleteSplitableHadron());
269 theTargets.clear();
270
274 }
275}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void setVect(const Hep3Vector &)
G4VSplitableHadron * GetSplitableHadron() const
Definition G4Nucleon.hh:97
virtual const G4LorentzVector & Get4Momentum() const
Definition G4Nucleon.hh:72
const G4String & GetParticleName() const
std::vector< G4InteractionContent * > theInteractions
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
std::vector< G4VSplitableHadron * > theTargets
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual G4int GetMassNumber()=0
virtual void SetProjectileNucleus(G4V3DNucleus *aNucleus)
G4V3DNucleus * theNucleus

◆ DeterminePartonMomenta()

G4bool G4QGSParticipants::DeterminePartonMomenta ( )
protected

Definition at line 1537 of file G4QGSParticipants.cc.

1538{
1539 if ( ! theProjectileSplitable ) return false;
1540
1541 const G4double aHugeValue = 1.0e+10;
1542
1543 #ifdef debugQGSParticipants
1544 G4cout<<G4endl<<"DeterminePartonMomenta()......"<<G4endl;
1545 G4cout<<"theProjectile status (0 -nondiffr, #0 diffr./reggeon): "<<theProjectileSplitable->GetStatus()<<G4endl;
1546 #endif
1547
1548 if (theProjectileSplitable->GetStatus() != 0) {return false;} // There were only diffractive interactions.
1549
1550 G4LorentzVector Projectile4Momentum = theProjectileSplitable->Get4Momentum();
1551 G4LorentzVector Psum = Projectile4Momentum;
1552
1553 G4double VqM_pr(0.), VaqM_pr(0.), VqM_tr(350.), VqqM_tr(700);
1554 if (std::abs(theProjectile.GetDefinition()->GetBaryonNumber()) != 0) {VqM_pr=350*MeV; VaqM_pr=700*MeV;}
1555
1556 #ifdef debugQGSParticipants
1557 G4cout<<"Projectile 4 momentum "<<Psum<<G4endl
1558 <<"Target nucleon momenta at start"<<G4endl;
1559 G4int NuclNo=0;
1560 #endif
1561
1562 std::vector<G4VSplitableHadron*>::iterator i;
1563
1564 for (i = theTargets.begin(); i != theTargets.end(); i++ )
1565 {
1566 Psum += (*i)->Get4Momentum();
1567 #ifdef debugQGSParticipants
1568 G4cout<<"Nusleus nucleon # and its 4Mom. "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1569 NuclNo++;
1570 #endif
1571 }
1572
1573 G4LorentzRotation toCms( -1*Psum.boostVector() );
1574
1575 G4LorentzVector Ptmp = toCms*Projectile4Momentum;
1576
1577 toCms.rotateZ( -1*Ptmp.phi() );
1578 toCms.rotateY( -1*Ptmp.theta() );
1579 G4LorentzRotation toLab(toCms.inverse());
1580 Projectile4Momentum.transform( toCms );
1581 // Ptarget.transform( toCms );
1582
1583 #ifdef debugQGSParticipants
1584 G4cout<<G4endl<<"In CMS---------------"<<G4endl;
1585 G4cout<<"Projectile 4 Mom "<<Projectile4Momentum<<G4endl;
1586 NuclNo=0;
1587 #endif
1588
1589 G4LorentzVector Target4Momentum(0.,0.,0.,0.);
1590 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1591 {
1592 G4LorentzVector tmp= (*i)->Get4Momentum(); tmp.transform( toCms );
1593 (*i)->Set4Momentum( tmp );
1594 #ifdef debugQGSParticipants
1595 G4cout<<"Target nucleon # and 4Mom "<<" "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1596 NuclNo++;
1597 #endif
1598 Target4Momentum += tmp;
1599 }
1600
1601 G4double S = Psum.mag2();
1602 G4double SqrtS = std::sqrt(S);
1603
1604 #ifdef debugQGSParticipants
1605 G4cout<<"Sum of target nucleons 4 momentum "<<Target4Momentum<<G4endl<<G4endl;
1606 G4cout<<"Target nucleons mom: px, py, z_1, m_i"<<G4endl;
1607 NuclNo=0;
1608 #endif
1609
1610 //G4double PplusProjectile = Projectile4Momentum.plus();
1611 G4double PminusTarget = Target4Momentum.minus();
1612
1613 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1614 {
1615 G4LorentzVector tmp = (*i)->Get4Momentum(); // tmp.boost(bstToCM);
1616
1617 //AR-19Jan2017 : the following line is causing a strange crash when Geant4
1618 // is built in optimized mode.
1619 // To fix it, I get the mass square instead of directly the
1620 // mass from the Lorentz vector, and then I take care of the
1621 // square root. If the mass square is negative, a JustWarning
1622 // exception is thrown, and the mass is set to 0.
1623 //G4double Mass = tmp.mag();
1624 G4double Mass2 = tmp.mag2();
1625 G4double Mass = 0.0;
1626 if ( Mass2 < 0.0 ) {
1628 ed << "Projectile " << theProjectile.GetDefinition()->GetParticleName()
1629 << "  4-momentum " << Psum << G4endl;
1630 ed << "LorentzVector tmp " << tmp << " with Mass2 " << Mass2 << G4endl;
1631 G4Exception( "G4QGSParticipants::DeterminePartonMomenta(): 4-momentum with negative mass!",
1632 "HAD_QGSPARTICIPANTS_001", JustWarning, ed );
1633 } else {
1634 Mass = std::sqrt( Mass2 );
1635 }
1636
1637 tmp.setPz(tmp.minus()/PminusTarget); tmp.setE(Mass);
1638 (*i)->Set4Momentum(tmp);
1639 #ifdef debugQGSParticipants
1640 G4cout<<"Target nucleons # and mom: "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1641 NuclNo++;
1642 #endif
1643 }
1644
1645 //+++++++++++++++++++++++++++++++++++++++++++
1646
1647 G4double SigPt = sigmaPt;
1648 G4Parton* aParton(0);
1649 G4ThreeVector aPtVector(0.,0.,0.);
1650 G4LorentzVector tmp(0.,0.,0.,0.);
1651
1652 G4double Mt(0.);
1653 G4double ProjSumMt(0.), ProjSumMt2perX(0.);
1654 G4double TargSumMt(0.), TargSumMt2perX(0.);
1655
1656
1657 G4double aBeta = beta; // Member of the class
1658
1659 const G4ParticleDefinition* theProjectileDefinition = theProjectileSplitable->GetDefinition();
1660 if (theProjectileDefinition == G4PionMinus::PionMinusDefinition()) aBeta = 1.;
1661 if (theProjectileDefinition == G4Gamma::GammaDefinition()) aBeta = 1.;
1662 if (theProjectileDefinition == G4PionPlus::PionPlusDefinition()) aBeta = 1.;
1663 if (theProjectileDefinition == G4PionZero::PionZeroDefinition()) aBeta = 1.;
1664 if (theProjectileDefinition == G4KaonPlus::KaonPlusDefinition()) aBeta = 0.;
1665 if (theProjectileDefinition == G4KaonMinus::KaonMinusDefinition()) aBeta = 0.;
1666
1667 G4double Xmin = 0.;
1668
1669 G4bool Success = true; G4int attempt = 0;
1670 const G4int maxNumberOfAttempts = 1000;
1671 do
1672 {
1673 attempt++; if( attempt == 100*(attempt/100) ) {SigPt/=2.;}
1674
1675 ProjSumMt=0.; ProjSumMt2perX=0.;
1676 TargSumMt=0.; TargSumMt2perX=0.;
1677
1678 Success = true;
1680 #ifdef debugQGSParticipants
1681 G4cout<<"attempt ------------------------ "<<attempt<<G4endl;
1682 G4cout<<"nSeaPair of proj "<<nSeaPair<<G4endl;
1683 #endif
1684
1685 G4double SumPx = 0.;
1686 G4double SumPy = 0.;
1687 G4double SumZ = 0.;
1688 G4int NumberOfUnsampledSeaQuarks = 2*nSeaPair;
1689
1690 G4double Qmass=0.;
1691 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1692 {
1693 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1694 #ifdef debugQGSParticipants
1695 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1696 #endif
1697 aPtVector = GaussianPt(SigPt, aHugeValue);
1698 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1699 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1700 Mt = std::sqrt(aPtVector.mag2()+sqr(Qmass));
1701 ProjSumMt += Mt;
1702
1703 // Sampling of Z fraction
1704 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1705 SumZ += tmp.z();
1706
1707 NumberOfUnsampledSeaQuarks--;
1708 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1709 tmp.setE(sqr(Mt));
1710 aParton->Set4Momentum(tmp);
1711
1712 aParton = theProjectileSplitable->GetNextAntiParton(); // for anti-quarks
1713 #ifdef debugQGSParticipants
1714 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1715 G4cout<<" "<<tmp<<" "<<SumZ<<G4endl;
1716 #endif
1717 aPtVector = GaussianPt(SigPt, aHugeValue);
1718 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1719 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1720 Mt = std::sqrt(aPtVector.mag2()+sqr(Qmass));
1721 ProjSumMt += Mt;
1722
1723 // Sampling of Z fraction
1724 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1725 SumZ += tmp.z();
1726
1727 NumberOfUnsampledSeaQuarks--;
1728 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1729 tmp.setE(sqr(Mt));
1730 aParton->Set4Momentum(tmp);
1731 #ifdef debugQGSParticipants
1732 G4cout<<" "<<tmp<<" "<<SumZ<<G4endl;
1733 #endif
1734 }
1735
1736 // For valence quark
1737 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1738 #ifdef debugQGSParticipants
1739 G4cout<<"Val quark of Pr"<<" "<<aParton->GetDefinition()->GetParticleName();
1740 #endif
1741 aPtVector = GaussianPt(SigPt, aHugeValue);
1742 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1743 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1744 Mt = std::sqrt(aPtVector.mag2()+sqr(VqM_pr));
1745 ProjSumMt += Mt;
1746
1747 // Sampling of Z fraction
1748 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1749 SumZ += tmp.z();
1750
1751 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1752 tmp.setE(sqr(Mt));
1753 aParton->Set4Momentum(tmp);
1754
1755 // For valence di-quark
1757 #ifdef debugQGSParticipants
1758 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1759 G4cout<<" "<<tmp<<" "<<SumZ<<" (z-fraction)"<<G4endl;
1760 #endif
1761 tmp.setPx(-SumPx); tmp.setPy(-SumPy);
1762 //Uzhi 2019 Mt = std::sqrt(aPtVector.mag2()+sqr(VaqM_pr));
1763 Mt = std::sqrt( sqr(SumPx) + sqr(SumPy) + sqr(VaqM_pr) ); //Uzhi 2019
1764 ProjSumMt += Mt;
1765 tmp.setPz(1.-SumZ);
1766
1767 ProjSumMt2perX +=sqr(Mt)/tmp.pz(); // QQmass=750 MeV
1768 tmp.setE(sqr(Mt));
1769 aParton->Set4Momentum(tmp);
1770 #ifdef debugQGSParticipants
1771 G4cout<<" "<<tmp<<" "<<SumZ+(1.-SumZ)<<" (z-fraction)"<<G4endl;
1772 NuclNo=0;
1773 #endif
1774
1775 // End of work with the projectile
1776
1777 // Work with target nucleons
1778
1779 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1780 {
1781 nSeaPair = (*i)->GetSoftCollisionCount()-1;
1782 #ifdef debugQGSParticipants
1783 G4cout<<"nSeaPair of target N "<<nSeaPair<<G4endl
1784 <<"Target nucleon 4Mom "<<(*i)->Get4Momentum()<<G4endl;
1785 #endif
1786
1787 SumPx = (*i)->Get4Momentum().px() * (-1.);
1788 SumPy = (*i)->Get4Momentum().py() * (-1.);
1789 SumZ = 0.;
1790
1791 NumberOfUnsampledSeaQuarks = 2*nSeaPair;
1792
1793 Qmass=0;
1794 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1795 {
1796 aParton = (*i)->GetNextParton(); // for quarks
1797 #ifdef debugQGSParticipants
1798 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1799 #endif
1800 aPtVector = GaussianPt(SigPt, aHugeValue);
1801 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1802 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1803 Mt=std::sqrt(aPtVector.mag2()+sqr(Qmass));
1804 TargSumMt += Mt;
1805
1806 // Sampling of Z fraction
1807 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1808 SumZ += tmp.z();
1809 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1810 NumberOfUnsampledSeaQuarks--;
1811 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1812 tmp.setE(sqr(Mt));
1813 aParton->Set4Momentum(tmp);
1814
1815 aParton = (*i)->GetNextAntiParton(); // for anti-quarks
1816 #ifdef debugQGSParticipants
1817 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1818 G4cout<<" "<<tmp<<" "<<SumZ<<G4endl;
1819 #endif
1820 aPtVector = GaussianPt(SigPt, aHugeValue);
1821 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1822 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1823 Mt=std::sqrt(aPtVector.mag2()+sqr(Qmass));
1824 TargSumMt += Mt;
1825
1826 // Sampling of Z fraction
1827 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1828 SumZ += tmp.z();
1829 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1830 NumberOfUnsampledSeaQuarks--;
1831 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1832 tmp.setE(sqr(Mt));
1833 aParton->Set4Momentum(tmp);
1834 #ifdef debugQGSParticipants
1835 G4cout<<" "<<tmp<<" "<<" "<<SumZ<<G4endl;
1836 #endif
1837 }
1838
1839 // Valence quark
1840 aParton = (*i)->GetNextParton(); // for quarks
1841 #ifdef debugQGSParticipants
1842 G4cout<<"Val quark of Tr"<<" "<<aParton->GetDefinition()->GetParticleName();
1843 #endif
1844 aPtVector = GaussianPt(SigPt, aHugeValue);
1845 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1846 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1847 Mt=std::sqrt(aPtVector.mag2()+sqr(VqM_tr));
1848 TargSumMt += Mt;
1849
1850 // Sampling of Z fraction
1851 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1852 SumZ += tmp.z();
1853 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1854 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1855 tmp.setE(sqr(Mt));
1856 aParton->Set4Momentum(tmp);
1857
1858 // Valence di-quark
1859 aParton = (*i)->GetNextAntiParton(); // for quarks
1860 #ifdef debugQGSParticipants
1861 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1862 G4cout<<" "<<tmp<<" "<<SumZ<<" (total z-sum) "<<G4endl;
1863 #endif
1864 tmp.setPx(-SumPx); tmp.setPy(-SumPy);
1865 //Uzhi 2019 Mt=std::sqrt(aPtVector.mag2()+sqr(VqqM_tr));
1866 Mt=std::sqrt( sqr(SumPx) + sqr(SumPy) + sqr(VqqM_tr) ); //Uzhi 2019
1867 TargSumMt += Mt;
1868
1869 tmp.setPz((*i)->Get4Momentum().pz()*(1.0 - SumZ));
1870 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1871 tmp.setE(sqr(Mt));
1872 aParton->Set4Momentum(tmp);
1873 #ifdef debugQGSParticipants
1874 G4cout<<" "<<tmp<<" "<<1.0<<" "<<(*i)->Get4Momentum().pz()<<G4endl;
1875 #endif
1876
1877 } // End of for(i = theTargets.begin(); i != theTargets.end(); i++ )
1878
1879 if( ProjSumMt + TargSumMt > SqrtS ) {
1880 Success = false; continue;}
1881 if( std::sqrt(ProjSumMt2perX) + std::sqrt(TargSumMt2perX) > SqrtS ) {
1882 Success = false; continue;}
1883
1884 } while( (!Success) &&
1885 attempt < maxNumberOfAttempts ); /* Loop checking, 07.08.2015, A.Ribon */
1886
1887 if ( attempt >= maxNumberOfAttempts ) {
1888 return false;
1889 }
1890
1891 //+++++++++++++++++++++++++++++++++++++++++++
1892
1893 G4double DecayMomentum2 = sqr(S) + sqr(ProjSumMt2perX) + sqr(TargSumMt2perX)
1894 - 2.0*S*ProjSumMt2perX - 2.0*S*TargSumMt2perX - 2.0*ProjSumMt2perX*TargSumMt2perX;
1895
1896 G4double targetWminus=( S - ProjSumMt2perX + TargSumMt2perX + std::sqrt( DecayMomentum2 ))/2.0/SqrtS;
1897 G4double projectileWplus = SqrtS - TargSumMt2perX/targetWminus;
1898
1899 G4LorentzVector Tmp(0.,0.,0.,0.);
1900 G4double z(0.);
1901
1903 #ifdef debugQGSParticipants
1904 G4cout<<"Backward transformation ===================="<<G4endl;
1905 G4cout<<"nSeaPair of proj "<<nSeaPair<<G4endl;
1906 #endif
1907
1908 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1909 {
1910 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1911 #ifdef debugQGSParticipants
1912 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1913 #endif
1914 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1915
1916 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1917 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1918 Tmp.transform( toLab );
1919
1920 aParton->Set4Momentum(Tmp);
1921
1922 aParton = theProjectileSplitable->GetNextAntiParton(); // for anti-quarks
1923 #ifdef debugQGSParticipants
1924 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1925 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1926 #endif
1927 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1928 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1929 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1930 Tmp.transform( toLab );
1931
1932 aParton->Set4Momentum(Tmp);
1933 #ifdef debugQGSParticipants
1934 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1935 #endif
1936 }
1937
1938 // For valence quark
1939 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1940 #ifdef debugQGSParticipants
1941 G4cout<<"Val quark of Pr"<<" "<<aParton->GetDefinition()->GetParticleName();
1942 #endif
1943 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1944 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1945 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1946 Tmp.transform( toLab );
1947
1948 aParton->Set4Momentum(Tmp);
1949
1950 // For valence di-quark
1952 #ifdef debugQGSParticipants
1953 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1954 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1955 #endif
1956 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1957 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1958 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1959 Tmp.transform( toLab );
1960
1961 aParton->Set4Momentum(Tmp);
1962
1963 #ifdef debugQGSParticipants
1964 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1965 NuclNo=0;
1966 #endif
1967
1968 // End of work with the projectile
1969
1970 // Work with target nucleons
1971 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1972 {
1973 nSeaPair = (*i)->GetSoftCollisionCount()-1;
1974 #ifdef debugQGSParticipants
1975 G4cout<<"nSeaPair of target and N# "<<nSeaPair<<" "<<NuclNo<<G4endl;
1976 NuclNo++;
1977 #endif
1978 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1979 {
1980 aParton = (*i)->GetNextParton(); // for quarks
1981 #ifdef debugQGSParticipants
1982 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1983 #endif
1984 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1985 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1986 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1987 Tmp.transform( toLab );
1988
1989 aParton->Set4Momentum(Tmp);
1990
1991 aParton = (*i)->GetNextAntiParton(); // for quarks
1992 #ifdef debugQGSParticipants
1993 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1994 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1995 #endif
1996 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1997 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1998 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1999 Tmp.transform( toLab );
2000
2001 aParton->Set4Momentum(Tmp);
2002 #ifdef debugQGSParticipants
2003 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
2004 #endif
2005 }
2006
2007 // Valence quark
2008
2009 aParton = (*i)->GetNextParton(); // for quarks
2010 #ifdef debugQGSParticipants
2011 G4cout<<"Val quark of Tr"<<" "<<aParton->GetDefinition()->GetParticleName();
2012 #endif
2013 Tmp =aParton->Get4Momentum(); z=Tmp.z();
2014 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
2015 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
2016 Tmp.transform( toLab );
2017
2018 aParton->Set4Momentum(Tmp);
2019
2020 // Valence di-quark
2021 aParton = (*i)->GetNextAntiParton(); // for quarks
2022 #ifdef debugQGSParticipants
2023 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
2024 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
2025 #endif
2026 Tmp =aParton->Get4Momentum(); z=Tmp.z();
2027 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
2028 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
2029 Tmp.transform( toLab );
2030
2031 aParton->Set4Momentum(Tmp);
2032 #ifdef debugQGSParticipants
2033 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
2034 NuclNo++;
2035 #endif
2036 } // End of for(i = theTargets.begin(); i != theTargets.end(); i++ )
2037
2038 return true;
2039}
G4double S(G4double temp)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
double theta() const
Hep3Vector boostVector() const
double minus() const
HepLorentzVector & rotateZ(double)
HepLorentzVector & rotateY(double)
HepLorentzVector & transform(const HepRotation &)
static G4Gamma * GammaDefinition()
Definition G4Gamma.cc:76
static G4KaonMinus * KaonMinusDefinition()
static G4KaonPlus * KaonPlusDefinition()
G4ParticleDefinition * GetDefinition()
Definition G4Parton.hh:161
static G4PionMinus * PionMinusDefinition()
static G4PionPlus * PionPlusDefinition()
Definition G4PionPlus.cc:88
static G4PionZero * PionZeroDefinition()
Definition G4PionZero.cc:97
virtual G4Parton * GetNextParton()
virtual G4Parton * GetNextAntiParton()
G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const

Referenced by BuildInteractions().

◆ DoLorentzBoost()

virtual void G4QGSParticipants::DoLorentzBoost ( G4ThreeVector aBoost)
inlinevirtual

Definition at line 55 of file G4QGSParticipants.hh.

56 {
57 theCurrentVelocity = -aBoost;
59 theBoost = aBoost;
60 }
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0

◆ GetNextPartonPair()

G4PartonPair * G4QGSParticipants::GetNextPartonPair ( )
inline

Definition at line 212 of file G4QGSParticipants.hh.

213{
214 if (thePartonPairs.empty()) return 0;
215 G4PartonPair * result = thePartonPairs.back();
216 thePartonPairs.pop_back();
217 return result;
218}
std::vector< G4PartonPair * > thePartonPairs

◆ operator!=()

G4bool G4QGSParticipants::operator!= ( const G4QGSParticipants & right) const

◆ operator=()

const G4QGSParticipants & G4QGSParticipants::operator= ( const G4QGSParticipants & right)

◆ operator==()

G4bool G4QGSParticipants::operator== ( const G4QGSParticipants & right) const

◆ PerformDiffractiveCollisions()

void G4QGSParticipants::PerformDiffractiveCollisions ( )
protected

Definition at line 1460 of file G4QGSParticipants.cc.

1461{
1462 #ifdef debugQGSParticipants
1463 G4cout<<G4endl<<"PerformDiffractiveCollisions()......"<<G4endl
1464 <<"theInteractions.size() "<<theInteractions.size()<<G4endl;
1465 #endif
1466
1467 unsigned int i;
1468 for (i = 0; i < theInteractions.size(); i++)
1469 {
1470 G4InteractionContent* anIniteraction = theInteractions[i];
1471 #ifdef debugQGSParticipants
1472 G4cout<<"Interaction # and its status "
1473 <<i<<" "<<theInteractions[i]->GetStatus()<<G4endl;
1474 #endif
1475
1476 G4int InterStatus = theInteractions[i]->GetStatus();
1477 if ( (InterStatus == PrD) || (InterStatus == TrD) || (InterStatus == DD))
1478 { // Selection of diffractive interactions
1479 #ifdef debugQGSParticipants
1480 G4cout<<"Simulation of diffractive interaction. "<<InterStatus <<" PrD/TrD/DD/ND/Qech - 0,1,2,3,4"<<G4endl;
1481 #endif
1482
1483 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1484
1485 #ifdef debugQGSParticipants
1486 G4cout<<"The proj. before inter "
1489 G4cout<<"The targ. before inter " <<aTarget->Get4Momentum()<<" "
1490 <<aTarget->Get4Momentum().mag()<<G4endl;
1491 #endif
1492
1493 if ( InterStatus == PrD )
1495
1496 if ( InterStatus == TrD )
1498
1499 if ( InterStatus == DD )
1501
1502 #ifdef debugQGSParticipants
1503 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1505 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1506 <<aTarget->Get4Momentum().mag()<<G4endl;
1507 #endif
1508 }
1509
1510 if ( InterStatus == Qexc )
1511 { // Quark exchange process
1512 #ifdef debugQGSParticipants
1513 G4cout<<"Simulation of interaction with quark exchange."<<G4endl;
1514 #endif
1515 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1516
1517 #ifdef debugQGSParticipants
1518 G4cout<<"The proj. before inter " <<theProjectileSplitable->Get4Momentum()<<" "
1520 G4cout<<"The targ. before inter "<<aTarget->Get4Momentum()<<" "
1521 <<aTarget->Get4Momentum().mag()<<G4endl;
1522 #endif
1523
1525
1526 #ifdef debugQGSParticipants
1527 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1529 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1530 <<aTarget->Get4Momentum().mag()<<G4endl;
1531 #endif
1532 }
1533 }
1534}
G4VSplitableHadron * GetTarget() const
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4bool ProjectileDiffraction=TRUE) const
G4QuarkExchange theQuarkExchange
G4SingleDiffractiveExcitation theSingleDiffExcitation
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4bool ProjectileDiffraction) const
#define TRUE
Definition globals.hh:41
#define FALSE
Definition globals.hh:38

Referenced by BuildInteractions().

◆ PerformSoftCollisions()

void G4QGSParticipants::PerformSoftCollisions ( )
protected

Definition at line 2331 of file G4QGSParticipants.cc.

2332{
2333 std::vector<G4InteractionContent*>::iterator i;
2334 G4LorentzVector str4Mom;
2335 i = theInteractions.begin();
2336 while ( i != theInteractions.end() ) /* Loop checking, 07.08.2015, A.Ribon */
2337 {
2338 G4InteractionContent* anIniteraction = *i;
2339 G4PartonPair * aPair = NULL;
2340 if (anIniteraction->GetNumberOfSoftCollisions())
2341 {
2342 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
2343 G4VSplitableHadron* pTarget = anIniteraction->GetTarget();
2344 for (G4int j = 0; j < anIniteraction->GetNumberOfSoftCollisions(); j++)
2345 {
2346 aPair = new G4PartonPair(pTarget->GetNextParton(), pProjectile->GetNextAntiParton(),
2348 #ifdef debugQGSParticipants
2349 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2350 << aPair->GetParton1()->Get4Momentum() << " "
2351 << aPair->GetParton1()->GetX() << " " << G4endl;
2352 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2353 << aPair->GetParton2()->Get4Momentum() << " "
2354 << aPair->GetParton2()->GetX() << " " << G4endl;
2355 #endif
2356 #ifdef debugQGSParticipants
2357 str4Mom += aPair->GetParton1()->Get4Momentum();
2358 str4Mom += aPair->GetParton2()->Get4Momentum();
2359 #endif
2360 thePartonPairs.push_back(aPair);
2361 aPair = new G4PartonPair(pProjectile->GetNextParton(), pTarget->GetNextAntiParton(),
2363 #ifdef debugQGSParticipants
2364 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2365 << aPair->GetParton1()->Get4Momentum() << " "
2366 << aPair->GetParton1()->GetX() << " " << G4endl;
2367 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2368 << aPair->GetParton2()->Get4Momentum() << " "
2369 << aPair->GetParton2()->GetX() << " " << G4endl;
2370 #endif
2371 #ifdef debugQGSParticipants
2372 str4Mom += aPair->GetParton1()->Get4Momentum();
2373 str4Mom += aPair->GetParton2()->Get4Momentum();
2374 #endif
2375 thePartonPairs.push_back(aPair);
2376 }
2377 delete *i;
2378 i=theInteractions.erase(i); // i now points to the next interaction
2379 } else {
2380 i++;
2381 }
2382 }
2383 #ifdef debugQGSParticipants
2384 G4cout << " string 4 mom " << str4Mom << G4endl;
2385 #endif
2386}
G4VSplitableHadron * GetProjectile() const
G4Parton * GetParton2(void)
G4Parton * GetParton1(void)
G4double GetX()
Definition G4Parton.hh:87
G4int GetPDGcode() const
Definition G4Parton.hh:127
const G4LorentzVector & Get4Momentum() const
Definition G4Parton.hh:143
virtual G4Parton * GetNextParton()=0
virtual G4Parton * GetNextAntiParton()=0

◆ SampleX()

G4double G4QGSParticipants::SampleX ( G4double anXmin,
G4int nSea,
G4int theTotalSea,
G4double aBeta )
protected

Definition at line 2042 of file G4QGSParticipants.cc.

2044{
2045 G4double Oalfa = 1./(alpha + 1.);
2046 G4double Obeta = 1./(aBeta + (alpha + 1.)*nSea + 1.); // ?
2047
2048 G4double Ksi1, Ksi2, r1, r2, r12;
2049 const G4int maxNumberOfLoops = 1000;
2050 G4int loopCounter = 0;
2051 do
2052 {
2053 Ksi1 = G4UniformRand(); r1 = G4Pow::GetInstance()->powA(Ksi1,Oalfa);
2054 Ksi2 = G4UniformRand(); r2 = G4Pow::GetInstance()->powA(Ksi2,Obeta);
2055 r12=r1+r2;
2056 } while( ( r12 > 1.) &&
2057 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
2058 if ( loopCounter >= maxNumberOfLoops ) {
2059 return 0.5; // Just an acceptable value, without any physics consideration.
2060 }
2061
2062 G4double result = r1/r12;
2063 return result;
2064}
#define G4UniformRand()
Definition Randomize.hh:52
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230

Referenced by DeterminePartonMomenta().

◆ SelectInteractions()

G4VSplitableHadron * G4QGSParticipants::SelectInteractions ( const G4ReactionProduct & thePrimary)
protectedvirtual

Definition at line 2390 of file G4QGSParticipants.cc.

2391{
2392 // Check reaction threshold - goes to CheckThreshold
2393
2396
2397 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
2398 G4LorentzVector aTargetNMomentum(0.,0.,0.,938.);
2399
2400 if ((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
2401 {
2402 throw G4HadronicException(__FILE__, __LINE__,
2403 "G4GammaParticipants::SelectInteractions: primary nan energy.");
2404 }
2405 G4double S = (aPrimaryMomentum + aTargetNMomentum).mag2();
2406 G4double ThresholdMass = thePrimary.GetMass() + 938.;
2407 ModelMode = SOFT;
2408
2409 #ifdef debugQGSParticipants
2410 G4cout << "Gamma projectile - SelectInteractions " << G4endl;
2411 G4cout<<"Energy and Nucleus "<<thePrimary.GetTotalEnergy()<<" "<<theNucleus->GetMassNumber()<<G4endl;
2412 G4cout << "SqrtS ThresholdMass ModelMode " <<std::sqrt(S)<<" "<<ThresholdMass<<" "<<ModelMode<< G4endl;
2413 #endif
2414
2415 if (sqr(ThresholdMass + ThresholdParameter) > S)
2416 {
2418 }
2419 if (sqr(ThresholdMass + QGSMThreshold) > S) // thus only diffractive in cascade!
2420 {
2422 }
2423
2424 // first find the collisions HPW
2425 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
2426 theInteractions.clear();
2427
2429 G4int NucleonNo=0;
2430
2432 G4Nucleon * pNucleon = 0;
2433
2434 while( (pNucleon = theNucleus->GetNextNucleon()) ) /* Loop checking, 07.08.2015, A.Ribon */
2435 { if(NucleonNo == theCurrent) break; NucleonNo++; }
2436
2437 if ( pNucleon ) {
2438 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(*pNucleon);
2439
2440 pNucleon->Hit(aTarget);
2441
2442 if ( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) )
2443 {
2446
2447 aInteraction->SetTarget(aTarget);
2448 aInteraction->SetTargetNucleon(pNucleon);
2449 aTarget->SetCollisionCount(0);
2450 aTarget->SetStatus(1);
2451
2452 aInteraction->SetNumberOfDiffractiveCollisions(1);
2453 aInteraction->SetNumberOfSoftCollisions(0);
2454 aInteraction->SetStatus(1);
2455
2456 theInteractions.push_back(aInteraction);
2457 }
2458 else
2459 {
2460 // nondiffractive soft interaction occurs
2461 aTarget->IncrementCollisionCount(1);
2462 aTarget->SetStatus(0);
2463 theTargets.push_back(aTarget);
2464
2467
2469 aInteraction->SetTarget(aTarget);
2470 aInteraction->SetTargetNucleon(pNucleon);
2471 aInteraction->SetNumberOfSoftCollisions(1);
2472 aInteraction->SetStatus(3);
2473 theInteractions.push_back(aInteraction);
2474 }
2475 }
2477}
void SetTargetNucleon(G4Nucleon *aNucleon)
void SetTarget(G4VSplitableHadron *aTarget)
void SetStatus(G4int aValue)
void Hit(G4VSplitableHadron *aHit)
Definition G4Nucleon.hh:91
G4double GetMass() const
void SetStatus(const G4int aStatus)
void SetCollisionCount(G4int aCount)
void IncrementCollisionCount(G4int aCount)

Referenced by BuildInteractions().

◆ SplitHadrons()

void G4QGSParticipants::SplitHadrons ( )
inlineprotected

Definition at line 220 of file G4QGSParticipants.hh.

221{
222 unsigned int i;
223 for(i = 0; i < theInteractions.size(); i++)
224 {
225 theInteractions[i]->SplitHadrons();
226 }
227}

Referenced by BuildInteractions().

◆ StartPartonPairLoop()

void G4QGSParticipants::StartPartonPairLoop ( )
inline

Definition at line 208 of file G4QGSParticipants.hh.

209{
210}

Member Data Documentation

◆ ModelMode

G4int G4QGSParticipants::ModelMode
protected

Definition at line 150 of file G4QGSParticipants.hh.

Referenced by SelectInteractions().

◆ nCutMax

const G4int G4QGSParticipants::nCutMax
protected

Definition at line 161 of file G4QGSParticipants.hh.

◆ QGSMThreshold

const G4double G4QGSParticipants::QGSMThreshold
protected

Definition at line 163 of file G4QGSParticipants.hh.

Referenced by SelectInteractions().

◆ theBoost

G4ThreeVector G4QGSParticipants::theBoost
protected

Definition at line 152 of file G4QGSParticipants.hh.

Referenced by DoLorentzBoost().

◆ theCurrentVelocity

G4ThreeVector G4QGSParticipants::theCurrentVelocity
protected

Definition at line 166 of file G4QGSParticipants.hh.

Referenced by DoLorentzBoost().

◆ theDiffExcitaton

G4QGSDiffractiveExcitation G4QGSParticipants::theDiffExcitaton
protected

Definition at line 149 of file G4QGSParticipants.hh.

Referenced by PerformDiffractiveCollisions().

◆ theInteractions

std::vector<G4InteractionContent*> G4QGSParticipants::theInteractions
protected

◆ theNucleonRadius

const G4double G4QGSParticipants::theNucleonRadius
protected

Definition at line 164 of file G4QGSParticipants.hh.

◆ thePartonPairs

std::vector<G4PartonPair*> G4QGSParticipants::thePartonPairs
protected

Definition at line 145 of file G4QGSParticipants.hh.

Referenced by GetNextPartonPair(), and PerformSoftCollisions().

◆ theProjectileSplitable

G4QGSMSplitableHadron* G4QGSParticipants::theProjectileSplitable
protected

◆ theQuarkExchange

G4QuarkExchange G4QGSParticipants::theQuarkExchange
protected

Definition at line 147 of file G4QGSParticipants.hh.

Referenced by PerformDiffractiveCollisions().

◆ theSingleDiffExcitation

G4SingleDiffractiveExcitation G4QGSParticipants::theSingleDiffExcitation
protected

Definition at line 148 of file G4QGSParticipants.hh.

Referenced by PerformDiffractiveCollisions().

◆ theTargets

std::vector<G4VSplitableHadron*> G4QGSParticipants::theTargets
protected

◆ ThresholdParameter

const G4double G4QGSParticipants::ThresholdParameter
protected

Definition at line 162 of file G4QGSParticipants.hh.

Referenced by SelectInteractions().


The documentation for this class was generated from the following files: