Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VPartonStringModel.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 //
27 // $Id: G4VPartonStringModel.cc 107318 2017-11-08 16:27:32Z gcosmo $
28 //
30 // GEANT 4 class implementation file
31 //
32 // ---------------- G4VPartonStringModel ----------------
33 // by Gunter Folger, May 1998.
34 // abstract class for all Parton String Models
35 // ------------------------------------------------------------
36 // debug switch
37 //#define debug_PartonStringModel
38 // ------------------------------------------------------------
39 
40 #include "G4VPartonStringModel.hh"
41 #include "G4ios.hh"
43 
44 #include "G4ParticleTable.hh"
45 #include "G4IonTable.hh"
46 
48  : G4VHighEnergyGenerator(modelName),
49  stringFragmentationModel(0),
50  theThis(0)
51 {
52 // Make shure Shotrylived particles are constructed.
53  G4ShortLivedConstructor ShortLived;
54  ShortLived.ConstructParticle();
55 }
56 
58 {
59 }
60 
62  const G4DynamicParticle &aPrimary)
63 {
64  G4ExcitedStringVector * strings = NULL;
65  G4DynamicParticle thePrimary=aPrimary;
66  G4LorentzVector SumStringMom(0.,0.,0.,0.);
67  G4KineticTrackVector * theResult = 0;
68  G4Nucleon * theNuclNucleon(0);
69 
70  #ifdef debug_PartonStringModel
71  G4cout<<G4endl;
72  G4cout<<"-----------------------Parton-String model is runnung ------------"<<G4endl;
73  G4cout<<"Projectile Name Mass "<<thePrimary.GetDefinition()->GetParticleName()<<" "
74  <<thePrimary.GetMass()<<G4endl;
75  G4cout<<" Momentum "<<thePrimary.Get4Momentum()<<G4endl;
76  G4cout<<"Target nucleus A Z "<<theNucleus.GetA_asInt()<<" "
77  <<theNucleus.GetZ_asInt()<<G4endl<<G4endl;
78  G4int Bsum=thePrimary.GetDefinition()->GetBaryonNumber() + theNucleus.GetA_asInt();
79  G4int Qsum=thePrimary.GetDefinition()->GetPDGCharge() + theNucleus.GetZ_asInt();
80  G4cout<<"Initial baryon number "<<Bsum<<G4endl;
81  G4cout<<"Initial charge "<<Qsum<<G4endl;
82  G4cout<<"-------------- Parton-String model: Generation of strings -------"<<G4endl<<G4endl;
83  Bsum -= theNucleus.GetA_asInt(); Qsum -= theNucleus.GetZ_asInt();
85  Bsum -= thePrimary.GetDefinition()->GetBaryonNumber();
86  Qsum -= thePrimary.GetDefinition()->GetPDGCharge();
87  }
88  G4int QsumSec(0), BsumSec(0);
89  G4LorentzVector SumPsecondr(0.,0.,0.,0.);
90  #endif
91 
93  G4LorentzVector Ptmp=thePrimary.Get4Momentum();
94  toZ.rotateZ(-1*Ptmp.phi());
95  toZ.rotateY(-1*Ptmp.theta());
96  thePrimary.Set4Momentum(toZ*Ptmp);
97  G4LorentzRotation toLab(toZ.inverse());
98 
99  G4bool Success=true;
100  G4int attempts = 0, maxAttempts=1000;
101  do
102  {
103  if (attempts++ > maxAttempts )
104  {
105  theThis->Init(theNucleus,thePrimary); // To put a nucleus into ground state
106  // But marks of hitted nucleons are left. They must be erased.
107  G4V3DNucleus * ResNucleus=theThis->GetWoundedNucleus();
108  theNuclNucleon = ResNucleus->StartLoop() ? ResNucleus->GetNextNucleon() : NULL;
109  while( theNuclNucleon )
110  {
111  if(theNuclNucleon->AreYouHit()) theNuclNucleon->Hit(nullptr);
112  theNuclNucleon = ResNucleus->GetNextNucleon();
113  }
114 
115  G4V3DNucleus * ProjResNucleus=theThis->GetProjectileNucleus();
116  if(ProjResNucleus != 0)
117  {
118  theNuclNucleon = ProjResNucleus->StartLoop() ? ProjResNucleus->GetNextNucleon() : NULL;
119  while( theNuclNucleon )
120  {
121  if(theNuclNucleon->AreYouHit()) theNuclNucleon->Hit(nullptr);
122  theNuclNucleon = ProjResNucleus->GetNextNucleon();
123  }
124  }
125 
127  ed << "Projectile Name Mass " <<thePrimary.GetDefinition()->GetParticleName()
128  << " " << thePrimary.GetMass()<< G4endl;
129  ed << "           Momentum  " << thePrimary.Get4Momentum() <<G4endl;
130  ed << "Target nucleus   A Z " << theNucleus.GetA_asInt() << " "
131  << theNucleus.GetZ_asInt() <<G4endl;
132  ed << "Initial states of projectile and target nucleus will be returned!"<<G4endl;
133  G4Exception( "G4VPartonStringModel::Scatter(): fails to generate or fragment strings ",
134  "HAD_PARTON_STRING_001", JustWarning, ed );
135 
136  G4ThreeVector Position(0.,0.,2*ResNucleus->GetOuterRadius());
137  G4KineticTrack* Hadron = new G4KineticTrack(aPrimary.GetParticleDefinition(), 0.,
138  Position, aPrimary.Get4Momentum());
139  if(theResult == 0) theResult = new G4KineticTrackVector();
140  theResult->push_back(Hadron);
141  return theResult;
142  }
143 
144  Success=true;
145 
146  theThis->Init(theNucleus,thePrimary);
147 
148  strings = GetStrings();
149 
150  if(strings == 0) {Success=false; continue;}
151 
152  G4double stringEnergy(0);
153  SumStringMom=G4LorentzVector(0.,0.,0.,0.);
154 
155  #ifdef debug_PartonStringModel
156  G4cout<<"------------ Parton-String model: Number of produced strings ---- "<<strings->size()<<G4endl;
157  #endif
158 
159  for ( unsigned int astring=0; astring < strings->size(); astring++)
160  {
161  // rotate string to lab frame, models have it aligned to z
162  if((*strings)[astring]->IsExcited())
163  {
164  stringEnergy += (*strings)[astring]->GetLeftParton()->Get4Momentum().t();
165  stringEnergy += (*strings)[astring]->GetRightParton()->Get4Momentum().t();
166  (*strings)[astring]->LorentzRotate(toLab);
167  SumStringMom+=(*strings)[astring]->Get4Momentum();
168  #ifdef debug_PartonStringModel
169  G4cout<<"String No "<<astring+1<<" "<<(*strings)[astring]->Get4Momentum()<<" "
170  <<(*strings)[astring]->Get4Momentum().mag()
171  <<" Partons "<<(*strings)[astring]->GetLeftParton()->GetDefinition()->GetPDGEncoding()
172  <<" "<<(*strings)[astring]->GetRightParton()->GetDefinition()->GetPDGEncoding()<<G4endl;
173  #endif
174  }
175  else
176  {
177  stringEnergy += (*strings)[astring]->GetKineticTrack()->Get4Momentum().t();
178  (*strings)[astring]->LorentzRotate(toLab);
179  SumStringMom+=(*strings)[astring]->GetKineticTrack()->Get4Momentum();
180  #ifdef debug_PartonStringModel
181  G4cout<<"A track No "<<astring+1<<" "
182  <<(*strings)[astring]->GetKineticTrack()->Get4Momentum()<<" "
183  <<(*strings)[astring]->GetKineticTrack()->Get4Momentum().mag()<<" "
184  <<(*strings)[astring]->GetKineticTrack()->GetDefinition()->GetParticleName()<<G4endl;
185  #endif
186  }
187  }
188 
189  #ifdef debug_PartonStringModel
190  G4cout<<G4endl<<"SumString4Mom "<<SumStringMom<<G4endl;
191  G4LorentzVector TargetResidual4Momentum(0.,0.,0.,0.);
192  G4LorentzVector ProjectileResidual4Momentum(0.,0.,0.,0.);
193  G4int hitsT(0), charged_hitsT(0);
194  G4int hitsP(0), charged_hitsP(0);
195  G4double ExcitationEt(0.), ExcitationEp(0.);
196  #endif
197 
198  G4V3DNucleus * ProjResNucleus=theThis->GetProjectileNucleus();
199 
200  if(ProjResNucleus != 0)
201  {
202  theNuclNucleon = ProjResNucleus->StartLoop() ? ProjResNucleus->GetNextNucleon() : NULL;
203  while( theNuclNucleon )
204  {
205  if(theNuclNucleon->AreYouHit())
206  {
207  G4LorentzVector tmp=toLab*theNuclNucleon->Get4Momentum();
208  #ifdef debug_PartonStringModel
209  ProjectileResidual4Momentum += tmp;
210  hitsP++;
211  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsP;
212  ExcitationEp +=theNuclNucleon->GetBindingEnergy();
213  #endif
214  theNuclNucleon->SetMomentum(tmp);
215  }
216  theNuclNucleon = ProjResNucleus->GetNextNucleon();
217  }
218  #ifdef debug_PartonStringModel
219  G4cout<<"Projectile residual A, Z and E* "
220  <<thePrimary.GetDefinition()->GetBaryonNumber() - hitsP<<" "
221  <<thePrimary.GetDefinition()->GetPDGCharge() - charged_hitsP<<" "
222  <<ExcitationEp<<G4endl;
223  G4cout<<"Projectile residual 4 momentum "<<ProjectileResidual4Momentum<<G4endl;
224  #endif
225  }
226 
227  G4V3DNucleus * ResNucleus=theThis->GetWoundedNucleus();
228 
229  // loop over wounded nucleus
230  theNuclNucleon = ResNucleus->StartLoop() ? ResNucleus->GetNextNucleon() : NULL;
231  while( theNuclNucleon )
232  {
233  if(theNuclNucleon->AreYouHit())
234  {
235  G4LorentzVector tmp=toLab*theNuclNucleon->Get4Momentum();
236  #ifdef debug_PartonStringModel
237  TargetResidual4Momentum += tmp;
238  hitsT++;
239  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsT;
240  ExcitationEt +=theNuclNucleon->GetBindingEnergy();
241  #endif
242  theNuclNucleon->SetMomentum(tmp);
243  }
244  theNuclNucleon = ResNucleus->GetNextNucleon();
245  }
246 
247  #ifdef debug_PartonStringModel
248  G4cout<<"Target residual A, Z and E* "
249  <<theNucleus.GetA_asInt() - hitsT<<" "
250  <<theNucleus.GetZ_asInt() - charged_hitsT<<" "
251  <<ExcitationEt<<G4endl;
252  G4cout<<"Target residual 4 momentum "<<TargetResidual4Momentum<<G4endl;
253  Bsum+=( hitsT + hitsP);
254  Qsum+=(charged_hitsT + charged_hitsP);
255  G4cout<<"Hitted # of nucleons of projectile and target "<<hitsP<<" "<<hitsT<<G4endl;
256  G4cout<<"Hitted # of protons of projectile and target "
257  <<charged_hitsP<<" "<<charged_hitsT<<G4endl<<G4endl;
258  G4cout<<"Bsum Qsum "<<Bsum<<" "<<Qsum<<G4endl<<G4endl;
259  #endif
260 
261  //=========================================================================================
262  // Fragment strings
263  #ifdef debug_PartonStringModel
264  G4cout<<"---------------- Attempt to fragment strings ------------- "<<attempts<<G4endl;
265  #endif
266 
267  G4double InvMass=SumStringMom.mag();
268  G4double SumMass(0.);
269 
270  #ifdef debug_PartonStringModel
271  QsumSec=0; BsumSec=0;
272  SumPsecondr=G4LorentzVector(0.,0.,0.,0.);
273  #endif
274 
275  if(theResult != 0)
276  {
277  std::for_each(theResult->begin(), theResult->end(), DeleteKineticTrack());
278  delete theResult;
279  }
280 
281  theResult = stringFragmentationModel->FragmentStrings(strings);
282 
283  #ifdef debug_PartonStringModel
284  G4cout<<"String fragmentation success (OK - #0, Not OK - 0) : "<<theResult<<G4endl;
285  #endif
286 
287  if(theResult == 0) {Success=false; continue;}
288 
289  #ifdef debug_PartonStringModel
290  G4cout<<"Parton-String model: Number of produced particles "<<theResult->size()<<G4endl;
291  SumPsecondr = G4LorentzVector(0.,0.,0.,0.);
292  QsumSec = 0; BsumSec = 0;
293  #endif
294 
295  SumMass=0.;
296 
297  for ( unsigned int i=0; i < theResult->size(); i++)
298  {
299  SumMass+=(*theResult)[i]->Get4Momentum().mag();
300  #ifdef debug_PartonStringModel
301  G4cout<<i<<" : "<<(*theResult)[i]->GetDefinition()->GetParticleName()<<" "
302  <<(*theResult)[i]->Get4Momentum()<<" "
303  <<(*theResult)[i]->Get4Momentum().mag()<<" "
304  <<(*theResult)[i]->GetDefinition()->GetPDGMass()<<G4endl;
305  SumPsecondr+=(*theResult)[i]->Get4Momentum();
306  BsumSec += (*theResult)[i]->GetDefinition()->GetBaryonNumber();
307  QsumSec += (*theResult)[i]->GetDefinition()->GetPDGCharge();
308  #endif
309  }
310 
311  #ifdef debug_PartonStringModel
312  G4cout<<G4endl<<"-----------------------Parton-String model: balances -------------"<<G4endl;
313  if(Qsum != QsumSec) {
314  G4cout<<"Charge is not conserved!!! ----"<<G4endl;
315  G4cout<<" Qsum != QsumSec "<<Qsum<<" "<<QsumSec<<G4endl;
316  }
317  if(Bsum != BsumSec) {
318  G4cout<<"Baryon number is not conserved!!!"<<G4endl;
319  G4cout<<" Bsum != BsumSec "<<Bsum<<" "<<BsumSec<<G4endl;
320  }
321  #endif
322 
323  if((SumMass > InvMass)||(SumMass == 0.)) {Success=false;}
324  std::for_each(strings->begin(), strings->end(), DeleteString() );
325  delete strings;
326 
327  } while(!Success);
328 
329  #ifdef debug_PartonStringModel
330  G4cout <<"Baryon number balance "<<Bsum-BsumSec<<G4endl;
331  G4cout <<"Charge balance "<<Qsum-QsumSec<<G4endl;
332  G4cout <<"4 momentum balance "<<SumStringMom-SumPsecondr<<G4endl;
333  G4cout<<"---------------------End of Parton-String model work -------------"<<G4endl<<G4endl;
334  #endif
335 
336  return theResult;
337 }
338 
339 void G4VPartonStringModel::ModelDescription(std::ostream& outFile) const
340 {
341  outFile << GetModelName() << " has no description yet.\n";
342 }
343 
345 { return 0; }
346 
virtual void Init(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
void SetMomentum(G4LorentzVector &aMomentum)
Definition: G4Nucleon.hh:70
std::vector< G4ExcitedString * > G4ExcitedStringVector
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
double phi() const
G4double GetPDGCharge() const
Float_t tmp
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual void ModelDescription(std::ostream &outFile) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
double theta() const
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
HepLorentzRotation & rotateZ(double delta)
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4String GetModelName() const
void Set4Momentum(const G4LorentzVector &momentum)
G4LorentzVector Get4Momentum() const
virtual G4V3DNucleus * GetProjectileNucleus() const
HepLorentzRotation inverse() const
virtual G4Nucleon * GetNextNucleon()=0
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:84
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)=0
G4VPartonStringModel(const G4String &modelName="Parton String Model")
G4VStringFragmentation * stringFragmentationModel
double mag() const
int G4int
Definition: G4Types.hh:78
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:89
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4VPartonStringModel * theThis
G4GLOB_DLL std::ostream G4cout
G4bool AreYouHit() const
Definition: G4Nucleon.hh:96
virtual G4ExcitedStringVector * GetStrings()=0
G4double GetMass() const
CLHEP::HepLorentzVector G4LorentzVector
HepLorentzRotation & rotateY(double delta)
const G4ParticleDefinition * GetParticleDefinition() const
virtual G4double GetOuterRadius()=0
virtual G4bool StartLoop()=0
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115