Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ExcitedStringDecay.cc
이 파일의 문서화 페이지로 가기
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 #include "G4SampleResonance.hh" // Uzhi July 2017
33 
34 //#define debug_G4ExcitedStringDecay
35 //#define debug_G4ExcitedStringCorr
36 
38 {}
39 
42  theStringDecay(aStringDecay)
43 {}
44 
47  theStringDecay(0)
48 {
49  throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::copy ctor not accessible");
50 }
51 
53 {
54 }
55 
57 {
58  throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::operator= meant to not be accessable");
59  return *this;
60 }
61 
63 {
64  return 0;
65 }
66 
68 {
69  return 1;
70 }
71 
73 {
75  return theStringDecay->FragmentString(theString);
76 }
77 
79 {
80  G4LorentzVector KTsum(0.,0.,0.,0.);
81 
82 #ifdef debug_G4ExcitedStringDecay
83  G4cout<<G4endl;
84  G4cout<<"--------------------------- G4ExcitedStringDecay ----------------------"<<G4endl;
85  G4cout<<"Hadronization of Excited Strings: theStrings->size() "<<theStrings->size()<<G4endl;
86 #endif
87 
88  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
89 // for ( unsigned int astring=0; astring < 1; astring++)
90  {
91 //G4cout<<"theStrings->operator[](astring)->IsExcited() "<<" "<<astring<<" "<<theStrings->operator[](astring)->IsExcited()<<G4endl;
92  if ( theStrings->operator[](astring)->IsExcited() )
93  {KTsum+= theStrings->operator[](astring)->Get4Momentum();}
94  else {KTsum+=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();}
95  }
96 
97  G4LorentzRotation toCms( -1 * KTsum.boostVector() );
98  G4LorentzRotation toLab(toCms.inverse());
99  G4LorentzVector Ptmp;
100  KTsum=G4LorentzVector(0.,0.,0.,0.);
101 
102  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
103 // for ( unsigned int astring=0; astring < 1; astring++)
104  {
105  if ( theStrings->operator[](astring)->IsExcited() )
106  {
107  Ptmp=toCms * theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
108  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
109 
110  Ptmp=toCms * theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
111  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
112 
113  KTsum+= theStrings->operator[](astring)->Get4Momentum();
114  }
115  else
116  {
117  Ptmp=toCms * theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
118  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
119  KTsum+= theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
120  }
121  }
122 
123  G4SampleResonance BrW; // Uzhi July 2017
124  const G4ParticleDefinition* TrackDefinition=0; // Uzhi July 2017
125 
126  G4KineticTrackVector * theResult = new G4KineticTrackVector;
127  G4int attempts(0);
128  G4bool success=false;
129  G4bool NeedEnergyCorrector=false;
130  do {
131 #ifdef debug_G4ExcitedStringDecay
132  G4cout<<"New try No "<<attempts<<" to hadronize strings"<<G4endl;
133 #endif
134 
135  std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack());
136  theResult->clear();
137 
138  attempts++;
139 
140  G4LorentzVector KTsecondaries(0.,0.,0.,0.);
141  NeedEnergyCorrector=false;
142 
143  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
144 // for ( unsigned int astring=0; astring < 1; astring++) // Uzhi For testing purposes
145  {
146 #ifdef debug_G4ExcitedStringDecay
147  G4cout<<"String No "<<astring+1<<" Excited? "<<theStrings->operator[](astring)->IsExcited()<<G4endl;
148 
149  G4cout<<"String No "<<astring+1<<" 4Momentum "<<theStrings->operator[](astring)->Get4Momentum()
150  <<" "<<theStrings->operator[](astring)->Get4Momentum().mag()<<G4endl;
151 #endif
152 
153  G4KineticTrackVector * generatedKineticTracks = NULL;
154  if ( theStrings->operator[](astring)->IsExcited() )
155  {
156 #ifdef debug_G4ExcitedStringDecay
157  G4cout<<"Fragment String with partons: "
158  <<theStrings->operator[](astring)->GetLeftParton()->GetPDGcode() <<" "
159  <<theStrings->operator[](astring)->GetRightParton()->GetPDGcode()<<" "
160  <<"Direction "<<theStrings->operator[](astring)->GetDirection()<<G4endl;
161 #endif
162  generatedKineticTracks=FragmentString(*theStrings->operator[](astring));
163 #ifdef debug_G4ExcitedStringDecay
164  G4cout<<"(G4ExcitedStringDecay) Number of produced hadrons = "
165  <<generatedKineticTracks->size()<<G4endl;
166 #endif
167  } else {
168 #ifdef debug_G4ExcitedStringDecay
169  G4cout<<" GetTrack from the String"<<G4endl;
170 #endif
171  G4LorentzVector Mom=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
172  G4KineticTrack * aTrack= new G4KineticTrack(
173  theStrings->operator[](astring)->GetKineticTrack()->GetDefinition(),
174  theStrings->operator[](astring)->GetKineticTrack()->GetFormationTime(),
175  G4ThreeVector(0), Mom);
176 
177  aTrack->SetPosition(theStrings->operator[](astring)->GetKineticTrack()->GetPosition());
178 
179 #ifdef debug_G4ExcitedStringDecay
180  G4cout<<" A particle stored in the track is "<<aTrack->GetDefinition()->GetParticleName()<<G4endl;
181 #endif
182 
183  generatedKineticTracks = new G4KineticTrackVector;
184  generatedKineticTracks->push_back(aTrack);
185  }
186 
187  if (generatedKineticTracks->size() == 0)
188  {
189 // G4cerr << "G4VPartonStringModel:No KineticTracks produced" << G4endl;
190 // continue;
191  success=false; NeedEnergyCorrector=false; break;
192  }
193 
194  G4LorentzVector KTsum1(0.,0.,0.,0.);
195  for ( unsigned int aTrack=0; aTrack<generatedKineticTracks->size();aTrack++)
196  {
197 #ifdef debug_G4ExcitedStringDecay
198  G4cout<<"Prod part No. "<<aTrack+1<<" "
199  <<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" "
200  <<(*generatedKineticTracks)[aTrack]->Get4Momentum()
201  <<(*generatedKineticTracks)[aTrack]->Get4Momentum().mag()<<G4endl;
202 #endif
203 // --------------- Sampling mass of unstable hadronic resonances ---------------- Uzhi July 2017
204  TrackDefinition = (*generatedKineticTracks)[aTrack]->GetDefinition();
205 
206  if(TrackDefinition->IsShortLived())
207  {
208  G4double NewTrackMass = BrW.SampleMass( TrackDefinition,
209  TrackDefinition->GetPDGMass() + 5.0*TrackDefinition->GetPDGWidth() );
210  G4LorentzVector Tmp=G4LorentzVector((*generatedKineticTracks)[aTrack]->Get4Momentum());
211  Tmp.setE(std::sqrt(sqr(NewTrackMass) + Tmp.vect().mag2()));
212 
213  (*generatedKineticTracks)[aTrack]->Set4Momentum(Tmp);
214 
215 #ifdef debug_G4ExcitedStringDecay
216  G4cout<<"Resonance *** "<<aTrack+1<<" "
217  <<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" "
218  <<(*generatedKineticTracks)[aTrack]->Get4Momentum()
219  <<(*generatedKineticTracks)[aTrack]->Get4Momentum().mag()<<G4endl;
220 #endif
221 
222  }
223 // //----------------------------------------------------------------------------------------------
224 
225  theResult->push_back(generatedKineticTracks->operator[](aTrack));
226  KTsum1+= (*generatedKineticTracks)[aTrack]->Get4Momentum();
227  }
228  KTsecondaries+=KTsum1;
229 
230 #ifdef debug_G4ExcitedStringDecay
231  G4cout << "String secondaries(" <<generatedKineticTracks->size()<< ")"<<G4endl
232  <<"Init string momentum: "<< theStrings->operator[](astring)->Get4Momentum()<<G4endl
233  <<"Final hadrons momentum: "<< KTsum1 << G4endl;
234 #endif
235 
236  if ( KTsum1.e() > 0 && std::abs((KTsum1.e()-theStrings->operator[](astring)->Get4Momentum().e()) / KTsum1.e()) > perMillion )
237  {
238  NeedEnergyCorrector=true;
239  }
240 
241 #ifdef debug_G4ExcitedStringDecay
242  G4cout<<"NeedEnergyCorrection yes/no "<<NeedEnergyCorrector<<G4endl;
243 #endif
244 
245 // clean up
246  delete generatedKineticTracks;
247  success=true;
248  }
249 
250  if ( NeedEnergyCorrector ) success=EnergyAndMomentumCorrector(theResult, KTsum);
251  } while(!success && (attempts < 100)); /* Loop checking, 07.08.2015, A.Ribon */
252 
253  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
254  {
255  Ptmp=(*theResult)[aTrack]->Get4Momentum();
256  Ptmp.transform( toLab);
257  (*theResult)[aTrack]->Set4Momentum(Ptmp);
258  }
259 
260 #ifdef debug_G4ExcitedStringDecay
261  G4cout<<"End of the strings fragmentation (G4ExcitedStringDecay)"<<G4endl;
262 
263  G4LorentzVector KTsum1(0.,0.,0.,0.);
264 
265  for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
266  {
267  G4cout << " corrected tracks .. " << (*theResult)[aTrack]->GetDefinition()->GetParticleName()
268  <<" " << (*theResult)[aTrack]->Get4Momentum()
269  <<" " << (*theResult)[aTrack]->Get4Momentum().mag()<< G4endl;
270  KTsum1+= (*theResult)[aTrack]->Get4Momentum();
271  }
272 
273  G4cout << "Needcorrector/success " << NeedEnergyCorrector << "/" << success
274  << ", Corrected total 4 momentum " << KTsum1 << G4endl;
275  if ( ! success ) G4cout << "failed to correct E/p" << G4endl;
276 
277  G4cout<<"End of the Hadronization (G4ExcitedStringDecay)"<<G4endl;
278 #endif
279 
280 if(!success)
281 {
282  if(theResult->size() != 0)
283  {std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack());
284  theResult->clear();
285  delete theResult; theResult=0;
286  }
287  for ( unsigned int astring=0; astring < theStrings->size(); astring++)
288 // for ( unsigned int astring=0; astring < 1; astring++) // Uzhi 2016 Need more correct. For testing purposes.
289  {
290  if ( theStrings->operator[](astring)->IsExcited() )
291  {
292  Ptmp=theStrings->operator[](astring)->GetLeftParton()->Get4Momentum();
293  Ptmp.transform( toLab);
294  theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp);
295 
296  Ptmp=theStrings->operator[](astring)->GetRightParton()->Get4Momentum();
297  Ptmp.transform( toLab);
298  theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp);
299  }
300  else
301  {
302  Ptmp=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();
303  Ptmp.transform( toLab);
304  theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp);
305  }
306  }
307 }
308  return theResult;
309 }
310 
312  (G4KineticTrackVector* Output, G4LorentzVector& TotalCollisionMom)
313  {
314  const int nAttemptScale = 500;
315  const double ErrLimit = 1.E-5;
316  if (Output->empty()) return TRUE;
317  G4LorentzVector SumMom;
318  G4double SumMass = 0;
319  G4double TotalCollisionMass = TotalCollisionMom.m();
320 
321  std::vector<G4double> HadronMass; G4double HadronM(0.);
322 
323 #ifdef debug_G4ExcitedStringCorr
324  G4cout<<G4endl<<"EnergyAndMomentumCorrector. Number of particles: "<<Output->size()<<G4endl;
325 #endif
326  // Calculate sum hadron 4-momenta and summing hadron mass
327  unsigned int cHadron;
328  for(cHadron = 0; cHadron < Output->size(); cHadron++)
329  {
330  SumMom += Output->operator[](cHadron)->Get4Momentum();
331  HadronM=Output->operator[](cHadron)->Get4Momentum().mag(); HadronMass.push_back(HadronM);
332  SumMass += Output->operator[](cHadron)->Get4Momentum().mag(); //GetDefinition()->GetPDGMass();
333  }
334 
335 #ifdef debug_G4ExcitedStringCorr
336  G4cout<<"Sum part mom "<<SumMom<<" "<<SumMom.mag()<<G4endl
337  <<"Sum str mom "<<TotalCollisionMom<<" "<<TotalCollisionMom.mag()<<G4endl;
338  G4cout<<"SumMass TotalCollisionMass "<<SumMass<<" "<<TotalCollisionMass<<G4endl;
339 #endif
340 
341  // Cannot correct a single particle
342  if (Output->size() < 2) return FALSE;
343 
344  if (SumMass > TotalCollisionMass) return FALSE;
345  SumMass = SumMom.m2();
346  if (SumMass < 0) return FALSE;
347  SumMass = std::sqrt(SumMass);
348 
349  // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s.
350 // G4ThreeVector Beta = -SumMom.boostVector();
351  G4ThreeVector Beta = -TotalCollisionMom.boostVector();
352  Output->Boost(Beta);
353 
354  // Scale total c.m.s. hadron energy (hadron system mass).
355  // It should be equal interaction mass
356  G4double Scale = 1;
357  G4int cAttempt = 0;
358  G4double Sum = 0;
359  G4bool success = false;
360  for(cAttempt = 0; cAttempt < nAttemptScale; cAttempt++)
361  {
362  Sum = 0;
363  for(cHadron = 0; cHadron < Output->size(); cHadron++)
364  {
365  HadronM = HadronMass.at(cHadron);
366  G4LorentzVector HadronMom = Output->operator[](cHadron)->Get4Momentum();
367  HadronMom.setVect(Scale*HadronMom.vect());
368  G4double E = std::sqrt(HadronMom.vect().mag2() + sqr(HadronM));
369  //sqr(Output->operator[](cHadron)->GetDefinition()->GetPDGMass()));
370  HadronMom.setE(E);
371  Output->operator[](cHadron)->Set4Momentum(HadronMom);
372  Sum += E;
373  }
374  Scale = TotalCollisionMass/Sum;
375 #ifdef debug_G4ExcitedStringCorr
376  G4cout << "Scale-1=" << Scale -1
377  << ", TotalCollisionMass=" << TotalCollisionMass
378  << ", Sum=" << Sum
379  << G4endl;
380 #endif
381  if (std::fabs(Scale - 1) <= ErrLimit)
382  {
383  success = true;
384  break;
385  }
386  }
387 
388 #ifdef debug_G4ExcitedStringCorr
389  if(!success)
390  {
391  G4cout << "G4ExcitedStringDecay::EnergyAndMomentumCorrector - Warning"<<G4endl;
392  G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl;
393  G4cout << " Number of secondaries: " << Output->size() << G4endl;
394  G4cout << " Wanted total energy: " << TotalCollisionMom.e() << G4endl;
395  G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl;
396 // throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay failed to correct...");
397  }
398 #endif
399  // Compute c.m.s. interaction velocity and KTV back boost
400  Beta = TotalCollisionMom.boostVector();
401  Output->Boost(Beta);
402 
403  return success;
404  }
const G4ExcitedStringDecay & operator=(const G4ExcitedStringDecay &right)
G4double GetPDGWidth() const
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4ExcitedString * > G4ExcitedStringVector
#define G4endl
Definition: G4ios.hh:61
void Boost(G4ThreeVector &Velocity)
const G4String & GetParticleName() const
const G4ParticleDefinition * GetDefinition() const
static constexpr double perMillion
Definition: G4SIunits.hh:334
void setVect(const Hep3Vector &)
G4double SampleMass(const G4double poleMass, const G4double gamma, const G4double minMass, const G4double maxMass) const
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
int operator!=(const G4ExcitedStringDecay &right) const
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)
G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)
#define FALSE
Definition: globals.hh:52
#define TRUE
Definition: globals.hh:55
double mag2() const
double m2() const
G4bool EnergyAndMomentumCorrector(G4KineticTrackVector *Output, G4LorentzVector &TotalCollisionMom)
double mag() const
int G4int
Definition: G4Types.hh:78
void SetPosition(const G4ThreeVector aPosition)
int operator==(const G4ExcitedStringDecay &right) const
G4VLongitudinalStringDecay * theStringDecay
G4GLOB_DLL std::ostream G4cout
Hep3Vector boostVector() const
Hep3Vector vect() const
T sqr(const T &x)
Definition: templates.hh:145
CLHEP::HepLorentzVector G4LorentzVector
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)=0
simulatedPeak Scale(1/simulationNormalisationFactor)