Geant4 10.7.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
 
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 SetVectorMesonProbability (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
 
G4double pspin_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 SetStringTensionParameter(1.*GeV/fermi);
59 SetStrangenessSuppression((1.0 - 0.12)/2.0);
61
62 // Check if charmed and bottom hadrons are enabled: if this is the case, then
63 // set the non-zero probabilities for c-cbar and b-bbar creation from the vacuum,
64 // else set them to 0.0. If these probabilities are/aren't zero then charmed or bottom
65 // hadrons can't/can be created during the string fragmentation of ordinary
66 // (i.e. not heavy) projectile hadron nuclear reactions.
67 if ( G4HadronicParameters::Instance()->EnableBCParticles() ) {
68 SetProbCCbar(0.005); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
69 SetProbBBbar(5.0e-5); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
70 } else {
71 SetProbCCbar(0.0);
72 SetProbBBbar(0.0);
73 }
74
75 SetMinMasses(); // For treating of small string decays
76}
static G4HadronicParameters * Instance()
virtual void SetMassCut(G4double aValue)
void SetDiquarkSuppression(G4double aValue)
void SetStrangenessSuppression(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)
void SetStringTensionParameter(G4double aValue)

◆ ~G4LundStringFragmentation()

G4LundStringFragmentation::~G4LundStringFragmentation ( )
virtual

Definition at line 1347 of file G4LundStringFragmentation.cc.

1348{}

Member Function Documentation

◆ FragmentString()

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

Implements G4VLongitudinalStringDecay.

Definition at line 80 of file G4LundStringFragmentation.cc.

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