Geant4 10.7.0
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)
 
virtual void SetProjectileNucleus (G4V3DNucleus *aNucleus)
 
virtual G4V3DNucleusGetProjectileNucleus () const
 

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 nCutMax(7),ThresholdParameter(0.000*GeV),
52 QGSMThreshold(3*GeV),theNucleonRadius(1.5*fermi),alpha(-0.5),beta(2.5)
53{
54 // Parameters setting
55 SetCofNuclearDestruction(1.);
56 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
57 SetDofNuclearDestruction( 0.3 );
58 SetPt2ofNuclearDestruction( 0.075*GeV*GeV );
59 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
60 SetExcitationEnergyPerWoundedNucleon( 40.0*MeV );
61
62 sigmaPt=0.25*sqr(GeV);
63}
const G4double ThresholdParameter
G4QGSDiffractiveExcitation theDiffExcitaton
const G4double QGSMThreshold
const G4double theNucleonRadius
T sqr(const T &x)
Definition: templates.hh:128

◆ G4QGSParticipants() [2/2]

G4QGSParticipants::G4QGSParticipants ( const G4QGSParticipants right)

◆ ~G4QGSParticipants()

G4QGSParticipants::~G4QGSParticipants ( )
virtual

Definition at line 71 of file G4QGSParticipants.cc.

71{}

Member Function Documentation

◆ BuildInteractions()

void G4QGSParticipants::BuildInteractions ( const G4ReactionProduct thePrimary)

Definition at line 73 of file G4QGSParticipants.cc.

