Geant4  10.00.p01
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 67999 2013-03-13 11:14: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 
47 //#define debug_PartonStringModel
48 
50  : G4VHighEnergyGenerator(modelName),
51  stringFragmentationModel(0),
52  theThis(0)
53 {
54 // Make shure Shotrylived partyicles are constructed.
55  G4ShortLivedConstructor ShortLived;
56  ShortLived.ConstructParticle();
57 }
58 
60 {
61 }
62 
64  const G4DynamicParticle &aPrimary)
65 {
66  G4ExcitedStringVector * strings = NULL;
67 
68  G4DynamicParticle thePrimary=aPrimary;
69 
70 #ifdef debug_PartonStringModel
71  G4cout<<"Parton-String model is runnung ------------"<<G4endl;
72  G4cout<<"Projectile Name Mass "<<thePrimary.GetDefinition()->GetParticleName()<<" "
73  <<thePrimary.GetMass()<<G4endl;
74  G4cout<<" Momentum "<<thePrimary.Get4Momentum()<<G4endl;
75  G4cout<<"Target nucleus A Z "<<theNucleus.GetA_asInt()<<" "
76  <<theNucleus.GetZ_asInt()<<G4endl<<G4endl;
77 #endif
78 
80  G4LorentzVector Ptmp=thePrimary.Get4Momentum();
81  toZ.rotateZ(-1*Ptmp.phi());
82  toZ.rotateY(-1*Ptmp.theta());
83  thePrimary.Set4Momentum(toZ*Ptmp);
84  G4LorentzRotation toLab(toZ.inverse());
85 
86  G4int attempts = 0, maxAttempts=20;
87  while ( strings == NULL )
88  {
89  if (attempts++ > maxAttempts )
90  {
91  throw G4HadronicException(__FILE__, __LINE__,
92  "G4VPartonStringModel::Scatter(): fails to generate strings");
93  }
94  theThis->Init(theNucleus,thePrimary);
95 
96  strings = GetStrings();
97  }
98 
99  G4double stringEnergy(0);
100  G4LorentzVector SumStringMom(0.,0.,0.,0.);
101 
102 #ifdef debug_PartonStringModel
103  G4cout<<"Parton-String model: Number of produced strings "<<strings->size()<<G4endl;
104 #endif
105 
106  for ( unsigned int astring=0; astring < strings->size(); astring++)
107  {
108 // rotate string to lab frame, models have it aligned to z
109  stringEnergy += (*strings)[astring]->GetLeftParton()->Get4Momentum().t();
110  stringEnergy += (*strings)[astring]->GetRightParton()->Get4Momentum().t();
111  (*strings)[astring]->LorentzRotate(toLab);
112  SumStringMom+=(*strings)[astring]->Get4Momentum();
113  }
114 
115  G4double InvMass=SumStringMom.mag();
116 
117 #ifdef debug_PartonStringModel
118  G4cout<<"Parton-String model: SumStringMom "<<SumStringMom<<G4endl;
119 
120  G4V3DNucleus * fancynucleus=theThis->GetWoundedNucleus();
121 
122  // loop over wounded nucleus
123  G4int hits(0), charged_hits(0);
124  G4ThreeVector hitNucleonMomentum(0.,0.,0.);
125  G4Nucleon * theCurrentNucleon = fancynucleus->StartLoop() ?
126  fancynucleus->GetNextNucleon() : NULL;
127  while( theCurrentNucleon )
128  {
129  if(theCurrentNucleon->AreYouHit())
130  {
131  hits++;
132  hitNucleonMomentum += theCurrentNucleon->Get4Momentum().vect();
133  if ( theCurrentNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hits;
134  }
135  theCurrentNucleon = fancynucleus->GetNextNucleon();
136  }
137 
138  G4int initialZ=fancynucleus->GetCharge();
139  G4int initialA=fancynucleus->GetMassNumber();
140  G4double initial_mass=
142  G4double final_mass =
144  initialZ-charged_hits, initialA-hits);
145  G4cout << "G4VPSM: " <<G4endl
146  << "strE "<<stringEnergy <<G4endl
147  << "hit nucleons "<<hits <<G4endl
148  << "Primary "<<Ptmp.e() <<G4endl
149  << "SumStringE "<<SumStringMom.e() <<G4endl
150  << "Nucleus intial "<<initial_mass <<G4endl
151  << "Nucleus final "<<final_mass <<G4endl
152  << "Excitation estimate "
153  <<Ptmp.e() + initial_mass - final_mass - stringEnergy << G4endl;
154  G4cout << "momentum balance "
155  << thePrimary.GetMomentum() + hitNucleonMomentum -
156  SumStringMom.vect()<< G4endl;
157 #endif
158 
159 // Fragment strings
160 
161  G4KineticTrackVector * theResult = 0;
162  G4double SumMass(0.);
163  attempts = 0;
164  maxAttempts=100;
165  do
166  {
167  attempts++;
168  if(theResult != 0)
169  {
170  std::for_each(theResult->begin(), theResult->end(), DeleteKineticTrack());
171  delete theResult;
172  }
173  theResult = stringFragmentationModel->FragmentStrings(strings);
174  if(attempts > maxAttempts ) break;
175 
176 #ifdef debug_PartonStringModel
177  G4cout<<"Parton-String model: Number of produced particles "<<theResult->size()<<G4endl;
178  G4LorentzVector SumPsecondr(0.,0.,0.,0.);
179 #endif
180 
181  SumMass=0.;
182 
183  for ( unsigned int i=0; i < theResult->size(); i++)
184  {
185  SumMass+=(*theResult)[i]->GetDefinition()->GetPDGMass();
186  //SumP+=(*theResult)[i]->Get4Momentum();
187 #ifdef debug_PartonStringModel
188  G4cout<<i<<" : "<<(*theResult)[i]->GetDefinition()->GetParticleName()<<" "
189  <<(*theResult)[i]->Get4Momentum()<<" "
190  <<(*theResult)[i]->Get4Momentum().mag()<<G4endl;
191  SumPsecondr+=(*theResult)[i]->Get4Momentum();
192 #endif
193  }
194 #ifdef debug_PartonStringModel
195  G4cout<<"SumP secondaries "<<SumPsecondr<<G4endl;
196 #endif
197  } while(SumMass > InvMass);
198 
199  std::for_each(strings->begin(), strings->end(), DeleteString() );
200  delete strings;
201 
202 #ifdef debug_PartonStringModel
203  G4cout<<"End of string model work ------------"<<G4endl<<G4endl;
204 #endif
205  return theResult;
206 }
207 
209 {
210  outFile << GetModelName() << " has no description yet.\n";
211 }
212 
214 { return 0;}
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
virtual G4int GetCharge()=0
CLHEP::Hep3Vector G4ThreeVector
std::ofstream outFile
Definition: GammaRayTel.cc:68
CLHEP::HepLorentzRotation G4LorentzRotation
virtual G4bool StartLoop()=0
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4int GetMassNumber()=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 void ModelDescription(std::ostream &outFile) const
G4IonTable * GetIonTable() const
G4GLOB_DLL std::ostream G4cout
G4VStringFragmentation * stringFragmentationModel
G4double GetMass() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1232
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
static G4ParticleTable * GetParticleTable()
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 G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
virtual G4Nucleon * GetNextNucleon()=0
double G4double
Definition: G4Types.hh:76
G4ThreeVector GetMomentum() const
CLHEP::HepLorentzVector G4LorentzVector