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

#include <G4LundStringFragmentation.hh>

Inheritance diagram for G4LundStringFragmentation:

Public Member Functions

 G4LundStringFragmentation ()
virtual ~G4LundStringFragmentation ()
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)
Public Member Functions inherited from G4VLongitudinalStringDecay
 G4VLongitudinalStringDecay (const G4String &name="StringDecay")
virtual ~G4VLongitudinalStringDecay ()
G4HadFinalStateApplyYourself (const G4HadProjectile &, G4Nucleus &) final
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 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)
G4ExcitedStringCopyExcited (const G4ExcitedString &string)
virtual G4ParticleDefinitionQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)
pDefPair CreatePartonPair (G4int NeedParticle, G4bool AllowDiquarks=true)
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 41 of file G4LundStringFragmentation.hh.

Constructor & Destructor Documentation

◆ G4LundStringFragmentation()

G4LundStringFragmentation::G4LundStringFragmentation ( )

Definition at line 49 of file G4LundStringFragmentation.cc.

50 : G4VLongitudinalStringDecay("LundStringFragmentation")
51{
52 SetMassCut(210.*MeV); // Mpi + Delta
53 // For ProduceOneHadron it is required
54 // that no one pi-meson can be produced.
55 SigmaQT = 0.435 * GeV;
56 Tmt = 190.0 * MeV;
57
58 SetStringTensionParameter(1.*GeV/fermi);
60
61 SetStrangenessSuppression((1.0 - 0.12)/2.0);
63
64 // Check if charmed and bottom hadrons are enabled: if this is the case, then
65 // set the non-zero probabilities for c-cbar and b-bbar creation from the vacuum,
66 // else set them to 0.0. If these probabilities are/aren't zero then charmed or bottom
67 // hadrons can't/can be created during the string fragmentation of ordinary
68 // (i.e. not heavy) projectile hadron nuclear reactions.
69 if ( G4HadronicParameters::Instance()->EnableBCParticles() ) {
70 SetProbCCbar(0.0002); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094; tuned by Uzhi Oct. 2022
71 SetProbBBbar(5.0e-5); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
72 } else {
73 SetProbCCbar(0.0);
74 SetProbBBbar(0.0);
75 }
76
77 SetMinMasses(); // For treating of small string decays
78}
static G4HadronicParameters * Instance()
virtual void SetMassCut(G4double aValue)
void SetStrangenessSuppression(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)
void SetStringTensionParameter(G4double aValue)
G4VLongitudinalStringDecay(const G4String &name="StringDecay")

◆ ~G4LundStringFragmentation()

G4LundStringFragmentation::~G4LundStringFragmentation ( )
virtual

Definition at line 1360 of file G4LundStringFragmentation.cc.

1361{}

Member Function Documentation

◆ FragmentString()

G4KineticTrackVector * G4LundStringFragmentation::FragmentString ( const G4ExcitedString & theString)
virtual

Implements G4VLongitudinalStringDecay.

Definition at line 82 of file G4LundStringFragmentation.cc.

83{
84 // Can no longer modify Parameters for Fragmentation.
85
86 PastInitPhase=true;
87
88 G4FragmentingString aString(theString);
89 SetMinimalStringMass(&aString);
90
91 #ifdef debug_LUNDfragmentation
92 G4cout<<G4endl<<"LUND StringFragmentation ------------------------------------"<<G4endl;
93 G4cout<<G4endl<<"LUND StringFragm: String Mass "
94 <<theString.Get4Momentum().mag()<<G4endl
95 <<"4Mom "<<theString.Get4Momentum()<<G4endl
96 <<"------------------------------------"<<G4endl;
97 G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
98 <<theString.GetRightParton()->GetPDGcode()<<" "
99 <<theString.GetDirection()<< G4endl;
100 G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
101 G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl<<G4endl;
102 G4cout<<"Check for Fragmentation "<<G4endl;
103 #endif
104
105 G4KineticTrackVector * LeftVector(0);
106
107 if (!aString.IsAFourQuarkString() && !IsItFragmentable(&aString))
108 {
109 #ifdef debug_LUNDfragmentation
110 G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
111 #endif
112 // SetMassCut(210.*MeV); // For ProduceOneHadron it is required
113 // that no one pi-meson can be produced.
114
115 G4double Mcut = GetMassCut();
116 SetMassCut(10000.*MeV);
117 LeftVector=ProduceOneHadron(&theString);
118 SetMassCut(Mcut);
119
120 if ( LeftVector )
121 {
122 if ( LeftVector->size() > 0)
123 {
124 LeftVector->operator[](0)->SetFormationTime(theString.GetTimeOfCreation());
125 LeftVector->operator[](0)->SetPosition(theString.GetPosition());
126 }
127 if (LeftVector->size() > 1)
128 {
129 // 2 hadrons created from qq-qqbar are stored
130 LeftVector->operator[](1)->SetFormationTime(theString.GetTimeOfCreation());
131 LeftVector->operator[](1)->SetPosition(theString.GetPosition());
132 }
133 }
134 return LeftVector;
135 }
136
137 #ifdef debug_LUNDfragmentation
138 G4cout<<"The string will be fragmented. "<<G4endl;
139 #endif
140
141 // The string can fragment. At least two particles can be produced.
142 LeftVector =new G4KineticTrackVector;
143 G4KineticTrackVector * RightVector=new G4KineticTrackVector;
144
145 G4bool success = Loop_toFragmentString(theString, LeftVector, RightVector);
146
147 if ( ! success )
148 {
149 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
150 LeftVector->clear();
151 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
152 delete RightVector;
153 return LeftVector;
154 }
155
156 // Join Left- and RightVector into LeftVector in correct order.
157 while (!RightVector->empty())
158 {
159 LeftVector->push_back(RightVector->back());
160 RightVector->erase(RightVector->end()-1);
161 }
162 delete RightVector;
163
164 return LeftVector;
165}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetTimeOfCreation() const
const G4ThreeVector & GetPosition() const
G4int GetDirection(void) const
G4Parton * GetRightParton(void) const
G4Parton * GetLeftParton(void) const
G4LorentzVector Get4Momentum() const
G4int GetPDGcode() const
Definition G4Parton.hh:127
const G4LorentzVector & Get4Momentum() const
Definition G4Parton.hh:143
G4KineticTrackVector * ProduceOneHadron(const G4ExcitedString *const theString)
void SetMinimalStringMass(const G4FragmentingString *const string)

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