Geant4  10.01.p01
G4ExcitedStringDecay.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // Historic fragment from M.Komogorov; clean-up still necessary @@@
27 
28 #include "G4ExcitedStringDecay.hh"
29 #include "G4SystemOfUnits.hh"
30 #include "G4KineticTrack.hh"
31 
32 //#define debug_G4ExcitedStringDecay
33 //#define debug_G4ExcitedStringCorr
34 
36 {}
37 
40  theStringDecay(aStringDecay)
41 {}
42 
45  theStringDecay(0)
46 {
47  throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::copy ctor not accessible");
48 }
49 
51 {
52 }
53 
55 {
56  throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::operator= meant to not be accessable");
57  return *this;
58 }
59 
61 {
62  return 0;
63 }
64 
66 {
67  return 1;
68 }
69 
71 {
73  return theStringDecay->FragmentString(theString);
74 }
75 
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 {KTsum+=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();}
91  }
92 
93  G4LorentzRotation toCms( -1 * KTsum.boostVector() ); // Uzhi 22 June 2014
94  G4LorentzRotation toLab(toCms.inverse()); // Uzhi 22 June 2014
95  G4LorentzVector Ptmp;
96  KTsum=G4LorentzVector(0.,0.,0.,0.);
97  for ( unsigned int astring=0; astring < theStrings->size(); astring++) // Uzhi 22 June 2014
98  {
99  if ( theStrings->operator[](astring)->IsExcited() )
100  {
101  Ptmp=toCms * theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
102  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
103 
104  Ptmp=toCms * theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
105  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
106 
107  KTsum+= theStrings->operator[](astring)->Get4Momentum();
108  }
109  else
110  {
111  Ptmp=toCms * theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
112  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
113  KTsum+= theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
114  }
115  } // Uzhi 22 June 2014
116 
117  G4KineticTrackVector * theResult = new G4KineticTrackVector;
118  G4int attempts(0);
119  G4bool success=false;
120  G4bool NeedEnergyCorrector=false;
121  do {
122 #ifdef debug_G4ExcitedStringDecay
123  G4cout<<"New try No "<<attempts<<" to hadronize strings"<<G4endl;
124 #endif
125 
126  std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack());
127  theResult->clear();
128 
129  attempts++;
130 
131  G4LorentzVector KTsecondaries(0.,0.,0.,0.);
132  NeedEnergyCorrector=false;
133 
134  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
135 // for ( unsigned int astring=0; astring < 1; astring++)
136  {
137 #ifdef debug_G4ExcitedStringDecay
138  G4cout<<"String No "<<astring+1<<" Excited? "<<theStrings->operator[](astring)->IsExcited()<<G4endl;
139 
140  G4cout<<"String No "<<astring+1<<" 4Momentum "<<theStrings->operator[](astring)->Get4Momentum()
141  <<" "<<theStrings->operator[](astring)->Get4Momentum().mag()<<G4endl;
142 #endif
143 
144  G4KineticTrackVector * generatedKineticTracks = NULL;
145  if ( theStrings->operator[](astring)->IsExcited() )
146  {
147 #ifdef debug_G4ExcitedStringDecay
148  G4cout<<"Fragment String with partons: "
149  <<theStrings->operator[](astring)->GetLeftParton()->GetPDGcode() <<" "
150  <<theStrings->operator[](astring)->GetRightParton()->GetPDGcode()<<" "
151  <<"Direction "<<theStrings->operator[](astring)->GetDirection()<<G4endl;
152 #endif
153  generatedKineticTracks=FragmentString(*theStrings->operator[](astring));
154 #ifdef debug_G4ExcitedStringDecay
155  G4cout<<"(G4ExcitedStringDecay) Number of produced hadrons = "<<generatedKineticTracks->size()<<G4endl;
156 #endif
157  } else {
158 #ifdef debug_G4ExcitedStringDecay
159  G4cout<<" GetTrack from the String"<<G4endl;
160 #endif
161  G4LorentzVector Mom=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
162  G4KineticTrack * aTrack= new G4KineticTrack(
163  theStrings->operator[](astring)->GetKineticTrack()->GetDefinition(),
164  theStrings->operator[](astring)->GetKineticTrack()->GetFormationTime(),
165  G4ThreeVector(0), Mom);
166 
167  aTrack->SetPosition(theStrings->operator[](astring)->GetKineticTrack()->GetPosition());
168 
169 #ifdef debug_G4ExcitedStringDecay
170  G4cout<<" A particle stored in the track is "<<aTrack->GetDefinition()->GetParticleName()<<G4endl;
171 #endif
172 
173  generatedKineticTracks = new G4KineticTrackVector;
174  generatedKineticTracks->push_back(aTrack);
175  }
176 
177  if (generatedKineticTracks == NULL)
178  {
179  G4cerr << "G4VPartonStringModel:No KineticTracks produced" << G4endl;
180  continue;
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  }
214 
215  success=true;
216 //NeedEnergyCorrector=false; // Vova
217  if ( NeedEnergyCorrector ) success=EnergyAndMomentumCorrector(theResult, KTsum);
218  } while(!success && (attempts < 10)); // It was 100 !!! Uzhi
219 
220  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++) // Uzhi 22 June 2014
221  {
222  Ptmp=(*theResult)[aTrack]->Get4Momentum();
223  Ptmp.transform( toLab);
224  (*theResult)[aTrack]->Set4Momentum(Ptmp);
225  } // Uzhi 22 June 2014
226 
227 #ifdef debug_G4ExcitedStringDecay
228  G4cout<<"End of the strings fragmentation (G4ExcitedStringDecay)"<<G4endl;
229 
230  G4LorentzVector KTsum1(0.,0.,0.,0.);
231 
232  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
233  {
234  G4cout << " corrected tracks .. " << (*theResult)[aTrack]->GetDefinition()->GetParticleName()
235  <<" " << (*theResult)[aTrack]->Get4Momentum()
236  <<" " << (*theResult)[aTrack]->Get4Momentum().mag()<< G4endl;
237  KTsum1+= (*theResult)[aTrack]->Get4Momentum();
238  }
239 
240  G4cout << "Needcorrector/success " << NeedEnergyCorrector << "/" << success << ", Corrected total 4 momentum " << KTsum1 << G4endl;
241  if ( ! success ) G4cout << "failed to correct E/p" << G4endl;
242 
243  G4cout<<"End of the Hadronization (G4ExcitedStringDecay)"<<G4endl;
244 #endif
245 
246  for ( unsigned int astring=0; astring < theStrings->size(); astring++) // Uzhi 24 Oct. 2014
247  {
248  if ( theStrings->operator[](astring)->IsExcited() )
249  {
250  Ptmp=theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
251  Ptmp.transform( toLab);
252  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
253 
254  Ptmp=theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
255  Ptmp.transform( toLab);
256  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
257  }
258  else
259  {
260  Ptmp=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
261  Ptmp.transform( toLab);
262  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
263  }
264  } // Uzhi 24 Oct. 2014
265 
266  return theResult;
267 }
268 
270  (G4KineticTrackVector* Output, G4LorentzVector& TotalCollisionMom)
271  {
272  const int nAttemptScale = 500;
273  const double ErrLimit = 1.E-5;
274  if (Output->empty()) return TRUE;
275  G4LorentzVector SumMom;
276  G4double SumMass = 0;
277  G4double TotalCollisionMass = TotalCollisionMom.m();
278 
279 #ifdef debug_G4ExcitedStringCorr
280  G4cout<<G4endl<<"EnergyAndMomentumCorrector. Number of particles: "<<Output->size()<<G4endl;
281 #endif
282  // Calculate sum hadron 4-momenta and summing hadron mass
283  unsigned int cHadron;
284  for(cHadron = 0; cHadron < Output->size(); cHadron++)
285  {
286  SumMom += Output->operator[](cHadron)->Get4Momentum();
287  SumMass += Output->operator[](cHadron)->GetDefinition()->GetPDGMass();
288  }
289 
290 #ifdef debug_G4ExcitedStringCorr
291  G4cout<<"Sum part mom "<<SumMom<<" "<<SumMom.mag()<<G4endl
292  <<"Sum str mom "<<TotalCollisionMom<<" "<<TotalCollisionMom.mag()<<G4endl;
293  G4cout<<"SumMass TotalCollisionMass "<<SumMass<<" "<<TotalCollisionMass<<G4endl;
294 #endif
295 
296  // Cannot correct a single particle
297  if (Output->size() < 2) return FALSE;
298 
299  if (SumMass > TotalCollisionMass) return FALSE;
300  SumMass = SumMom.m2();
301  if (SumMass < 0) return FALSE;
302  SumMass = std::sqrt(SumMass);
303 
304  // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s.
305 // G4ThreeVector Beta = -SumMom.boostVector(); // Uzhi 22 June 2014
306  G4ThreeVector Beta = -TotalCollisionMom.boostVector(); // Uzhi 22 June 2014
307  Output->Boost(Beta);
308 
309  // Scale total c.m.s. hadron energy (hadron system mass).
310  // It should be equal interaction mass
311  G4double Scale = 1;
312  G4int cAttempt = 0;
313  G4double Sum = 0;
314  G4bool success = false;
315  for(cAttempt = 0; cAttempt < nAttemptScale; cAttempt++)
316  {
317  Sum = 0;
318  for(cHadron = 0; cHadron < Output->size(); cHadron++)
319  {
320  G4LorentzVector HadronMom = Output->operator[](cHadron)->Get4Momentum();
321  HadronMom.setVect(Scale*HadronMom.vect());
322  G4double E = std::sqrt(HadronMom.vect().mag2() + sqr(Output->operator[](cHadron)->GetDefinition()->GetPDGMass()));
323  HadronMom.setE(E);
324  Output->operator[](cHadron)->Set4Momentum(HadronMom);
325  Sum += E;
326  }
327  Scale = TotalCollisionMass/Sum;
328 #ifdef debug_G4ExcitedStringCorr
329  G4cout << "Scale-1=" << Scale -1
330  << ", TotalCollisionMass=" << TotalCollisionMass
331  << ", Sum=" << Sum
332  << G4endl;
333 #endif
334  if (std::fabs(Scale - 1) <= ErrLimit)
335  {
336  success = true;
337  break;
338  }
339  }
340 #ifdef debug_G4ExcitedStringCorr
341  if(!success)
342  {
343  G4cout << "G4ExcitedStringDecay::EnergyAndMomentumCorrector - Warning"<<G4endl;
344  G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl;
345  G4cout << " Number of secondaries: " << Output->size() << G4endl;
346  G4cout << " Wanted total energy: " << TotalCollisionMom.e() << G4endl;
347  G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl;
348 // throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay failed to correct...");
349  }
350 #endif
351  // Compute c.m.s. interaction velocity and KTV back boost
352  Beta = TotalCollisionMom.boostVector();
353  Output->Boost(Beta);
354 
355  return success;
356  }
G4bool EnergyAndMomentumCorrector(G4KineticTrackVector *Output, G4LorentzVector &TotalCollisionMom)
std::vector< G4ExcitedString * > G4ExcitedStringVector
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepLorentzRotation G4LorentzRotation
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)=0
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)
const G4ExcitedStringDecay & operator=(const G4ExcitedStringDecay &right)
G4VLongitudinalStringDecay * theStringDecay
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetPosition(const G4ThreeVector aPosition)
int operator!=(const G4ExcitedStringDecay &right) const
#define FALSE
Definition: globals.hh:52
#define TRUE
Definition: globals.hh:55
void Boost(G4ThreeVector &Velocity)
static const double perMillion
Definition: G4SIunits.hh:298
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
int operator==(const G4ExcitedStringDecay &right) const
const G4ParticleDefinition * GetDefinition() const
G4GLOB_DLL std::ostream G4cerr
CLHEP::HepLorentzVector G4LorentzVector