Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4QGSMFragmentation Class Reference

#include <G4QGSMFragmentation.hh>

Inheritance diagram for G4QGSMFragmentation:
Collaboration diagram for G4QGSMFragmentation:

Public Member Functions

 G4QGSMFragmentation ()
 
 ~G4QGSMFragmentation ()
 
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)
 
- Public Member Functions inherited from G4VLongitudinalStringDecay
 G4VLongitudinalStringDecay ()
 
virtual ~G4VLongitudinalStringDecay ()
 
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)
 
G4ExcitedStringCPExcited (const G4ExcitedString &string)
 
G4KineticTrackSplitup (G4FragmentingString *string, G4FragmentingString *&newString)
 
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 ()
 
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 ( )

Definition at line 49 of file G4QGSMFragmentation.cc.

49  :
50 arho(0.5), aphi(0.), an(-0.5), ala(-0.75), aksi(-1.), alft(0.5)
51 {
55 }
void SetStrangenessSuppression(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)
void SetDiquarkSuppression(G4double aValue)

Here is the call graph for this function:

G4QGSMFragmentation::~G4QGSMFragmentation ( )

Definition at line 57 of file G4QGSMFragmentation.cc.

58 {
59 }

Member Function Documentation

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

Implements G4VLongitudinalStringDecay.

Definition at line 63 of file G4QGSMFragmentation.cc.

64 {
65  #ifdef debug_QGSMfragmentation
66  G4cout<<G4endl<<"QGSM StringFragm: String Mass " <<theString.Get4Momentum().mag()<<" Pz "
67  <<theString.Get4Momentum().pz()
68  <<"------------------------------------"<<G4endl;
69  G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
70  <<theString.GetRightParton()->GetPDGcode()<<" "
71  <<theString.GetDirection()<< G4endl;
72  G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
73  G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl;
74  G4cout<<"Check for Fragmentation "<<G4endl;
75  #endif
76 
77  // Can no longer modify Parameters for Fragmentation.
78  PastInitPhase=true;
79 
80  // check if string has enough mass to fragment...
81 
82  G4KineticTrackVector * LeftVector=LightFragmentationTest(&theString);
83 
84  #ifdef debug_QGSMfragmentation
85  if ( LeftVector != 0 ) G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
86  #endif
87 
88  if ( LeftVector != 0 ) return LeftVector;
89 
90  #ifdef debug_QGSMfragmentation
91  G4cout<<"The string will be fragmented. "<<G4endl;
92  #endif
93 
94  LeftVector = new G4KineticTrackVector;
96 
97  // this should work but its only a semi deep copy. %GF G4ExcitedString theStringInCMS(theString);
98  G4ExcitedString *theStringInCMS=CPExcited(theString);
99  G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
100 
101  G4bool success=false, inner_sucess=true;
102  G4int attempt=0;
103  while ( !success && attempt++ < StringLoopInterrupt ) /* Loop checking, 07.08.2015, A.Ribon */
104  {
105  #ifdef debug_QGSMfragmentation
106  G4cout<<"Loop_toFrag "<<theStringInCMS->GetLeftParton()->GetPDGcode()<<" "
107  <<theStringInCMS->GetRightParton()->GetPDGcode()<<" "
108  <<theStringInCMS->GetDirection()<< G4endl;
109  #endif
110 
111  G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
112 
113  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
114  LeftVector->clear();
115  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
116  RightVector->clear();
117 
118  inner_sucess=true; // set false on failure..
119  const G4int maxNumberOfLoops = 1000;
120  G4int loopCounter = -1;
121  while (! StopFragmenting(currentString) && ++loopCounter < maxNumberOfLoops ) /* Loop checking, 07.08.2015, A.Ribon */
122  { // Split current string into hadron + new string
123 
124  #ifdef debug_QGSMfragmentation
125  G4cout<<"The string can fragment. "<<G4endl;;
126  #endif
127 
128  G4FragmentingString *newString=0; // used as output from SplitUp...
129  G4KineticTrack * Hadron=Splitup(currentString,newString);
130 
131  if ( Hadron != 0 ) // && IsFragmentable(newString))
132  {
133  #ifdef debug_QGSMfragmentation
134  G4cout<<"Hadron prod at fragm. "<<Hadron->GetDefinition()->GetParticleName()<<G4endl;
135  #endif
136 
137  if ( currentString->GetDecayDirection() > 0 ) LeftVector->push_back(Hadron);
138  else RightVector->push_back(Hadron);
139 
140  delete currentString;
141  currentString=newString;
142 
143  } else {
144 
145  #ifdef debug_QGSMfragmentation
146  G4cout<<"abandon ... start from the beginning ---------------"<<G4endl;
147  #endif
148 
149  // abandon ... start from the beginning
150  if (newString) delete newString;
151  inner_sucess=false;
152  break;
153  }
154  }
155  if ( loopCounter >= maxNumberOfLoops ) inner_sucess=false;
156 
157  // Split current string into 2 final Hadrons
158  #ifdef debug_QGSMfragmentation
159  G4cout<<"Split remaining string into 2 final hadrons."<<G4endl;
160  #endif
161 
162  if ( inner_sucess && SplitLast(currentString,LeftVector, RightVector) )
163  {
164  success=true;
165  }
166  delete currentString;
167  } // End of while loop
168 
169  delete theStringInCMS;
170 
171  if ( ! success )
172  {
173  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
174  LeftVector->clear();
175  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
176  delete RightVector;
177  return LeftVector;
178  }
179 
180  // Join Left- and RightVector into LeftVector in correct order.
181  while(!RightVector->empty()) /* Loop checking, 07.08.2015, A.Ribon */
182  {
183  LeftVector->push_back(RightVector->back());
184  RightVector->erase(RightVector->end()-1);
185  }
186  delete RightVector;
187 
188  CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
189 
190  G4LorentzRotation toObserverFrame(toCms.inverse());
191 
192  for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
193  {
194  G4KineticTrack* Hadron = LeftVector->operator[](C1);
195  G4LorentzVector Momentum = Hadron->Get4Momentum();
196  Momentum = toObserverFrame*Momentum;
197  Hadron->Set4Momentum(Momentum);
198  G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
199  Momentum = toObserverFrame*Coordinate;
200  Hadron->SetFormationTime(Momentum.e());
201  G4ThreeVector aPosition(Momentum.vect());
202  Hadron->SetPosition(theString.GetPosition()+aPosition);
203  }
204  return LeftVector;
205 }
G4int GetPDGcode() const
Definition: G4Parton.hh:127
const double C1
const G4ThreeVector & GetPosition() const
G4ExcitedString * CPExcited(const G4ExcitedString &string)
void SetFormationTime(G4double aFormationTime)
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:143
G4Parton * GetLeftParton(void) const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
G4GLOB_DLL std::ostream G4cout
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
G4LorentzVector Get4Momentum() const
G4double GetFormationTime() const
double mag() const
bool G4bool
Definition: G4Types.hh:79
G4LorentzRotation TransformToAlignedCms()
void SetPosition(const G4ThreeVector aPosition)
G4int GetDecayDirection() const
void Set4Momentum(const G4LorentzVector &a4Momentum)
G4Parton * GetRightParton(void) const
double pz() const
G4int GetDirection(void) const
#define G4endl
Definition: G4ios.hh:61
HepLorentzRotation inverse() const
const G4LorentzVector & Get4Momentum() const
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
const G4ParticleDefinition * GetDefinition() const
G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)

Here is the call graph for this function:


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