Geant4  10.02.p02
G4VPartonStringModel.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 //
27 // $Id: G4VPartonStringModel.cc 91858 2015-08-07 13:57:22Z 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 partyicles are constructed.
53  G4ShortLivedConstructor ShortLived;
54  ShortLived.ConstructParticle();
55 }
56 
58 {
59 }
60 
62  const G4DynamicParticle &aPrimary)
63 {
64  G4ExcitedStringVector * strings = NULL;
65 
66  G4DynamicParticle thePrimary=aPrimary;
67 
68 #ifdef debug_PartonStringModel
69  G4cout<<G4endl;
70  G4cout<<"-----------------------Parton-String model is runnung ------------"<<G4endl;
71  G4cout<<"Projectile Name Mass "<<thePrimary.GetDefinition()->GetParticleName()<<" "
72  <<thePrimary.GetMass()<<G4endl;
73  G4cout<<" Momentum "<<thePrimary.Get4Momentum()<<G4endl;
74  G4cout<<"Target nucleus A Z "<<theNucleus.GetA_asInt()<<" "
75  <<theNucleus.GetZ_asInt()<<G4endl<<G4endl;
76  G4int Bsum=thePrimary.GetDefinition()->GetBaryonNumber() + theNucleus.GetA_asInt();
77  G4int Qsum=thePrimary.GetDefinition()->GetPDGCharge() + theNucleus.GetZ_asInt();
78  G4cout<<"Initial baryon number "<<Bsum<<G4endl;
79  G4cout<<"Initial charge "<<Qsum<<G4endl;
80 
81  Bsum -= theNucleus.GetA_asInt(); Qsum -= theNucleus.GetZ_asInt();
83  Bsum -= thePrimary.GetDefinition()->GetBaryonNumber();
84  Qsum -= thePrimary.GetDefinition()->GetPDGCharge();
85  }
86 #endif
87 
89  G4LorentzVector Ptmp=thePrimary.Get4Momentum();
90  toZ.rotateZ(-1*Ptmp.phi());
91  toZ.rotateY(-1*Ptmp.theta());
92  thePrimary.Set4Momentum(toZ*Ptmp);
93  G4LorentzRotation toLab(toZ.inverse());
94 
95  G4int attempts = 0, maxAttempts=20;
96  while ( strings == NULL ) /* Loop checking, 07.08.2015, A.Ribon */
97  {
98  if (attempts++ > maxAttempts )
99  {
100  throw G4HadronicException(__FILE__, __LINE__,
101  "G4VPartonStringModel::Scatter(): fails to generate strings");
102  }
103  theThis->Init(theNucleus,thePrimary);
104 
105  strings = GetStrings();
106  }
107 
108  G4double stringEnergy(0);
109  G4LorentzVector SumStringMom(0.,0.,0.,0.);
110 
111 #ifdef debug_PartonStringModel
112  G4cout<<"Parton-String model: Number of produced strings "<<strings->size()<<G4endl;
113 #endif
114 
115  for ( unsigned int astring=0; astring < strings->size(); astring++)
116  {
117 // rotate string to lab frame, models have it aligned to z
118  if((*strings)[astring]->IsExcited())
119  {
120  stringEnergy += (*strings)[astring]->GetLeftParton()->Get4Momentum().t();
121  stringEnergy += (*strings)[astring]->GetRightParton()->Get4Momentum().t();
122  (*strings)[astring]->LorentzRotate(toLab);
123  SumStringMom+=(*strings)[astring]->Get4Momentum();
124 #ifdef debug_PartonStringModel
125 G4cout<<"String No "<<astring+1<<" "<<(*strings)[astring]->Get4Momentum()<<" "
126  <<(*strings)[astring]->Get4Momentum().mag()
127  <<" Partons "<<(*strings)[astring]->GetLeftParton()->GetDefinition()->GetPDGEncoding()
128  <<" "<<(*strings)[astring]->GetRightParton()->GetDefinition()->GetPDGEncoding()<<G4endl;
129 #endif
130  }
131  else
132  {
133  stringEnergy += (*strings)[astring]->GetKineticTrack()->Get4Momentum().t();
134  (*strings)[astring]->LorentzRotate(toLab);
135  SumStringMom+=(*strings)[astring]->GetKineticTrack()->Get4Momentum();
136 #ifdef debug_PartonStringModel
137  G4cout<<"A track No "<<astring+1<<" "
138  <<(*strings)[astring]->GetKineticTrack()->Get4Momentum()<<" "
139  <<(*strings)[astring]->GetKineticTrack()->Get4Momentum().mag()<<" "
140  <<(*strings)[astring]->GetKineticTrack()->GetDefinition()->GetParticleName()<<G4endl;
141 #endif
142  }
143  }
144 #ifdef debug_PartonStringModel
145  G4cout<<G4endl<<"SumString4Mom "<<SumStringMom<<G4endl;
146 
147  G4LorentzVector TargetResidual4Momentum(0.,0.,0.,0.);
148  G4LorentzVector ProjectileResidual4Momentum(0.,0.,0.,0.);
149  G4int hitsT(0), charged_hitsT(0);
150  G4int hitsP(0), charged_hitsP(0);
151  G4double ExcitationEt(0.), ExcitationEp(0.);
152 #endif
153 
154  G4V3DNucleus * ProjResNucleus=theThis->GetProjectileNucleus();
155  G4Nucleon * theNuclNucleon(0);
156 
157  if(ProjResNucleus != 0)
158  {
159  theNuclNucleon = ProjResNucleus->StartLoop() ?
160  ProjResNucleus->GetNextNucleon() : NULL;
161  while( theNuclNucleon ) /* Loop checking, 07.08.2015, A.Ribon */
162  {
163  if(theNuclNucleon->AreYouHit())
164  {
165  G4LorentzVector tmp=toLab*theNuclNucleon->Get4Momentum();
166  #ifdef debug_PartonStringModel
167  ProjectileResidual4Momentum += tmp;
168  hitsP++;
169  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsP;
170  ExcitationEp +=theNuclNucleon->GetBindingEnergy();
171  #endif
172  theNuclNucleon->SetMomentum(tmp);
173  }
174  theNuclNucleon = ProjResNucleus->GetNextNucleon();
175  }
176 #ifdef debug_PartonStringModel
177  G4cout<<"Projectile residual A, Z and E* "
178  <<thePrimary.GetDefinition()->GetBaryonNumber() - hitsP<<" "
179  <<thePrimary.GetDefinition()->GetPDGCharge() - charged_hitsP<<" "
180  <<ExcitationEp<<G4endl;
181  G4cout<<"Projectile residual 4 momentum "<<ProjectileResidual4Momentum<<G4endl;
182 #endif
183  }
184 
185  G4V3DNucleus * ResNucleus=theThis->GetWoundedNucleus();
186 
187  // loop over wounded nucleus
188  theNuclNucleon = ResNucleus->StartLoop() ?
189  ResNucleus->GetNextNucleon() : NULL;
190  while( theNuclNucleon ) /* Loop checking, 07.08.2015, A.Ribon */
191  {
192  if(theNuclNucleon->AreYouHit())
193  {
194  G4LorentzVector tmp=toLab*theNuclNucleon->Get4Momentum();
195  #ifdef debug_PartonStringModel
196  TargetResidual4Momentum += tmp;
197  hitsT++;
198  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsT;
199  ExcitationEt +=theNuclNucleon->GetBindingEnergy();
200  #endif
201  theNuclNucleon->SetMomentum(tmp);
202  }
203  theNuclNucleon = ResNucleus->GetNextNucleon();
204  }
205 #ifdef debug_PartonStringModel
206  G4cout<<"Target residual A, Z and E* "
207  <<theNucleus.GetA_asInt() - hitsT<<" "
208  <<theNucleus.GetZ_asInt() - charged_hitsT<<" "
209  <<ExcitationEt<<G4endl;
210  G4cout<<"Target residual 4 momentum "<<TargetResidual4Momentum<<G4endl;
211 
212  Bsum+=( hitsT + hitsP);
213  Qsum+=(charged_hitsT + charged_hitsP);
214 
215  G4cout<<"Hitted # of nucleons of projectile and target "<<hitsP<<" "<<hitsT<<G4endl;
216  G4cout<<"Hitted # of protons of projectile and target "
217  <<charged_hitsP<<" "<<charged_hitsT<<G4endl<<G4endl;
218 G4cout<<"Bsum Qsum "<<Bsum<<" "<<Qsum<<G4endl;
219 #endif
220 
221 //=========================================================================================
222 // Fragment strings
223 
224  G4KineticTrackVector * theResult = 0;
225  G4double InvMass=SumStringMom.mag();
226  G4double SumMass(0.);
227  attempts = 0;
228  maxAttempts=100;
229 
230  #ifdef debug_PartonStringModel
231  G4LorentzVector SumPsecondr(0.,0.,0.,0.);
232  G4int QsumSec(0), BsumSec(0);
233  #endif
234 
235  do
236  {
237  attempts++;
238  if(theResult != 0)
239  {
240  std::for_each(theResult->begin(), theResult->end(), DeleteKineticTrack());
241  delete theResult;
242  }
243 
244  theResult = stringFragmentationModel->FragmentStrings(strings);
245 
246  if(attempts > maxAttempts ) break;
247 
248  #ifdef debug_PartonStringModel
249  G4cout<<"Attempt to fragment the strings "<<attempts<<G4endl;
250  G4cout<<"Parton-String model: Number of produced particles "<<theResult->size()<<G4endl;
251  SumPsecondr = G4LorentzVector(0.,0.,0.,0.);
252  QsumSec = 0; BsumSec = 0;
253  #endif
254 
255  SumMass=0.;
256 
257  for ( unsigned int i=0; i < theResult->size(); i++)
258  {
259  SumMass+=(*theResult)[i]->Get4Momentum().mag();
260  #ifdef debug_PartonStringModel
261  G4cout<<i<<" : "<<(*theResult)[i]->GetDefinition()->GetParticleName()<<" "
262  <<(*theResult)[i]->Get4Momentum()<<" "
263  <<(*theResult)[i]->Get4Momentum().mag()<<" "
264  <<(*theResult)[i]->GetDefinition()->GetPDGMass()<<G4endl;
265  SumPsecondr+=(*theResult)[i]->Get4Momentum();
266  BsumSec += (*theResult)[i]->GetDefinition()->GetBaryonNumber();
267  QsumSec += (*theResult)[i]->GetDefinition()->GetPDGCharge();
268  #endif
269  }
270 
271 #ifdef debug_PartonStringModel
272  if(Qsum != QsumSec) {
273  G4cout<<"Charge is not conserved!!! ----"<<G4endl;
274  G4cout<<" Qsum != QsumSec "<<Qsum<<" "<<QsumSec<<G4endl;
275  }
276  if(Bsum != BsumSec) {
277  G4cout<<"Baryon number is not conserved!!!"<<G4endl;
278  G4cout<<" Bsum != BsumSec "<<Bsum<<" "<<BsumSec<<G4endl;
279  }
280 #endif
281  } while(SumMass > InvMass); /* Loop checking, 07.08.2015, A.Ribon */
282 
283  std::for_each(strings->begin(), strings->end(), DeleteString() );
284  delete strings;
285 
286 #ifdef debug_PartonStringModel
287  G4cout<<G4endl<<"Baryon number balance "<<Bsum-BsumSec<<G4endl;
288  G4cout <<"Charge balance "<<Qsum-QsumSec<<G4endl;
289  G4cout <<"4 momentum balance "<<SumStringMom-SumPsecondr<<G4endl;
290  G4cout<<"End of string model work ------------"<<G4endl<<G4endl;
291 #endif
292 
293  return theResult;
294 }
295 
296 void G4VPartonStringModel::ModelDescription(std::ostream& outFile) const
297 {
298  outFile << GetModelName() << " has no description yet.\n";
299 }
300 
302 { return 0;}
303 
G4VPartonStringModel(const G4String &modelName="Parton String Model")
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
std::vector< G4ExcitedString * > G4ExcitedStringVector
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)=0
G4VPartonStringModel * theThis
void SetMomentum(G4LorentzVector &aMomentum)
Definition: G4Nucleon.hh:70
CLHEP::HepLorentzRotation G4LorentzRotation
virtual G4bool StartLoop()=0
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4V3DNucleus * GetProjectileNucleus() const
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
virtual void ModelDescription(std::ostream &outFile) const
G4GLOB_DLL std::ostream G4cout
G4VStringFragmentation * stringFragmentationModel
G4double GetMass() const
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
virtual G4String GetModelName() const
virtual void Init(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
virtual G4ExcitedStringVector * GetStrings()=0
#define G4endl
Definition: G4ios.hh:61
virtual G4Nucleon * GetNextNucleon()=0
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
CLHEP::HepLorentzVector G4LorentzVector