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

#include <G4QGSMFragmentation.hh>

+ Inheritance diagram for G4QGSMFragmentation:

Public Member Functions

 G4QGSMFragmentation ()
 
 ~G4QGSMFragmentation ()
 
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)
 
- Public Member Functions inherited from G4VLongitudinalStringDecay
 G4VLongitudinalStringDecay (const G4String &name="StringDecay")
 
virtual ~G4VLongitudinalStringDecay ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &, G4Nucleus &) final
 
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)=0
 
void AddNewParticles ()
 
void EraseNewParticles ()
 
void SetMinMasses ()
 
void SetMinimalStringMass (const G4FragmentingString *const string)
 
void SetMinimalStringMass2 (const G4double aValue)
 
G4int SampleQuarkFlavor (void)
 
G4ThreeVector SampleQuarkPt (G4double ptMax=-1.)
 
void SetSigmaTransverseMomentum (G4double aQT)
 
void SetStrangenessSuppression (G4double aValue)
 
void SetDiquarkSuppression (G4double aValue)
 
void SetDiquarkBreakProbability (G4double aValue)
 
void SetSpinThreeHalfBarionProbability (G4double aValue)
 
void SetScalarMesonMixings (std::vector< G4double > aVector)
 
void SetVectorMesonMixings (std::vector< G4double > aVector)
 
void SetStringTensionParameter (G4double aValue)
 
void SetProbCCbar (G4double aValue)
 
void SetProbEta_c (G4double aValue)
 
void SetProbBBbar (G4double aValue)
 
void SetProbEta_b (G4double aValue)
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Public Attributes inherited from G4VLongitudinalStringDecay
G4double Mass_of_light_quark
 
G4double Mass_of_s_quark
 
G4double Mass_of_c_quark
 
G4double Mass_of_b_quark
 
G4double Mass_of_string_junction
 
G4double minMassQQbarStr [5][5]
 
G4double minMassQDiQStr [5][5][5]
 
G4double MinimalStringMass
 
G4double MinimalStringMass2
 
G4int Qcharge [5]
 
G4int Meson [5][5][7]
 
G4double MesonWeight [5][5][7]
 
G4int Baryon [5][5][5][4]
 
G4double BaryonWeight [5][5][5][4]
 
G4double Prob_QQbar [5]
 
G4int DecayQuark
 
G4int NewQuark
 
G4ParticleDefinitionFS_LeftHadron [350]
 
G4ParticleDefinitionFS_RightHadron [350]
 
G4double FS_Weight [350]
 
G4int NumberOf_FS
 
- Protected Types inherited from G4VLongitudinalStringDecay
typedef std::pair< G4ParticleDefinition *, G4ParticleDefinition * > pDefPair
 
typedef G4ParticleDefinition *(G4HadronBuilder::* Pcreate) (G4ParticleDefinition *, G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VLongitudinalStringDecay
virtual void SetMassCut (G4double aValue)
 
G4double GetMassCut ()
 
G4KineticTrackVectorProduceOneHadron (const G4ExcitedString *const theString)
 
G4double PossibleHadronMass (const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
 
G4ParticleDefinitionFindParticle (G4int Encoding)
 
virtual G4bool StopFragmenting (const G4FragmentingString *const string)=0
 
virtual G4bool IsItFragmentable (const G4FragmentingString *const string)=0
 
virtual G4bool SplitLast (G4FragmentingString *string, G4KineticTrackVector *LeftVector, G4KineticTrackVector *RightVector)=0
 
virtual void Sample4Momentum (G4LorentzVector *Mom, G4double Mass, G4LorentzVector *AntiMom, G4double AntiMass, G4double InitialMass)=0
 
G4ExcitedStringCopyExcited (const G4ExcitedString &string)
 
virtual G4KineticTrackSplitup (G4FragmentingString *string, G4FragmentingString *&newString)=0
 
virtual G4ParticleDefinitionQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)
 
virtual G4ParticleDefinitionDiQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)=0
 