74{
75 theProjectile = thePrimary;
76
77 Regge = new G4Reggeons(theProjectile.GetDefinition());
78
80
81 NumberOfInvolvedNucleonsOfProjectile= 0;
82 G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
83 ProjectileResidualMassNumber = 0;
84 ProjectileResidualCharge = 0;
85 ProjectileResidualExcitationEnergy = 0.0;
86 ProjectileResidual4Momentum = tmp;
87
88 NumberOfInvolvedNucleonsOfTarget= 0;
89 TargetResidualMassNumber = theNucleus->GetMassNumber();
90 TargetResidualCharge = theNucleus->GetCharge();
91 TargetResidualExcitationEnergy = 0.0;
92
94 G4Nucleon* NuclearNucleon;
95 while ( ( NuclearNucleon = theNucleus->GetNextNucleon() ) ) {
96 tmp+=NuclearNucleon->Get4Momentum();
97 }
98 TargetResidual4Momentum = tmp;
99
100 if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
101 // Projectile is a hadron : meson or baryon
102 ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
103 ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
104 ProjectileResidualExcitationEnergy = 0.0;
105 ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
106 ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
107 }
108
109 #ifdef debugQGSParticipants
110 G4cout <<G4endl<< "G4QGSParticipants::BuildInteractions---------" << G4endl
111 << "thePrimary " << thePrimary.GetDefinition()->GetParticleName()<<" "
112 <<ProjectileResidual4Momentum<<G4endl;
113 G4cout << "Target A and Z " << theNucleus->GetMassNumber() <<" "<<theNucleus->GetCharge()<<" "
114 << TargetResidual4Momentum<<G4endl;
115 #endif
116
117 G4bool Success( true );
118
119 const G4int maxNumberOfLoops = 1000;
120 G4int loopCounter = 0;
121 do
122 {
123 const G4int maxNumberOfInternalLoops = 1000;
124 G4int internalLoopCounter = 0;
125 do
126 {
127 if(std::abs(theProjectile.GetDefinition()->GetPDGEncoding()) < 100.0)
128 {
129 SelectInteractions(theProjectile); // for lepton projectile
130 }
131 else
132 {
133 GetList( theProjectile ); // Get list of participating nucleons for hadron projectile
134 }
135
136 if ( theInteractions.size() == 0 ) return;
137
138 StoreInvolvedNucleon(); // Store participating nucleon
139
140 ReggeonCascade(); // Make reggeon cascading. Involve nucleons in the cascading.
141
142 Success = PutOnMassShell(); // Ascribe momenta to the involved and participating nucleons
143
144 if(!Success) PrepareInitialState( thePrimary );
145
146 } while( (!Success) && ++internalLoopCounter < maxNumberOfInternalLoops );
147
148 if ( internalLoopCounter >= maxNumberOfInternalLoops ) {
149 Success = false;
150 }
151
152 if ( Success ) {
153 #ifdef debugQGSParticipants
154 G4cout<<G4endl<<"PerformDiffractiveCollisions(); if they happend." <<G4endl;
155 #endif
156
158
159 #ifdef debugQGSParticipants
160 G4cout<<G4endl<<"SplitHadrons();" <<G4endl;
161 #endif
162
163 SplitHadrons();
164
166 #ifdef debugQGSParticipants
167 G4cout<<"Perform non-Diffractive Collisions if they happend. Determine Parton Momenta and so on." <<G4endl;
168 #endif
169 Success = DeterminePartonMomenta();
170 }
171
172 if(!Success) PrepareInitialState( thePrimary );
173 }
174 } while( (!Success) && ++loopCounter < maxNumberOfLoops );
175
176 if ( loopCounter >= maxNumberOfLoops ) {
177 Success = false;
178 #ifdef debugQGSParticipants
179 G4cout<<"NOT Successful ======" <<G4endl;
180 #endif
181 }
182
183 if ( Success ) {
184 CreateStrings(); // To create strings
185
186 GetResiduals(); // To calculate residual nucleus properties
187
188 #ifdef debugQGSParticipants
189 G4cout<<"All O.K. ======" <<G4endl;
190 #endif
191 }
192
193 // clean-up, if necessary
194 #ifdef debugQGSParticipants
195 G4cout<<"Clearing "<<G4endl;
196 G4cout<<"theInteractions.size() "<<theInteractions.size()<<G4endl;
197 #endif
198
199 if( Regge ) delete Regge;
200
201 std::for_each( theInteractions.begin(), theInteractions.end(), DeleteInteractionContent() );
202 theInteractions.clear();
203
204 // Erasing of target involved nucleons.
205 #ifdef debugQGSParticipants
206 G4cout<<"Erasing of involved target nucleons "<<NumberOfInvolvedNucleonsOfTarget<<G4endl;
207 #endif
208
209 if ( NumberOfInvolvedNucleonsOfTarget != 0 ) {
210 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; i++ ) {
211 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfTarget[i]->GetSplitableHadron();
212 if ( (aNucleon != 0 ) && (aNucleon->GetStatus() >= 1) ) delete aNucleon;
213 }
214 }
215
216 // Erasing of projectile involved nucleons.
217 if ( NumberOfInvolvedNucleonsOfProjectile != 0 ) {
218 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfProjectile; i++ ) {
219 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfProjectile[i]->GetSplitableHadron();
220 if ( aNucleon ) delete aNucleon;
221 }
222 }
223
224 #ifdef debugQGSParticipants
225 G4cout<<"Delition of target nucleons from soft interactions "<<theTargets.size()
226 <<G4endl<<G4endl;
227 #endif
228 std::for_each(theTargets.begin(), theTargets.end(), DeleteSplitableHadron());
229 theTargets.clear();
230
234 }
235}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void setVect(const Hep3Vector &)
G4VSplitableHadron * GetSplitableHadron() const
Definition: G4Nucleon.hh:95
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
G4double GetPDGCharge() const
const G4String & GetParticleName() const
std::vector< G4InteractionContent * > theInteractions
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
G4QGSMSplitableHadron * theProjectileSplitable
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 1497 of file G4QGSParticipants.cc.

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

