Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VPartonStringModel Class Referenceabstract

#include <G4VPartonStringModel.hh>

Inheritance diagram for G4VPartonStringModel:
Collaboration diagram for G4VPartonStringModel:

Public Member Functions

 G4VPartonStringModel (const G4String &modelName="Parton String Model")
 
virtual ~G4VPartonStringModel ()
 
void SetFragmentationModel (G4VStringFragmentation *aModel)
 
G4KineticTrackVectorScatter (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
virtual G4V3DNucleusGetWoundedNucleus () const =0
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual G4V3DNucleusGetProjectileNucleus () const
 
- Public Member Functions inherited from G4VHighEnergyGenerator
 G4VHighEnergyGenerator (const G4String &modelName="High Energy Generator")
 
virtual ~G4VHighEnergyGenerator ()
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double AbsoluteLevel)
 
virtual G4String GetModelName () const
 

Protected Member Functions

virtual void Init (const G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)=0
 
virtual G4ExcitedStringVectorGetStrings ()=0
 
void SetThisPointer (G4VPartonStringModel *aPointer)
 
G4bool EnergyAndMomentumCorrector (G4KineticTrackVector *Output, G4LorentzVector &TotalCollisionMomentum)
 

Detailed Description

Definition at line 51 of file G4VPartonStringModel.hh.

Constructor & Destructor Documentation

G4VPartonStringModel::G4VPartonStringModel ( const G4String modelName = "Parton String Model")

Definition at line 47 of file G4VPartonStringModel.cc.

48 : G4VHighEnergyGenerator(modelName), stringFragmentationModel(0), theThis(0)
49 {
50  // Make sure that Shortlived particles are constructed.
51  G4ShortLivedConstructor ShortLived;
52  ShortLived.ConstructParticle();
53 }
G4VHighEnergyGenerator(const G4String &modelName="High Energy Generator")

Here is the call graph for this function:

G4VPartonStringModel::~G4VPartonStringModel ( )
virtual

Definition at line 56 of file G4VPartonStringModel.cc.

57 {
58 }

Member Function Documentation

G4bool G4VPartonStringModel::EnergyAndMomentumCorrector ( G4KineticTrackVector Output,
G4LorentzVector TotalCollisionMomentum 
)
protected
G4V3DNucleus * G4VPartonStringModel::GetProjectileNucleus ( ) const
virtual

Reimplemented from G4VHighEnergyGenerator.

Reimplemented in G4FTFModel, G4QGSModel< ParticipantType >, G4QGSModel< G4GammaParticipants >, and G4QGSModel< G4QGSParticipants >.

Definition at line 322 of file G4VPartonStringModel.cc.

323 {
324  return 0;
325 }

Here is the caller graph for this function:

virtual G4ExcitedStringVector* G4VPartonStringModel::GetStrings ( )
protectedpure virtual

Implemented in G4FTFModel, G4QGSModel< ParticipantType >, G4QGSModel< G4GammaParticipants >, and G4QGSModel< G4QGSParticipants >.

Here is the caller graph for this function:

virtual G4V3DNucleus* G4VPartonStringModel::GetWoundedNucleus ( ) const
pure virtual

Implements G4VHighEnergyGenerator.

Implemented in G4FTFModel, G4QGSModel< ParticipantType >, G4QGSModel< G4GammaParticipants >, and G4QGSModel< G4QGSParticipants >.

Here is the caller graph for this function:

virtual void G4VPartonStringModel::Init ( const G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)
protectedpure virtual

Implemented in G4FTFModel, G4QGSModel< ParticipantType >, G4QGSModel< G4GammaParticipants >, and G4QGSModel< G4QGSParticipants >.

Here is the caller graph for this function:

void G4VPartonStringModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VHighEnergyGenerator.

Reimplemented in G4FTFModel, G4QGSModel< ParticipantType >, G4QGSModel< G4GammaParticipants >, and G4QGSModel< G4QGSParticipants >.

Definition at line 316 of file G4VPartonStringModel.cc.

317 {
318  outFile << GetModelName() << " has no description yet.\n";
319 }
virtual G4String GetModelName() const

Here is the call graph for this function:

G4KineticTrackVector * G4VPartonStringModel::Scatter ( const G4Nucleus theNucleus,
const G4DynamicParticle thePrimary 
)
virtual

Implements G4VHighEnergyGenerator.

