Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
54 const G4ExcitedStringDecay & G4ExcitedStringDecay::operator=(const G4ExcitedStringDecay &)
55 {
56  throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::operator= meant to not be accessable");
57  return *this;
58 }
59 
60 int G4ExcitedStringDecay::operator==(const G4ExcitedStringDecay &) const
61 {
62  return 0;
63 }
64 
65 int G4ExcitedStringDecay::operator!=(const G4ExcitedStringDecay &) const
66 {
67  return 1;
68 }
69 
70 G4KineticTrackVector *G4ExcitedStringDecay::FragmentString(const G4ExcitedString &theString)
71 {
72  if ( theStringDecay == NULL ) theStringDecay=new G4LundStringFragmentation();
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 {
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 }
267 
268 G4bool G4ExcitedStringDecay::EnergyAndMomentumCorrector
269  (G4KineticTrackVector* Output, G4LorentzVector& TotalCollisionMom)
270 {
271  const int nAttemptScale = 500;
272  const double ErrLimit = 1.E-5;
273  if (Output->empty()) return TRUE;
274  G4LorentzVector SumMom;
275  G4double SumMass = 0;
276  G4double TotalCollisionMass = TotalCollisionMom.m();
277 
278  #ifdef debug_G4ExcitedStringCorr
279  G4cout<<G4endl<<"EnergyAndMomentumCorrector. Number of particles: "<<Output->size()<<G4endl;
280  #endif
281  // Calculate sum hadron 4-momenta and summing hadron mass
282  unsigned int cHadron;
283  for(cHadron = 0; cHadron < Output->size(); cHadron++)
284  {
285  SumMom += Output->operator[](cHadron)->Get4Momentum();
286  SumMass += Output->operator[](cHadron)->GetDefinition()->GetPDGMass();
287  }
288 
289  #ifdef debug_G4ExcitedStringCorr
290  G4cout<<"Sum part mom "<<SumMom<<" "<<SumMom.mag()<<G4endl
291  <<"Sum str mom "<<TotalCollisionMom<<" "<<TotalCollisionMom.mag()<<G4endl;
292  G4cout<<"SumMass TotalCollisionMass "<<SumMass<<" "<<TotalCollisionMass<<G4endl;
293  #endif
294 
295  // Cannot correct a single particle
296  if (Output->size() < 2) return FALSE;
297 
298  if (SumMass > TotalCollisionMass) return FALSE;
299  SumMass = SumMom.m2();
300  if (SumMass < 0) return FALSE;
301  SumMass = std::sqrt(SumMass);
302 
303  // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s.
304  //G4ThreeVector Beta = -SumMom.boostVector();
305  G4ThreeVector Beta = -TotalCollisionMom.boostVector();
306  Output->Boost(Beta);
307 
308  // Scale total c.m.s. hadron energy (hadron system mass).
309  // It should be equal interaction mass
310  G4double Scale = 1;
311  G4int cAttempt = 0;
312  G4double Sum = 0;
313  G4bool success = false;
314  for(cAttempt = 0; cAttempt < nAttemptScale; cAttempt++)
315  {
316  Sum = 0;
317  for(cHadron = 0; cHadron < Output->size(); cHadron++)
318  {
319  G4LorentzVector HadronMom = Output->operator[](cHadron)->Get4Momentum();
320  HadronMom.setVect(Scale*HadronMom.vect());
321  G4double E = std::sqrt(HadronMom.vect().mag2() + sqr(Output->operator[](cHadron)->GetDefinition()->GetPDGMass()));
322  HadronMom.setE(E);
323  Output->operator[](cHadron)->Set4Momentum(HadronMom);
324  Sum += E;
325  }
326  Scale = TotalCollisionMass/Sum;
327  #ifdef debug_G4ExcitedStringCorr
328  G4cout << "Scale-1=" << Scale -1
329  << ", TotalCollisionMass=" << TotalCollisionMass
330  << ", Sum=" << Sum
331  << G4endl;
332  #endif
333  if (std::fabs(Scale - 1) <= ErrLimit)
334  {
335  success = true;
336  break;
337  }
338  }
339  #ifdef debug_G4ExcitedStringCorr
340  if(!success)
341  {
342  G4cout << "G4ExcitedStringDecay::EnergyAndMomentumCorrector - Warning"<<G4endl;
343  G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl;
344  G4cout << " Number of secondaries: " << Output->size() << G4endl;
345  G4cout << " Wanted total energy: " << TotalCollisionMom.e() << G4endl;
346  G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl;
347  //throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay failed to correct...");
348  }
349  #endif
350  // Compute c.m.s. interaction velocity and KTV back boost
351  Beta = TotalCollisionMom.boostVector();
352  Output->Boost(Beta);
353 
354  return success;
355 }
356 
static constexpr double perMillion
Definition: G4SIunits.hh:334
Hep3Vector boostVector() const
std::vector< G4ExcitedString * > G4ExcitedStringVector
CLHEP::Hep3Vector G4ThreeVector
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)=0
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
Hep3Vector vect() const
G4GLOB_DLL std::ostream G4cout
double mag() const
bool G4bool
Definition: G4Types.hh:79
void SetPosition(const G4ThreeVector aPosition)
#define FALSE
Definition: globals.hh:52
#define TRUE
Definition: globals.hh:55
void Boost(G4ThreeVector &Velocity)
double m2() const
double mag2() const
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
void setVect(const Hep3Vector &)
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * GetDefinition() const
CLHEP::HepLorentzVector G4LorentzVector