pDefPair CreatePartonPair (G4int NeedParticle, G4bool AllowDiquarks=true)
 
virtual G4LorentzVectorSplitEandP (G4ParticleDefinition *pHadron, G4FragmentingString *string, G4FragmentingString *newString)=0
 
virtual G4double GetLightConeZ (G4double zmin, G4double zmax, G4int PartonEncoding, G4ParticleDefinition *pHadron, G4double Px, G4double Py)=0
 
void CalculateHadronTimePosition (G4double theInitialStringMass, G4KineticTrackVector *)
 
void ConstructParticle ()
 
G4ParticleDefinitionCreateHadron (G4int id1, G4int id2, G4bool theGivenSpin, G4int theSpin)
 
G4double GetDiquarkSuppress ()
 
G4double GetDiquarkBreakProb ()
 
G4double GetStrangeSuppress ()
 
G4int GetClusterLoopInterrupt ()
 
G4double GetProbCCbar ()
 
G4double GetProbEta_c ()
 
G4double GetProbBBbar ()
 
G4double GetProbEta_b ()
 
G4double GetStringTensionParameter ()
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4VLongitudinalStringDecay
G4double MassCut
 
G4double SigmaQT
 
G4double DiquarkSuppress
 
G4double DiquarkBreakProb
 
G4double StrangeSuppress
 
G4int StringLoopInterrupt
 
G4int ClusterLoopInterrupt
 
G4HadronBuilderhadronizer
 
std::vector< G4doublepspin_meson
 
G4double pspin_barion
 
std::vector< G4doublevectorMesonMix
 
std::vector< G4doublescalarMesonMix
 
G4double ProbCCbar
 
G4double ProbEta_c
 
G4double ProbBBbar
 
G4double ProbEta_b
 
G4double ProbCB
 
G4double MaxMass
 
G4bool PastInitPhase
 
G4double Kappa
 
std::vector< G4ParticleDefinition * > NewParticles
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 39 of file G4QGSMFragmentation.hh.

Constructor & Destructor Documentation

◆ G4QGSMFragmentation()

G4QGSMFragmentation::G4QGSMFragmentation ( )

Definition at line 49 of file G4QGSMFragmentation.cc.

50{
51 SigmaQT = 0.45 * GeV;
52
53 MassCut = 0.35*GeV;
54
55 SetStrangenessSuppression((1.0 - 0.16)/2.);
56
57 // Check if charmed and bottom hadrons are enabled: if this is the case, then
58 // set the non-zero probabilities for c-cbar and b-bbar creation from the vacuum,
59 // else set them to 0.0. If these probabilities are/aren't zero then charmed or bottom
60 // hadrons can't/can be created during the string fragmentation of ordinary
61 // (i.e. not heavy) projectile hadron nuclear reactions.
62 if ( G4HadronicParameters::Instance()->EnableBCParticles() ) {
63 SetProbCCbar(0.0002); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094; tuned by Uzhi Oct. 2022
64 SetProbBBbar(5.0e-5); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
65 } else {
66 SetProbCCbar(0.0);
67 SetProbBBbar(0.0);
68 }
69
72
74
75 arho = 0.5; // alpha_rho0
76 aphi = 0.0; // alpha_fi
77 aJPs =-2.2; // alpha_J/Psi
78 aUps =-8.0; // alpha_Y ??? O. Piskunova Yad. Phys. 56 (1993) 1094.
79
80 aksi =-1.0;
81 alft = 0.5; // 2 * alpha'_R *<Pt^2>
82
83 an = -0.5 ;
84 ala = -0.75; // an - arho/2 + aphi/2
85 alaC = an - arho/2.0 + aJPs/2.0;
86 alaB = an - arho/2.0 + aUps/2.0;
87 aXi = 0.0; // ??
88 aXiC = 0.0; // ??
89 aXiB = 0.0; // ??
90 aXiCC = 0.0; // ??
91 aXiCB = 0.0; // ??
92 aXiBB = 0.0; // ??
93
94 SetFFq2q();
95 SetFFq2qq();
96 SetFFqq2q();
97 SetFFqq2qq();
98 // d u s c b
99 G4int Index[5][5] = { { 0, 1, 2, 3, 4 }, // d
100 { 1, 5, 6, 7, 8 }, // u
101 { 2, 6, 9, 10, 11 }, // s
102 { 3, 7, 10, 12, 13 }, // c
103 { 4, 8, 11, 13, 14 } }; // b
104 for (G4int i = 0; i < 5; i++ ) {
105 for ( G4int j = 0; j < 5; j++ ) {
106 IndexDiQ[i][j] = Index[i][j];
107 }
108 };
109}
int G4int
Definition: G4Types.hh:85
static G4HadronicParameters * Instance()
void SetDiquarkSuppression(G4double aValue)
void SetStrangenessSuppression(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)