Definition at line 61 of file G4VPartonStringModel.cc.

63 {
64  G4ExcitedStringVector * strings = NULL;
65  G4DynamicParticle thePrimary = aPrimary;
66  G4LorentzVector SumStringMom( 0.0, 0.0, 0.0, 0.0 );
67  G4KineticTrackVector * theResult = 0;
68 
69  #ifdef debug_PartonStringModel
70  G4cout << G4endl;
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  G4int Bsum = thePrimary.GetDefinition()->GetBaryonNumber() + theNucleus.GetA_asInt();
78  G4int Qsum = thePrimary.GetDefinition()->GetPDGCharge() + theNucleus.GetZ_asInt();
79  G4cout << "Initial baryon number " << Bsum << G4endl;
80  G4cout << "Initial charge " << Qsum << G4endl;
81  G4cout << "-------------- Parton-String model: Generation of strings -------" << G4endl
82  << G4endl;
83  Bsum -= theNucleus.GetA_asInt(); Qsum -= theNucleus.GetZ_asInt();
84  if ( theThis->GetProjectileNucleus() ) {
85  Bsum -= thePrimary.GetDefinition()->GetBaryonNumber();
86  Qsum -= thePrimary.GetDefinition()->GetPDGCharge();
87  }
88  G4int QsumSec( 0 ), BsumSec( 0 );
89  G4LorentzVector SumPsecondr( 0.0, 0.0, 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 ) {
105  ed << " Projectile Name Mass " << thePrimary.GetDefinition()->GetParticleName()
106  << " " << thePrimary.GetMass()<< G4endl;
107  ed << " Momentum " << thePrimary.Get4Momentum() << G4endl;
108  ed << " Target nucleus A Z " << theNucleus.GetA_asInt() << " "
109  << theNucleus.GetZ_asInt() << G4endl;
110  G4Exception( "G4VPartonStringModel::Scatter(): fails to generate or fragment strings ",
111  "HAD_PARTON_STRING_001", JustWarning, ed );
112  }
113 
114  Success = true;
115 
116  theThis->Init( theNucleus, thePrimary );
117  strings = GetStrings();
118 
119  if ( strings == 0 ) {
120  Success = false;
121  continue;
122  }
123 
124  G4double stringEnergy( 0.0 );
125  SumStringMom = G4LorentzVector( 0.0, 0.0, 0.0, 0.0 );
126 
127  #ifdef debug_PartonStringModel
128  G4cout << "------------ Parton-String model: Number of produced strings ---- " << strings->size() << G4endl;
129  #endif
130 
131  for ( unsigned int astring = 0; astring < strings->size(); astring++ ) {
132  // rotate string to lab frame, models have it aligned to z
133  if ( (*strings)[astring]->IsExcited() ) {
134  stringEnergy += (*strings)[astring]->GetLeftParton()->Get4Momentum().t();
135  stringEnergy += (*strings)[astring]->GetRightParton()->Get4Momentum().t();
136  (*strings)[astring]->LorentzRotate( toLab );
137  SumStringMom += (*strings)[astring]->Get4Momentum();
138 
139  #ifdef debug_PartonStringModel
140  G4cout << "String No " << astring + 1 << " " << (*strings)[astring]->Get4Momentum()
141  << " " << (*strings)[astring]->Get4Momentum().mag()
142  << " Partons " << (*strings)[astring]->GetLeftParton()->GetDefinition()->GetPDGEncoding()
143  << " " << (*strings)[astring]->GetRightParton()->GetDefinition()->GetPDGEncoding()
144  <<G4endl;
145  #endif
146 
147  } else {
148  stringEnergy += (*strings)[astring]->GetKineticTrack()->Get4Momentum().t();
149  (*strings)[astring]->LorentzRotate( toLab );
150  SumStringMom += (*strings)[astring]->GetKineticTrack()->Get4Momentum();
151 
152  #ifdef debug_PartonStringModel
153  G4cout << "A track No " << astring + 1 << " "
154  << (*strings)[astring]->GetKineticTrack()->Get4Momentum() << " "
155  << (*strings)[astring]->GetKineticTrack()->Get4Momentum().mag() << " "
156  << (*strings)[astring]->GetKineticTrack()->GetDefinition()->GetParticleName()
157  << G4endl;
158  #endif
159 
160  }
161  }
162 
163  #ifdef debug_PartonStringModel
164  G4cout << G4endl << "SumString4Mom " << SumStringMom << G4endl;
165  G4LorentzVector TargetResidual4Momentum( 0.0, 0.0, 0.0, 0.0 );
166  G4LorentzVector ProjectileResidual4Momentum( 0.0, 0.0, 0.0, 0.0 );
167  G4int hitsT( 0 ), charged_hitsT( 0 );
168  G4int hitsP( 0 ), charged_hitsP( 0 );
169  G4double ExcitationEt( 0.0 ), ExcitationEp( 0.0 );
170  #endif
171 
172  G4V3DNucleus * ProjResNucleus = theThis->GetProjectileNucleus();
173  G4Nucleon * theNuclNucleon( 0 );
174 
175  if ( ProjResNucleus != 0 ) {
176  theNuclNucleon = ProjResNucleus->StartLoop() ? ProjResNucleus->GetNextNucleon() : NULL;
177  while ( theNuclNucleon ) { /* Loop checking, 07.08.2015, A.Ribon */
178  if ( theNuclNucleon->AreYouHit() ) {
179  G4LorentzVector tmp = toLab * theNuclNucleon->Get4Momentum();
180 
181  #ifdef debug_PartonStringModel
182  ProjectileResidual4Momentum += tmp;
183  hitsP++;
184  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsP;
185  ExcitationEp += theNuclNucleon->GetBindingEnergy();
186  #endif
187 
188  theNuclNucleon->SetMomentum( tmp );
189  }
190  theNuclNucleon = ProjResNucleus->GetNextNucleon();
191  }
192 
193  #ifdef debug_PartonStringModel
194  G4cout << "Projectile residual A, Z and E* "
195  << thePrimary.GetDefinition()->GetBaryonNumber() - hitsP << " "
196  << thePrimary.GetDefinition()->GetPDGCharge() - charged_hitsP << " "
197  << ExcitationEp << G4endl;
198  G4cout << "Projectile residual 4 momentum " << ProjectileResidual4Momentum << G4endl;
199  #endif
200  }
201 
202  G4V3DNucleus * ResNucleus = theThis->GetWoundedNucleus();
203 
204  // loop over wounded nucleus
205  theNuclNucleon = ResNucleus->StartLoop() ? ResNucleus->GetNextNucleon() : NULL;
206  while ( theNuclNucleon ) { /* Loop checking, 07.08.2015, A.Ribon */
207  if ( theNuclNucleon->AreYouHit() ) {
208  G4LorentzVector tmp = toLab * theNuclNucleon->Get4Momentum();
209 
210  #ifdef debug_PartonStringModel
211  TargetResidual4Momentum += tmp;
212  hitsT++;
213  if ( theNuclNucleon->GetDefinition() == G4Proton::Proton() ) ++charged_hitsT;
214  ExcitationEt += theNuclNucleon->GetBindingEnergy();
215  #endif
216 
217  theNuclNucleon->SetMomentum( tmp );
218  }
219  theNuclNucleon = ResNucleus->GetNextNucleon();
220  }
221 
222  #ifdef debug_PartonStringModel
223  G4cout << "Target residual A, Z and E* "
224  << theNucleus.GetA_asInt() - hitsT << " "
225  << theNucleus.GetZ_asInt() - charged_hitsT << " "
226  << ExcitationEt << G4endl;
227  G4cout << "Target residual 4 momentum " << TargetResidual4Momentum << G4endl;
228  Bsum += ( hitsT + hitsP );
229  Qsum += ( charged_hitsT + charged_hitsP );
230  G4cout << "Hitted # of nucleons of projectile and target "
231  << hitsP << " " << hitsT << G4endl;
232  G4cout << "Hitted # of protons of projectile and target "
233  << charged_hitsP << " " << charged_hitsT << G4endl << G4endl;
234  G4cout << "Bsum Qsum " << Bsum << " " << Qsum << G4endl << G4endl;
235  #endif
236 
237  //--- Fragment strings ---
238 
239  #ifdef debug_PartonStringModel
240  G4cout << "------Parton-String model: Number of produced strings ----------- " << strings->size() << G4endl;
241  #endif
242 
243  G4double InvMass = SumStringMom.mag();
244  G4double SumMass( 0.0 );
245 
246  #ifdef debug_PartonStringModel
247  QsumSec = 0; BsumSec = 0;
248  SumPsecondr = G4LorentzVector( 0.0, 0.0, 0.0, 0.0 );
249  #endif
250 
251  if ( theResult != 0 ) {
252  std::for_each( theResult->begin(), theResult->end(), DeleteKineticTrack() );
253  delete theResult;
254  }
255 
256  theResult = stringFragmentationModel->FragmentStrings( strings );
257 
258  if ( theResult == 0 ) {
259  Success = false;
260  continue;
261  }
262 
263  #ifdef debug_PartonStringModel
264  G4cout << "Attempt to fragment the strings " << attempts << G4endl;
265  G4cout << "Parton-String model: Number of produced particles " << theResult->size() << G4endl;
266  SumPsecondr = G4LorentzVector( 0.0, 0.0, 0.0, 0.0 );
267  QsumSec = 0; BsumSec = 0;
268  #endif
269 
270  SumMass = 0.0;
271 
272  for ( unsigned int i = 0; i < theResult->size(); i++ ) {
273  SumMass += (*theResult)[i]->Get4Momentum().mag();
274 
275  #ifdef debug_PartonStringModel
276  G4cout << i << " : " << (*theResult)[i]->GetDefinition()->GetParticleName() << " "
277  << (*theResult)[i]->Get4Momentum() << " "
278  << (*theResult)[i]->Get4Momentum().mag() << " "
279  << (*theResult)[i]->GetDefinition()->GetPDGMass() << G4endl;
280  SumPsecondr += (*theResult)[i]->Get4Momentum();
281  BsumSec += (*theResult)[i]->GetDefinition()->GetBaryonNumber();
282  QsumSec += (*theResult)[i]->GetDefinition()->GetPDGCharge();
283  #endif
284  }
285 
286  #ifdef debug_PartonStringModel
287  G4cout << G4endl << "-----------------------Parton-String model: balances -------------" << G4endl;
288  if ( Qsum != QsumSec ) {
289  G4cout << "Charge is not conserved!!! ----" << G4endl;
290  G4cout << " Qsum != QsumSec " << Qsum << " " << QsumSec << G4endl;
291  }
292  if ( Bsum != BsumSec ) {
293  G4cout << "Baryon number is not conserved!!!" << G4endl;
294  G4cout << " Bsum != BsumSec " << Bsum << " " << BsumSec << G4endl;
295  }
296  #endif
297 
298  if ( SumMass > InvMass || SumMass == 0.0 ) Success = false;
299 
300  std::for_each( strings->begin(), strings->end(), DeleteString() );
301  delete strings;
302 
303  } while ( ! Success ); /* Loop checking, 07.08.2015, A.Ribon */
304 
305  #ifdef debug_PartonStringModel
306  G4cout << "Baryon number balance " << Bsum - BsumSec << G4endl;
307  G4cout << "Charge balance " << Qsum - QsumSec << G4endl;
308  G4cout << "4 momentum balance " << SumStringMom - SumPsecondr << G4endl;
309  G4cout << "---------------------End of Parton-String model work -------------" << G4endl << G4endl;
310  #endif
311 
312  return theResult;
313 }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
std::vector< G4ExcitedString * > G4ExcitedStringVector
virtual G4KineticTrackVector * FragmentStrings(const G4ExcitedStringVector *theStrings)=0
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
virtual G4bool StartLoop()=0
virtual G4V3DNucleus * GetWoundedNucleus() const =0
virtual G4V3DNucleus * GetProjectileNucleus() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
HepLorentzRotation & rotateY(double delta)
double phi() const
G4GLOB_DLL std::ostream G4cout
double theta() const
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4LorentzVector Get4Momentum() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
HepLorentzRotation & rotateZ(double delta)
void Set4Momentum(const G4LorentzVector &momentum)
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
HepLorentzRotation inverse() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:

void G4VPartonStringModel::SetFragmentationModel ( G4VStringFragmentation aModel)
inline

Definition at line 82 of file G4VPartonStringModel.hh.

83 {
84  stringFragmentationModel = aModel;
85 }

Here is the caller graph for this function:

void G4VPartonStringModel::SetThisPointer ( G4VPartonStringModel aPointer)
inlineprotected

Definition at line 87 of file G4VPartonStringModel.hh.

88 {
89  theThis=aPointer;
90 }

Here is the caller graph for this function:


The documentation for this class was generated from the following files: