Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
28//
29//// ------------------------------------------------------------
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
41#include "G4ios.hh"
43
44#include "G4ParticleTable.hh"
45#include "G4IonTable.hh"
46
47
49 : G4VHighEnergyGenerator(modelName),
50 stringFragmentationModel(0),
51 theThis(0)
52{
53// Make shure Shotrylived partyicles are constructed.
54 G4ShortLivedConstructor ShortLived;
55 ShortLived.ConstructParticle();
56}
57
59{
60}
61
63 const G4DynamicParticle &aPrimary)
64{
65 G4ExcitedStringVector * strings = NULL;
66
67 G4DynamicParticle thePrimary=aPrimary;
68
70 G4LorentzVector Ptmp=thePrimary.Get4Momentum();
71 toZ.rotateZ(-1*Ptmp.phi());
72 toZ.rotateY(-1*Ptmp.theta());
73 thePrimary.Set4Momentum(toZ*Ptmp);
74 G4LorentzRotation toLab(toZ.inverse());
75
76 G4int attempts = 0, maxAttempts=20;
77 while ( strings == NULL )
78 {
79 if (attempts++ > maxAttempts )
80 {
81 throw G4HadronicException(__FILE__, __LINE__, "G4VPartonStringModel::Scatter(): fails to generate strings");
82 }
83
84 theThis->Init(theNucleus,thePrimary);
85
86 strings = GetStrings();
87 }
88
89 G4double stringEnergy(0);
90 G4LorentzVector SumStringMom(0.,0.,0.,0.);
91
92 for ( unsigned int astring=0; astring < strings->size(); astring++)
93 {
94// rotate string to lab frame, models have it aligned to z
95 stringEnergy += (*strings)[astring]->GetLeftParton()->Get4Momentum().t();
96 stringEnergy += (*strings)[astring]->GetRightParton()->Get4Momentum().t();
97 (*strings)[astring]->LorentzRotate(toLab);
98 SumStringMom+=(*strings)[astring]->Get4Momentum();
99 }
100
101 G4double InvMass=SumStringMom.mag();
102
103//#define debug_PartonStringModel
104#ifdef debug_PartonStringModel
105
106 G4V3DNucleus * fancynucleus=theThis->GetWoundedNucleus();
107
108 // loop over wounded nucleus
109 G4int hits(0), charged_hits(0);
110 G4ThreeVector hitNucleonMomentum(0.,0.,0.);
111 G4Nucleon * theCurrentNucleon = fancynucleus->StartLoop() ? fancynucleus->GetNextNucleon() : NULL;
112 while( theCurrentNucleon )
113 {
114 if(theCurrentNucleon->AreYouHit())
115 {
116 hits++;
117 hitNucleonMomentum += theCurrentNucleon->Get4Momentum().vect();
118 if ( theCurrentNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hits;
119 }
120 theCurrentNucleon = fancynucleus->GetNextNucleon();
121 }
122
123 G4int initialZ=fancynucleus->GetCharge();
124 G4int initialA=fancynucleus->GetMassNumber();
125 G4double initial_mass=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(initialZ,initialA);
126 G4double final_mass = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(initialZ-charged_hits, initialA-hits);
127 G4cout << "G4VPSM: strE, hit nucleons, Primary, SumStringE, nucleus intial, nucleus final, excitation estimate "
128 << stringEnergy << " " << hits << ", " << Ptmp.e() << ", "<< SumStringMom.e() << ", "
129 << initial_mass<< ", " << final_mass<< ", " << Ptmp.e() + initial_mass - final_mass - stringEnergy << G4endl;
130 G4cout << "momentum balance " << thePrimary.GetMomentum() + hitNucleonMomentum - SumStringMom.vect()<< G4endl;
131#endif
132
133// Fragment strings
134
135 G4KineticTrackVector * theResult = 0;
136 G4double SumMass(0.);
137 attempts = 0;
138 maxAttempts=100;
139 do
140 {
141 attempts++;
142 if(theResult != 0)
143 {
144 std::for_each(theResult->begin(), theResult->end(), DeleteKineticTrack());
145 delete theResult;
146 }
147 theResult = stringFragmentationModel->FragmentStrings(strings);
148 if(attempts > maxAttempts ) break;
149
150 //G4cout<<"G4endl<<"G4VPartonStringModel:: Final Result, Size "<<theResult->size()<<G4endl;
151
152 SumMass=0.;
153 //G4LorentzVector SumP(0.,0.,0.,0.);
154 for ( unsigned int i=0; i < theResult->size(); i++)
155 {
156 SumMass+=(*theResult)[i]->GetDefinition()->GetPDGMass();
157 //SumP+=(*theResult)[i]->Get4Momentum();
158 //G4cout<<i<<" : "<<(*theResult)[i]->GetDefinition()->GetParticleName();
159 //G4cout<<"p= " << (*theResult)[i]->Get4Momentum()<<" m= "<<(*theResult)[i]->Get4Momentum().mag()<<G4endl;
160 }
161
162 //G4cout<<"SumP "<<SumP<<G4endl;
163 } while(SumMass > InvMass);
164
165 std::for_each(strings->begin(), strings->end(), DeleteString() );
166 delete strings;
167
168 return theResult;
169}
170
171void G4VPartonStringModel::ModelDescription(std::ostream& outFile) const
172{
173 outFile << GetModelName() << " has no description yet.\n";
174}
std::vector< G4ExcitedString * > G4ExcitedStringVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
double theta() const
Hep3Vector vect() const
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4ThreeVector GetMomentum() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
virtual G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual G4int GetMassNumber()=0
virtual G4String GetModelName() const
G4KineticTrackVector * Scatter(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
virtual G4V3DNucleus * GetWoundedNucleus() const =0
G4VPartonStringModel(const G4String &modelName="Parton String Model")
virtual void Init(const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
virtual G4ExcitedStringVector * GetStrings()=0
virtual void ModelDescription(std::ostream &outFile) const
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)=0