1421{
1422 #ifdef debugQGSParticipants
1423 G4cout<<G4endl<<"PerformDiffractiveCollisions()......"<<G4endl
1424 <<"theInteractions.size() "<<theInteractions.size()<<G4endl;
1425 #endif
1426
1427 unsigned int i;
1428 for (i = 0; i < theInteractions.size(); i++)
1429 {
1430 G4InteractionContent* anIniteraction = theInteractions[i];
1431 #ifdef debugQGSParticipants
1432 G4cout<<"Interaction # and its status "
1433 <<i<<" "<<theInteractions[i]->GetStatus()<<G4endl;
1434 #endif
1435
1436 G4int InterStatus = theInteractions[i]->GetStatus();
1437 if ( (InterStatus == PrD) || (InterStatus == TrD) || (InterStatus == DD))
1438 { // Selection of diffractive interactions
1439 #ifdef debugQGSParticipants
1440 G4cout<<"Simulation of diffractive interaction. "<<InterStatus <<" PrD/TrD/DD/ND/Qech - 0,1,2,3,4"<<G4endl;
1441 #endif
1442
1443 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1444
1445 #ifdef debugQGSParticipants
1446 G4cout<<"The proj. before inter "
1449 G4cout<<"The targ. before inter " <<aTarget->Get4Momentum()<<" "
1450 <<aTarget->Get4Momentum().mag()<<G4endl;
1451 #endif
1452
1453 if ( InterStatus == PrD )
1455
1456 if ( InterStatus == TrD )
1458
1459 if ( InterStatus == DD )
1461
1462 #ifdef debugQGSParticipants
1463 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1465 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1466 <<aTarget->Get4Momentum().mag()<<G4endl;
1467 #endif
1468 }
1469
1470 if ( InterStatus == Qexc )
1471 { // Quark exchange process
1472 #ifdef debugQGSParticipants
1473 G4cout<<"Simulation of interaction with quark exchange."<<G4endl;
1474 #endif
1475 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1476
1477 #ifdef debugQGSParticipants
1478 G4cout<<"The proj. before inter " <<theProjectileSplitable->Get4Momentum()<<" "
1480 G4cout<<"The targ. before inter "<<aTarget->Get4Momentum()<<" "
1481 <<aTarget->Get4Momentum().mag()<<G4endl;
1482 #endif
1483
1485
1486 #ifdef debugQGSParticipants
1487 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1489 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1490 <<aTarget->Get4Momentum().mag()<<G4endl;
1491 #endif
1492 }
1493 }
1494}
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
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

Referenced by BuildInteractions().

◆ PerformSoftCollisions()

void G4QGSParticipants::PerformSoftCollisions ( )
protected

Definition at line 2297 of file G4QGSParticipants.cc.