◆ ~G4QGSMFragmentation()

G4QGSMFragmentation::~G4QGSMFragmentation ( )

Definition at line 111 of file G4QGSMFragmentation.cc.

112{}

Member Function Documentation

◆ FragmentString()

G4KineticTrackVector * G4QGSMFragmentation::FragmentString ( const G4ExcitedString theString)
virtual

Implements G4VLongitudinalStringDecay.

Definition at line 116 of file G4QGSMFragmentation.cc.

117{
118
119 G4FragmentingString aString(theString);
120 SetMinimalStringMass(&aString);
121
122 #ifdef debug_QGSMfragmentation
123 G4cout<<G4endl<<"QGSM StringFragm: String Mass "
124 <<theString.Get4Momentum().mag()<<" Pz "
125 <<theString.Get4Momentum().pz()
126 <<"------------------------------------"<<G4endl;
127 G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
128 <<theString.GetRightParton()->GetPDGcode()<<" "
129 <<theString.GetDirection()<< G4endl;
130 G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
131 G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl;
132 G4cout<<"Check for Fragmentation "<<G4endl;
133 #endif
134
135 // Can no longer modify Parameters for Fragmentation.
136 PastInitPhase=true;
137
138 // Check if string has enough mass to fragment...
139 G4KineticTrackVector * LeftVector=NULL;
140
141 if ( !IsItFragmentable(&aString) ) {
142 LeftVector=ProduceOneHadron(&theString);
143
144 #ifdef debug_QGSMfragmentation
145 if ( LeftVector != 0 ) G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
146 #endif
147
148 if ( LeftVector == nullptr ) LeftVector = new G4KineticTrackVector;
149 return LeftVector;
150 }
151
152 #ifdef debug_QGSMfragmentation
153 G4cout<<"The string will be fragmented. "<<G4endl;
154 #endif
155
156 LeftVector = new G4KineticTrackVector;
158
159 G4ExcitedString *theStringInCMS=CopyExcited(theString);
160 G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
161
162 G4bool success=false, inner_sucess=true;
163 G4int attempt=0;
164 while ( !success && attempt++ < StringLoopInterrupt ) /* Loop checking, 07.08.2015, A.Ribon */
165 {
166 #ifdef debug_QGSMfragmentation
167 G4cout<<"Loop_toFrag "<<theStringInCMS->GetLeftParton()->GetPDGcode()<<" "
168 <<theStringInCMS->GetRightParton()->GetPDGcode()<<" "
169 <<theStringInCMS->GetDirection()<< G4endl;
170 #endif
171
172 G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
173
174 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
175 LeftVector->clear();
176 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
177 RightVector->clear();
178
179 inner_sucess=true; // set false on failure..
180 const G4int maxNumberOfLoops = 1000;
181 G4int loopCounter = -1;
182 while (! StopFragmenting(currentString) && ++loopCounter < maxNumberOfLoops ) /* Loop checking, 07.08.2015, A.Ribon */
183 { // Split current string into hadron + new string
184
185 #ifdef debug_QGSMfragmentation
186 G4cout<<"The string can fragment. "<<G4endl;;
187 #endif
188 G4FragmentingString *newString=0; // used as output from SplitUp...
189 G4KineticTrack * Hadron=Splitup(currentString,newString);
190
191 if ( Hadron != 0 )
192 {
193 #ifdef debug_QGSMfragmentation
194 G4cout<<"Hadron prod at fragm. "<<Hadron->GetDefinition()->GetParticleName()<<G4endl;
195 #endif
196 // To close the production of hadrons at fragmentation stage
197 if ( currentString->GetDecayDirection() > 0 )
198 LeftVector->push_back(Hadron);
199 else
200 RightVector->push_back(Hadron);
201
202 delete currentString;
203 currentString=newString;
204
205 } else {
206
207 #ifdef debug_QGSMfragmentation
208 G4cout<<"abandon ... start from the beginning ---------------"<<G4endl;
209 #endif
210
211 // Abandon ... start from the beginning
212 if (newString) delete newString;
213 inner_sucess=false;
214 break;
215 }
216 }
217 if ( loopCounter >= maxNumberOfLoops ) {
218 inner_sucess=false;
219 }
220
221 // Split current string into 2 final Hadrons
222 #ifdef debug_QGSMfragmentation
223 if( inner_sucess ) {
224 G4cout<<"Split remaining string into 2 final hadrons."<<G4endl;
225 } else {
226 G4cout<<" New attempt to fragment string"<<G4endl;
227 }
228 #endif
229 // To the close production of hadrons at last string decay
230 if ( inner_sucess &&
231 SplitLast(currentString,LeftVector, RightVector) )
232 {
233 success=true;
234 }
235 delete currentString;
236 }
237
238 delete theStringInCMS;
239
240 if ( ! success )
241 {
242 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
243 LeftVector->clear();
244 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
245 delete RightVector;
246 return LeftVector;
247 }
248
249 // Join Left- and RightVector into LeftVector in correct order.
250 while(!RightVector->empty()) /* Loop checking, 07.08.2015, A.Ribon */
251 {
252 LeftVector->push_back(RightVector->back());
253 RightVector->erase(RightVector->end()-1);
254 }
255 delete RightVector;
256
257 CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
258
259 G4LorentzRotation toObserverFrame(toCms.inverse());
260
261 for (size_t C1 = 0; C1 < LeftVector->size(); C1++)
262 {
263 G4KineticTrack* Hadron = LeftVector->operator[](C1);
265 Momentum = toObserverFrame*Momentum;
266 Hadron->Set4Momentum(Momentum);
267 G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
268 Momentum = toObserverFrame*Coordinate;
269 Hadron->SetFormationTime(Momentum.e());
270 G4ThreeVector aPosition(Momentum.vect());
271 Hadron->SetPosition(theString.GetPosition()+aPosition);
272 }
273 return LeftVector;
274}
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define C1
HepLorentzRotation inverse() const
const G4ThreeVector & GetPosition() const
G4int GetDirection(void) const
G4LorentzRotation TransformToAlignedCms()
G4Parton * GetRightParton(void) const
G4Parton * GetLeftParton(void) const
G4LorentzVector Get4Momentum() const
G4int GetDecayDirection() const
G4double GetFormationTime() const
void SetPosition(const G4ThreeVector aPosition)
void Set4Momentum(const G4LorentzVector &a4Momentum)
const G4ThreeVector & GetPosition() const
const G4ParticleDefinition * GetDefinition() const
void SetFormationTime(G4double aFormationTime)
const G4LorentzVector & Get4Momentum() const
const G4String & GetParticleName() const
G4int GetPDGcode() const
Definition: G4Parton.hh:127
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:143
G4KineticTrackVector * ProduceOneHadron(const G4ExcitedString *const theString)
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
void SetMinimalStringMass(const G4FragmentingString *const string)
G4ExcitedString * CopyExcited(const G4ExcitedString &string)

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