Geant4 9.6.0
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 ()
 
virtual ~G4VLongitudinalStringDecay ()
 
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)=0
 
G4int SampleQuarkFlavor (void)
 
G4ThreeVector SampleQuarkPt (G4double ptMax=-1.)
 
G4KineticTrackVectorDecayResonans (G4KineticTrackVector *aHadrons)
 
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)
 

Additional Inherited Members

- 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)
 
G4KineticTrackVectorLightFragmentationTest (const G4ExcitedString *const theString)
 
G4double FragmentationMass (const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
 
G4ParticleDefinitionFindParticle (G4int Encoding)
 
virtual void Sample4Momentum (G4LorentzVector *Mom, G4double Mass, G4LorentzVector *AntiMom, G4double AntiMass, G4double InitialMass)=0
 
virtual G4bool StopFragmenting (const G4FragmentingString *const string)=0
 
virtual G4bool IsFragmentable (const G4FragmentingString *const string)=0
 
virtual G4bool SplitLast (G4FragmentingString *string, G4KineticTrackVector *LeftVector, G4KineticTrackVector *RightVector)=0
 
G4ExcitedStringCPExcited (const G4ExcitedString &string)
 
G4KineticTrackSplitup (G4FragmentingString *string, G4FragmentingString *&newString)
 
G4ParticleDefinitionQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)
 
G4ParticleDefinitionDiQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)
 
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 ()
 
G4double GetClusterMass ()
 
G4int GetClusterLoopInterrupt ()
 
G4double GetStringTensionParameter ()
 
- Protected Attributes inherited from G4VLongitudinalStringDecay
G4double MassCut
 
G4double ClusterMass
 
G4double SigmaQT
 
G4double DiquarkSuppress
 
G4double DiquarkBreakProb
 
G4double SmoothParam
 
G4double StrangeSuppress
 
G4int StringLoopInterrupt
 
G4int ClusterLoopInterrupt
 
G4HadronBuilderhadronizer
 
G4double pspin_meson
 
G4double pspin_barion
 
std::vector< G4doublevectorMesonMix
 
std::vector< G4doublescalarMesonMix
 
G4bool PastInitPhase
 
G4double Kappa
 

Detailed Description

Definition at line 40 of file G4QGSMFragmentation.hh.

Constructor & Destructor Documentation

◆ G4QGSMFragmentation()

G4QGSMFragmentation::G4QGSMFragmentation ( )

Definition at line 45 of file G4QGSMFragmentation.cc.

45 :
46arho(0.5), aphi(0.), an(-0.5), ala(-0.75), aksi(-1.), alft(0.5)
47 {
48 }

◆ ~G4QGSMFragmentation()

G4QGSMFragmentation::~G4QGSMFragmentation ( )

Definition at line 50 of file G4QGSMFragmentation.cc.

51 {
52 }

Member Function Documentation

◆ FragmentString()

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

Implements G4VLongitudinalStringDecay.

Definition at line 56 of file G4QGSMFragmentation.cc.

57{
58// Can no longer modify Parameters for Fragmentation.
59 PastInitPhase=true;
60
61// check if string has enough mass to fragment...
62 G4KineticTrackVector * LeftVector=LightFragmentationTest(&theString);
63 if ( LeftVector != 0 ) return LeftVector;
64
65 LeftVector = new G4KineticTrackVector;
67
68// this should work but its only a semi deep copy. %GF G4ExcitedString theStringInCMS(theString);
69 G4ExcitedString *theStringInCMS=CPExcited(theString);
70 G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
71
72 G4bool success=false, inner_sucess=true;
73 G4int attempt=0;
74 while ( !success && attempt++ < StringLoopInterrupt )
75 {
76 G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
77
78 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
79 LeftVector->clear();
80 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
81 RightVector->clear();
82
83 inner_sucess=true; // set false on failure..
84 while (! StopFragmenting(currentString) )
85 { // Split current string into hadron + new string
86 G4FragmentingString *newString=0; // used as output from SplitUp...
87 G4KineticTrack * Hadron=Splitup(currentString,newString);
88 if ( Hadron != 0 && IsFragmentable(newString))
89 {
90 if ( currentString->GetDecayDirection() > 0 )
91 LeftVector->push_back(Hadron);
92 else
93 RightVector->push_back(Hadron);
94 delete currentString;
95 currentString=newString;
96 } else {
97 // abandon ... start from the beginning
98 if (newString) delete newString; // Uzhi restore 20.06.08
99 if (Hadron) delete Hadron;
100 inner_sucess=false;
101 break;
102 }
103 }
104 // Split current string into 2 final Hadrons
105 if ( inner_sucess &&
106 SplitLast(currentString,LeftVector, RightVector) )
107 {
108 success=true;
109 }
110 delete currentString;
111 }
112
113 delete theStringInCMS;
114
115 if ( ! success )
116 {
117 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
118 LeftVector->clear();
119 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
120 delete RightVector;
121 return LeftVector;
122 }
123
124 // Join Left- and RightVector into LeftVector in correct order.
125 while(!RightVector->empty())
126 {
127 LeftVector->push_back(RightVector->back());
128 RightVector->erase(RightVector->end()-1);
129 }
130 delete RightVector;
131
132 CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
133
134 G4LorentzRotation toObserverFrame(toCms.inverse());
135
136 for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
137 {
138 G4KineticTrack* Hadron = LeftVector->operator[](C1);
139 G4LorentzVector Momentum = Hadron->Get4Momentum();
140 Momentum = toObserverFrame*Momentum;
141 Hadron->Set4Momentum(Momentum);
142 G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
143 Momentum = toObserverFrame*Coordinate;
144 Hadron->SetFormationTime(Momentum.e());
145 G4ThreeVector aPosition(Momentum.vect());
146 Hadron->SetPosition(theString.GetPosition()+aPosition);
147 }
148 return LeftVector;
149
150
151
152}
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define C1
HepLorentzRotation inverse() const
const G4ThreeVector & GetPosition() const
G4LorentzRotation TransformToAlignedCms()
G4LorentzVector Get4Momentum() const
G4int GetDecayDirection() const
G4double GetFormationTime() const
void SetPosition(const G4ThreeVector aPosition)
void Set4Momentum(const G4LorentzVector &a4Momentum)
const G4ThreeVector & GetPosition() const
void SetFormationTime(G4double aFormationTime)
const G4LorentzVector & Get4Momentum() const
G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
G4ExcitedString * CPExcited(const G4ExcitedString &string)
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)

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