2298{
2299 std::vector<G4InteractionContent*>::iterator i;
2300 G4LorentzVector str4Mom;
2301 i = theInteractions.begin();
2302 while ( i != theInteractions.end() ) /* Loop checking, 07.08.2015, A.Ribon */
2303 {
2304 G4InteractionContent* anIniteraction = *i;
2305 G4PartonPair * aPair = NULL;
2306 if (anIniteraction->GetNumberOfSoftCollisions())
2307 {
2308 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
2309 G4VSplitableHadron* pTarget = anIniteraction->GetTarget();
2310 for (G4int j = 0; j < anIniteraction->GetNumberOfSoftCollisions(); j++)
2311 {
2312 aPair = new G4PartonPair(pTarget->GetNextParton(), pProjectile->GetNextAntiParton(),
2314 #ifdef debugQGSParticipants
2315 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2316 << aPair->GetParton1()->Get4Momentum() << " "
2317 << aPair->GetParton1()->GetX() << " " << G4endl;
2318 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2319 << aPair->GetParton2()->Get4Momentum() << " "
2320 << aPair->GetParton2()->GetX() << " " << G4endl;
2321 #endif
2322 #ifdef debugQGSParticipants
2323 str4Mom += aPair->GetParton1()->Get4Momentum();
2324 str4Mom += aPair->GetParton2()->Get4Momentum();
2325 #endif
2326 thePartonPairs.push_back(aPair);
2327 aPair = new G4PartonPair(pProjectile->GetNextParton(), pTarget->GetNextAntiParton(),
2329 #ifdef debugQGSParticipants
2330 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2331 << aPair->GetParton1()->Get4Momentum() << " "
2332 << aPair->GetParton1()->GetX() << " " << G4endl;
2333 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2334 << aPair->GetParton2()->Get4Momentum() << " "
2335 << aPair->GetParton2()->GetX() << " " << G4endl;
2336 #endif
2337 #ifdef debugQGSParticipants
2338 str4Mom += aPair->GetParton1()->Get4Momentum();
2339 str4Mom += aPair->GetParton2()->Get4Momentum();
2340 #endif
2341 thePartonPairs.push_back(aPair);
2342 }
2343 delete *i;
2344 i=theInteractions.erase(i); // i now points to the next interaction
2345 } else {
2346 i++;
2347 }
2348 }
2349 #ifdef debugQGSParticipants
2350 G4cout << " string 4 mom " << str4Mom << G4endl;
2351 #endif
2352}
G4VSplitableHadron * GetProjectile() const
G4Parton * GetParton2(void)
Definition: G4PartonPair.hh:76
G4Parton * GetParton1(void)
Definition: G4PartonPair.hh:71
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 2007 of file G4QGSParticipants.cc.

2009{
2010 G4double Xmin=anXmin; G4int Nsea=totalSea; Xmin*=1.; Nsea++; // Must be erased
2011 G4double Oalfa = 1./(alpha + 1.);
2012 G4double Obeta = 1./(aBeta + (alpha + 1.)*nSea + 1.); // ?
2013
2014 G4double Ksi1, Ksi2, r1, r2, r12;
2015 const G4int maxNumberOfLoops = 1000;
2016 G4int loopCounter = 0;
2017 do
2018 {
2019 Ksi1 = G4UniformRand(); r1 = G4Pow::GetInstance()->powA(Ksi1,Oalfa);
2020 Ksi2 = G4UniformRand(); r2 = G4Pow::GetInstance()->powA(Ksi2,Obeta);
2021 r12=r1+r2;
2022 } while( ( r12 > 1.) &&
2023 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
2024 if ( loopCounter >= maxNumberOfLoops ) {
2025 return 0.5; // Just an acceptable value, without any physics consideration.
2026 }
2027
2028 G4double result = r1/r12;
2029 return result;
2030}
#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 2356 of file G4QGSParticipants.cc.

2357{
2358 // Check reaction threshold - goes to CheckThreshold
2359
2362
2363 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
2364 G4LorentzVector aTargetNMomentum(0.,0.,0.,938.);
2365
2366 if ((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
2367 {
2368 throw G4HadronicException(__FILE__, __LINE__,
2369 "G4GammaParticipants::SelectInteractions: primary nan energy.");
2370 }
2371 G4double S = (aPrimaryMomentum + aTargetNMomentum).mag2();
2372 G4double ThresholdMass = thePrimary.GetMass() + 938.;
2373 ModelMode = SOFT;
2374
2375 #ifdef debugQGSParticipants
2376 G4cout << "Gamma projectile - SelectInteractions " << G4endl;
2377 G4cout<<"Energy and Nucleus "<<thePrimary.GetTotalEnergy()<<" "<<theNucleus->GetMassNumber()<<G4endl;
2378 G4cout << "SqrtS ThresholdMass ModelMode " <<std::sqrt(S)<<" "<<ThresholdMass<<" "<<ModelMode<< G4endl;
2379 #endif
2380
2381 if (sqr(ThresholdMass + ThresholdParameter) > S)
2382 {
2384 }
2385 if (sqr(ThresholdMass + QGSMThreshold) > S) // thus only diffractive in cascade!
2386 {
2388 }
2389
2390 // first find the collisions HPW
2391 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
2392 theInteractions.clear();
2393 G4int totalCuts = 0;
2394
2396 G4int NucleonNo=0;
2397
2399 G4Nucleon * pNucleon = 0;
2400
2401 while( (pNucleon = theNucleus->GetNextNucleon()) ) /* Loop checking, 07.08.2015, A.Ribon */
2402 { if(NucleonNo == theCurrent) break; NucleonNo++; }
2403
2404 if ( pNucleon ) {
2405 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(*pNucleon);
2406
2407 pNucleon->Hit(aTarget);
2408
2409 if ( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) )
2410 {
2413
2414 aInteraction->SetTarget(aTarget);
2415 aInteraction->SetTargetNucleon(pNucleon);
2416 aTarget->SetCollisionCount(0);
2417 aTarget->SetStatus(1);
2418
2419 aInteraction->SetNumberOfDiffractiveCollisions(1);
2420 aInteraction->SetNumberOfSoftCollisions(0);
2421 aInteraction->SetStatus(1);
2422
2423 theInteractions.push_back(aInteraction);
2424 totalCuts += 1;
2425 }
2426 else
2427 {
2428 // nondiffractive soft interaction occurs
2429 aTarget->IncrementCollisionCount(1);
2430 aTarget->SetStatus(0);
2431 theTargets.push_back(aTarget);
2432
2435
2437 aInteraction->SetTarget(aTarget);
2438 aInteraction->SetTargetNucleon(pNucleon);
2439 aInteraction->SetNumberOfSoftCollisions(1);
2440 aInteraction->SetStatus(3);
2441 theInteractions.push_back(aInteraction);
2442 totalCuts += 1;
2443 }
2444 }
2446}
void SetNumberOfDiffractiveCollisions(int)
void SetTargetNucleon(G4Nucleon *aNucleon)
void SetTarget(G4VSplitableHadron *aTarget)
void SetStatus(G4int aValue)
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:89
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: