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

#include <G4ExcitedStringDecay.hh>

Inheritance diagram for G4ExcitedStringDecay:
Collaboration diagram for G4ExcitedStringDecay:

Public Member Functions

 G4ExcitedStringDecay ()
 
 G4ExcitedStringDecay (G4VLongitudinalStringDecay *aStringDecay)
 
virtual ~G4ExcitedStringDecay ()
 
virtual G4KineticTrackVectorFragmentStrings (const G4ExcitedStringVector *theStrings)
 
- Public Member Functions inherited from G4VStringFragmentation
 G4VStringFragmentation ()
 
virtual ~G4VStringFragmentation ()
 

Detailed Description

Definition at line 38 of file G4ExcitedStringDecay.hh.

Constructor & Destructor Documentation

G4ExcitedStringDecay::G4ExcitedStringDecay ( )

Definition at line 35 of file G4ExcitedStringDecay.cc.

35  : G4VStringFragmentation(),theStringDecay(0)
36 {}
G4ExcitedStringDecay::G4ExcitedStringDecay ( G4VLongitudinalStringDecay aStringDecay)

Definition at line 38 of file G4ExcitedStringDecay.cc.

40  theStringDecay(aStringDecay)
41 {}
G4ExcitedStringDecay::~G4ExcitedStringDecay ( )
virtual

Definition at line 50 of file G4ExcitedStringDecay.cc.

51 {
52 }

Member Function Documentation

G4KineticTrackVector * G4ExcitedStringDecay::FragmentStrings ( const G4ExcitedStringVector theStrings)
virtual

Implements G4VStringFragmentation.

Definition at line 76 of file G4ExcitedStringDecay.cc.

77 {
78  G4LorentzVector KTsum(0.,0.,0.,0.);
79 
80  #ifdef debug_G4ExcitedStringDecay
81  G4cout<<G4endl;
82  G4cout<<"--------------------------- G4ExcitedStringDecay ----------------------"<<G4endl;
83  G4cout<<"Hadronization of Excited Strings: theStrings->size() "<<theStrings->size()<<G4endl;
84  #endif
85 
86  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
87  {
88  if ( theStrings->operator[](astring)->IsExcited() ) {
89  KTsum+= theStrings->operator[](astring)->Get4Momentum();
90  } else {
91  KTsum+=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
92  }
93  }
94 
95  G4LorentzRotation toCms( -1 * KTsum.boostVector() );
96  G4LorentzRotation toLab(toCms.inverse());
97  G4LorentzVector Ptmp;
98  KTsum=G4LorentzVector(0.,0.,0.,0.);
99  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
100  {
101  if ( theStrings->operator[](astring)->IsExcited() ) {
102  Ptmp=toCms * theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
103  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
104 
105  Ptmp=toCms * theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
106  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
107 
108  KTsum+= theStrings->operator[](astring)->Get4Momentum();
109  } else {
110  Ptmp=toCms * theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
111  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
112  KTsum+= theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
113  }
114  }
115 
116  G4KineticTrackVector * theResult = new G4KineticTrackVector;
117  G4int attempts(0);
118  G4bool success=false;
119  G4bool NeedEnergyCorrector=false;
120  do {
121  #ifdef debug_G4ExcitedStringDecay
122  G4cout<<"New try No "<<attempts<<" to hadronize strings"<<G4endl;
123  #endif
124 
125  std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack());
126  theResult->clear();
127 
128  attempts++;
129 
130  G4LorentzVector KTsecondaries(0.,0.,0.,0.);
131  NeedEnergyCorrector=false;
132 
133  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
134  {
135  #ifdef debug_G4ExcitedStringDecay
136  G4cout<<"String No "<<astring+1<<" Excited? "<<theStrings->operator[](astring)->IsExcited()<<G4endl;
137 
138  G4cout<<"String No "<<astring+1<<" 4Momentum "<<theStrings->operator[](astring)->Get4Momentum()
139  <<" "<<theStrings->operator[](astring)->Get4Momentum().mag()<<G4endl;
140  #endif
141 
142  G4KineticTrackVector * generatedKineticTracks = NULL;
143  if ( theStrings->operator[](astring)->IsExcited() )
144  {
145  #ifdef debug_G4ExcitedStringDecay
146  G4cout<<"Fragment String with partons: "
147  <<theStrings->operator[](astring)->GetLeftParton()->GetPDGcode() <<" "
148  <<theStrings->operator[](astring)->GetRightParton()->GetPDGcode()<<" "
149  <<"Direction "<<theStrings->operator[](astring)->GetDirection()<<G4endl;
150  #endif
151  generatedKineticTracks=FragmentString(*theStrings->operator[](astring));
152  #ifdef debug_G4ExcitedStringDecay
153  G4cout<<"(G4ExcitedStringDecay) Number of produced hadrons = "<<generatedKineticTracks->size()<<G4endl;
154  #endif
155  } else {
156  #ifdef debug_G4ExcitedStringDecay
157  G4cout<<" GetTrack from the String"<<G4endl;
158  #endif
159  G4LorentzVector Mom=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
160  G4KineticTrack * aTrack= new G4KineticTrack(
161  theStrings->operator[](astring)->GetKineticTrack()->GetDefinition(),
162  theStrings->operator[](astring)->GetKineticTrack()->GetFormationTime(),
163  G4ThreeVector(0), Mom);
164 
165  aTrack->SetPosition(theStrings->operator[](astring)->GetKineticTrack()->GetPosition());
166 
167  #ifdef debug_G4ExcitedStringDecay
168  G4cout<<" A particle stored in the track is "<<aTrack->GetDefinition()->GetParticleName()<<G4endl;
169  #endif
170 
171  generatedKineticTracks = new G4KineticTrackVector;
172  generatedKineticTracks->push_back(aTrack);
173  }
174 
175  //if (generatedKineticTracks == NULL)
176  if (generatedKineticTracks->size() == 0)
177  {
178  //G4cerr << "G4VPartonStringModel:No KineticTracks produced" << G4endl;
179  //continue;
180  success=false; NeedEnergyCorrector=false; break;
181  }
182 
183  G4LorentzVector KTsum1(0.,0.,0.,0.);
184  for ( unsigned int aTrack=0; aTrack<generatedKineticTracks->size();aTrack++)
185  {
186  #ifdef debug_G4ExcitedStringDecay
187  G4cout<<"Prod part No. "<<aTrack+1<<" "
188  <<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" "
189  <<(*generatedKineticTracks)[aTrack]->Get4Momentum()<<G4endl;
190  #endif
191  theResult->push_back(generatedKineticTracks->operator[](aTrack));
192  KTsum1+= (*generatedKineticTracks)[aTrack]->Get4Momentum();
193  }
194  KTsecondaries+=KTsum1;
195 
196  #ifdef debug_G4ExcitedStringDecay
197  G4cout << "String secondaries(" <<generatedKineticTracks->size()<< ")"<<G4endl
198  <<"Init string momentum: "<< theStrings->operator[](astring)->Get4Momentum()<<G4endl
199  <<"Final hadrons momentum: "<< KTsum1 << G4endl;
200  #endif
201 
202  if ( KTsum1.e() > 0 && std::abs((KTsum1.e()-theStrings->operator[](astring)->Get4Momentum().e()) / KTsum1.e()) > perMillion )
203  {
204  NeedEnergyCorrector=true;
205  }
206 
207  #ifdef debug_G4ExcitedStringDecay
208  G4cout<<"NeedEnergyCorrection yes/no "<<NeedEnergyCorrector<<G4endl;
209  #endif
210 
211  // clean up
212  delete generatedKineticTracks;
213  success=true;
214  }
215 
216  //success=true;
217 
218  if ( NeedEnergyCorrector ) success=EnergyAndMomentumCorrector(theResult, KTsum);
219  } while(!success && (attempts < 100)); /* Loop checking, 07.08.2015, A.Ribon */
220 
221  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
222  {
223  Ptmp=(*theResult)[aTrack]->Get4Momentum();
224  Ptmp.transform( toLab);
225  (*theResult)[aTrack]->Set4Momentum(Ptmp);
226  }
227 
228  #ifdef debug_G4ExcitedStringDecay
229  G4cout<<"End of the strings fragmentation (G4ExcitedStringDecay)"<<G4endl;
230 
231  G4LorentzVector KTsum1(0.,0.,0.,0.);
232 
233  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
234  {
235  G4cout << " corrected tracks .. " << (*theResult)[aTrack]->GetDefinition()->GetParticleName()
236  <<" " << (*theResult)[aTrack]->Get4Momentum()
237  <<" " << (*theResult)[aTrack]->Get4Momentum().mag()<< G4endl;
238  KTsum1+= (*theResult)[aTrack]->Get4Momentum();
239  }
240 
241  G4cout << "Needcorrector/success " << NeedEnergyCorrector << "/" << success
242  << ", Corrected total 4 momentum " << KTsum1 << G4endl;
243  if ( ! success ) G4cout << "failed to correct E/p" << G4endl;
244 
245  G4cout<<"End of the Hadronization (G4ExcitedStringDecay)"<<G4endl;
246  #endif
247 
248  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
249  {
250  if ( theStrings->operator[](astring)->IsExcited() ) {
251  Ptmp=theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
252  Ptmp.transform( toLab);
253  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
254 
255  Ptmp=theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
256  Ptmp.transform( toLab);
257  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
258  } else {
259  Ptmp=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
260  Ptmp.transform( toLab);
261  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
262  }
263  }
264 
265  return theResult;
266 }
static constexpr double perMillion
Definition: G4SIunits.hh:334
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetPosition(const G4ThreeVector aPosition)
#define G4endl
Definition: G4ios.hh:61
const G4ParticleDefinition * GetDefinition() const
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:


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