Geant4  10.02.p03
G4LMsdGenerator Class Reference

#include <G4LMsdGenerator.hh>

Inheritance diagram for G4LMsdGenerator:
Collaboration diagram for G4LMsdGenerator:

Public Member Functions

 G4LMsdGenerator (const G4String &name="LMsdGenerator")
 
 ~G4LMsdGenerator ()
 
G4bool IsApplicable (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
G4double SampleMx (const G4HadProjectile *aParticle)
 
G4double SampleT (const G4HadProjectile *aParticle, G4double Mx)
 
void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

 G4LMsdGenerator (const G4LMsdGenerator &right)
 
const G4LMsdGeneratoroperator= (const G4LMsdGenerator &right)
 
int operator== (const G4LMsdGenerator &right) const
 
int operator!= (const G4LMsdGenerator &right) const
 

Private Attributes

G4int fPDGencoding
 

Static Private Attributes

static const G4double fMxBdata [23][2]
 
static const G4double fProbMx [60][2]
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 60 of file G4LMsdGenerator.hh.

Constructor & Destructor Documentation

◆ G4LMsdGenerator() [1/2]

G4LMsdGenerator::G4LMsdGenerator ( const G4String name = "LMsdGenerator")

Definition at line 46 of file G4LMsdGenerator.cc.

47  : G4HadronicInteraction(name)
48 
49 {
50  fPDGencoding = 0;
51 
52  // theParticleChange = new G4HadFinalState;
53 }
G4HadronicInteraction(const G4String &modelName="HadronicModel")

◆ ~G4LMsdGenerator()

G4LMsdGenerator::~G4LMsdGenerator ( )

Definition at line 55 of file G4LMsdGenerator.cc.

56 {
57  // delete theParticleChange;
58 }

◆ G4LMsdGenerator() [2/2]

G4LMsdGenerator::G4LMsdGenerator ( const G4LMsdGenerator right)
private

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4LMsdGenerator::ApplyYourself ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 111 of file G4LMsdGenerator.cc.

113 {
115 
116  const G4HadProjectile* aParticle = &aTrack;
117  G4double eTkin = aParticle->GetKineticEnergy();
118 
119  if( eTkin <= 1.*CLHEP::GeV )
120  {
122  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
123  return &theParticleChange;
124  }
125 
126  G4int A = targetNucleus.GetA_asInt();
127  G4int Z = targetNucleus.GetZ_asInt();
128 
129  G4double plab = aParticle->GetTotalMomentum();
130  G4double plab2 = plab*plab;
131 
132  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
133  G4double partMass = theParticle->GetPDGMass();
134 
135  G4double oldE = partMass + eTkin;
136 
138  G4double targMass2 = targMass*targMass;
139 
140  G4LorentzVector partLV = aParticle->Get4Momentum();
141 
142  G4double sumE = oldE + targMass;
143  G4double sumE2 = sumE*sumE;
144 
145  G4ThreeVector p1 = partLV.vect();
146  // G4cout<<"p1 = "<<p1<<G4endl;
147  G4ParticleMomentum p1unit = p1.unit();
148 
149  G4double Mx = SampleMx(aParticle); // in GeV
150  G4double t = 0.; // SampleT( aParticle, Mx); // in GeV
151 
152  Mx *= CLHEP::GeV;
153 
154  G4double Mx2 = Mx*Mx;
155 
156  // equation for q|| based on sum-E-P and new invariant mass
157 
158  G4double B = sumE2 + targMass2 - Mx2 - plab2;
159 
160  G4double a = 4*(plab2 - sumE2);
161  G4double b = 4*plab*B;
162  G4double c = B*B - 4*sumE2*targMass2;
163  G4double det2 = b*b - 4*a*c;
164  G4double qLong, det, eRetard; // , x2, x3, e2;
165 
166  if( det2 >= 0.)
167  {
168  det = std::sqrt(det2);
169  qLong = (-b - det)/2./a;
170  eRetard = std::sqrt((plab-qLong)*(plab-qLong)+Mx2);
171  }
172  else
173  {
175  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
176  return &theParticleChange;
177  }
179 
180  plab -= qLong;
181 
182  G4ThreeVector pRetard = plab*p1unit;
183 
184  G4ThreeVector pTarg = p1 - pRetard;
185 
186  G4double eTarg = std::sqrt( targMass2 + pTarg.mag2()); // std::sqrt( targMass*targMass + pTarg.mag2() );
187 
188  G4LorentzVector lvRetard(pRetard, eRetard);
189  G4LorentzVector lvTarg(pTarg, eTarg);
190 
191  lvTarg += lvRetard; // sum LV
192 
193  G4ThreeVector bst = lvTarg.boostVector();
194 
195  lvRetard.boost(-bst); // to CNS
196 
197  G4ThreeVector pCMS = lvRetard.vect();
198  G4double momentumCMS = pCMS.mag();
199  G4double tMax = 4.0*momentumCMS*momentumCMS;
200 
201  if( t > tMax ) t = tMax*G4UniformRand();
202 
203  G4double cost = 1. - 2.0*t/tMax;
204 
205 
207  G4double sint;
208 
209  if( cost > 1.0 || cost < -1.0 ) //
210  {
211  cost = 1.0;
212  sint = 0.0;
213  }
214  else // normal situation
215  {
216  sint = std::sqrt( (1.0-cost)*(1.0+cost) );
217  }
218  G4ThreeVector v1( sint*std::cos(phi), sint*std::sin(phi), cost);
219 
220  v1 *= momentumCMS;
221 
222  G4LorentzVector lvRes( v1.x(),v1.y(),v1.z(), std::sqrt( momentumCMS*momentumCMS + Mx2));
223 
224  lvRes.boost(bst); // to LS
225 
226  lvTarg -= lvRes;
227 
228  G4double eRecoil = lvTarg.e() - targMass;
229 
230  if( eRecoil > 100.*CLHEP::MeV ) // add recoil nucleus
231  {
232  G4ParticleDefinition * recoilDef = 0;
233 
234  if ( Z == 1 && A == 1 ) { recoilDef = G4Proton::Proton(); }
235  else if ( Z == 1 && A == 2 ) { recoilDef = G4Deuteron::Deuteron(); }
236  else if ( Z == 1 && A == 3 ) { recoilDef = G4Triton::Triton(); }
237  else if ( Z == 2 && A == 3 ) { recoilDef = G4He3::He3(); }
238  else if ( Z == 2 && A == 4 ) { recoilDef = G4Alpha::Alpha(); }
239  else
240  {
241  recoilDef =
243  }
244  G4DynamicParticle * aSec = new G4DynamicParticle( recoilDef, lvTarg);
246  }
247  else if( eRecoil > 0.0 )
248  {
250  }
251 
253  FindParticle(fPDGencoding);
254 
255  // G4cout<<fPDGencoding<<", "<<ddPart->GetParticleName()<<", "<<ddPart->GetPDGMass()<<" MeV; lvRes = "<<lvRes<<G4endl;
256 
257  G4DynamicParticle * aRes = new G4DynamicParticle( ddPart, lvRes);
258  theParticleChange.AddSecondary(aRes); // simply return resonance
259 
260 
261  /*
262  // Recursive decay using methods of G4KineticTrack
263 
264  G4KineticTrack ddkt( ddPart, 0., G4ThreeVector(0.,0.,0.), lvRes);
265  G4KineticTrackVector* ddktv = ddkt.Decay();
266 
267  G4DecayKineticTracks decay( ddktv );
268 
269  for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
270  {
271  G4DynamicParticle * aNew =
272  new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
273  ddktv->operator[](i)->Get4Momentum());
274  // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
275 
276  theParticleChange.AddSecondary(aNew);
277  delete ddktv->operator[](i);
278  }
279  delete ddktv;
280  */
281  return &theParticleChange;
282 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4LorentzVector & Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
int G4int
Definition: G4Types.hh:78
static const double GeV
double mag2() const
void SetStatusChange(G4HadFinalStateStatus aS)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
G4IonTable * GetIonTable() const
Float_t Z
G4double SampleMx(const G4HadProjectile *aParticle)
double mag() const
HepLorentzVector & boost(double, double, double)
Hep3Vector unit() const
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double GetTotalMomentum() const
const G4ParticleDefinition * GetDefinition() const
void SetEnergyChange(G4double anEnergy)
G4double GetKineticEnergy() const
static G4ParticleTable * GetParticleTable()
void SetLocalEnergyDeposit(G4double aE)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
static const double twopi
Definition: SystemOfUnits.h:54
static G4He3 * He3()
Definition: G4He3.cc:94
void SetMomentumChange(const G4ThreeVector &aV)
static const double MeV
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4LMsdGenerator::IsApplicable ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 77 of file G4LMsdGenerator.cc.

79 {
80  G4bool applied = false;
81 
82  if( ( aTrack.GetDefinition() == G4Proton::Proton() ||
83  aTrack.GetDefinition() == G4Neutron::Neutron() ) &&
84  targetNucleus.GetA_asInt() >= 1 &&
85  aTrack.GetKineticEnergy() > 1800*CLHEP::MeV ) // 750*CLHEP::MeV )
86  {
87  applied = true;
88  }
89  else if( ( aTrack.GetDefinition() == G4PionPlus::PionPlus() ||
90  aTrack.GetDefinition() == G4PionMinus::PionMinus() ) &&
91  targetNucleus.GetA_asInt() >= 1 &&
92  aTrack.GetKineticEnergy() > 2340*CLHEP::MeV )
93  {
94  applied = true;
95  }
96  else if( ( aTrack.GetDefinition() == G4KaonPlus::KaonPlus() ||
97  aTrack.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
98  targetNucleus.GetA_asInt() >= 1 &&
99  aTrack.GetKineticEnergy() > 1980*CLHEP::MeV )
100  {
101  applied = true;
102  }
103  return applied;
104 }
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
bool G4bool
Definition: G4Types.hh:79
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static const double MeV
Here is the call graph for this function:

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 60 of file G4LMsdGenerator.cc.

61 {
62  outFile << GetModelName() <<" consists of a "
63  << " string model and a stage to de-excite the excited nuclear fragment."
64  << "\n<p>"
65  << "The string model simulates the interaction of\n"
66  << "an incident hadron with a nucleus, forming \n"
67  << "excited strings, decays these strings into hadrons,\n"
68  << "and leaves an excited nucleus. \n"
69  << "<p>The string model:\n";
70 }
const G4String & GetModelName() const
Here is the call graph for this function:

◆ operator!=()

int G4LMsdGenerator::operator!= ( const G4LMsdGenerator right) const
private

◆ operator=()

const G4LMsdGenerator& G4LMsdGenerator::operator= ( const G4LMsdGenerator right)
private

◆ operator==()

int G4LMsdGenerator::operator== ( const G4LMsdGenerator right) const
private

◆ SampleMx()

G4double G4LMsdGenerator::SampleMx ( const G4HadProjectile aParticle)

Definition at line 288 of file G4LMsdGenerator.cc.

289 {
290  G4double Mx = 0.;
291  G4int i;
292  G4double rand = G4UniformRand();
293 
294  for( i = 0; i < 60; i++)
295  {
296  if( rand >= fProbMx[i][1] ) break;
297  }
298  if(i <= 0) Mx = fProbMx[0][0];
299  else if(i >= 59) Mx = fProbMx[59][0];
300  else Mx = fProbMx[i][0];
301 
302  fPDGencoding = 0;
303 
304  if ( Mx <= 1.45 )
305  {
306  if( aParticle->GetDefinition() == G4Proton::Proton() )
307  {
308  Mx = 1.44;
309  fPDGencoding = 12212;
310  }
311  else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
312  {
313  Mx = 1.44;
314  fPDGencoding = 12112;
315  }
316  else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
317  {
318  // Mx = 1.3;
319  // fPDGencoding = 100211;
320  Mx = 1.26;
321  fPDGencoding = 20213; // a1(1260)+
322  }
323  else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
324  {
325  // Mx = 1.3;
326  // fPDGencoding = -100211;
327  Mx = 1.26;
328  fPDGencoding = -20213; // a1(1260)-
329  }
330  else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
331  {
332  Mx = 1.27;
333  fPDGencoding = 10323;
334  }
335  else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
336  {
337  Mx = 1.27;
338  fPDGencoding = -10323;
339  }
340  }
341  else if ( Mx <= 1.55 )
342  {
343  if( aParticle->GetDefinition() == G4Proton::Proton() )
344  {
345  Mx = 1.52;
346  fPDGencoding = 2124;
347  }
348  else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
349  {
350  Mx = 1.52;
351  fPDGencoding = 1214;
352  }
353  else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
354  {
355  // Mx = 1.45;
356  // fPDGencoding = 10211;
357  Mx = 1.32;
358  fPDGencoding = 215; // a2(1320)+
359  }
360  else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
361  {
362  // Mx = 1.45;
363  // fPDGencoding = -10211;
364  Mx = 1.32;
365  fPDGencoding = -215; // a2(1320)-
366  }
367  else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
368  {
369  Mx = 1.46;
370  fPDGencoding = 100321;
371  }
372  else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
373  {
374  Mx = 1.46;
375  fPDGencoding = -100321;
376  }
377  }
378  else
379  {
380  if( aParticle->GetDefinition() == G4Proton::Proton() )
381  {
382  Mx = 1.68;
383  fPDGencoding = 12216;
384  }
385  else if( aParticle->GetDefinition() == G4Neutron::Neutron() )
386  {
387  Mx = 1.68;
388  fPDGencoding = 12116;
389  }
390  else if( aParticle->GetDefinition() == G4PionPlus::PionPlus() )
391  {
392  Mx = 1.67;
393  fPDGencoding = 10215; // pi2(1670)+
394  // Mx = 1.45;
395  // fPDGencoding = 10211;
396  }
397  else if( aParticle->GetDefinition() == G4PionMinus::PionMinus() )
398  {
399  Mx = 1.67; // f0 problems->4pi vmg 20.11.14
400  fPDGencoding = -10215; // pi2(1670)-
401  // Mx = 1.45;
402  // fPDGencoding = -10211;
403  }
404  else if( aParticle->GetDefinition() == G4KaonPlus::KaonPlus() )
405  {
406  Mx = 1.68;
407  fPDGencoding = 30323;
408  }
409  else if( aParticle->GetDefinition() == G4KaonMinus::KaonMinus() )
410  {
411  Mx = 1.68;
412  fPDGencoding = -30323;
413  }
414  }
415  if(fPDGencoding == 0)
416  {
417  Mx = 1.44;
418  fPDGencoding = 12212;
419  }
420  G4ParticleDefinition* myResonance =
422 
423  if ( myResonance ) Mx = myResonance->GetPDGMass();
424 
425  // G4cout<<"PDG-ID = "<<fPDGencoding<<"; with mass = "<<Mx/CLHEP::GeV<<" GeV"<<G4endl;
426 
427  return Mx/CLHEP::GeV;
428 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static const G4double fProbMx[60][2]
int G4int
Definition: G4Types.hh:78
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static const double GeV
#define G4UniformRand()
Definition: Randomize.hh:97
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4ParticleDefinition * GetDefinition() const
static G4ParticleTable * GetParticleTable()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleT()

G4double G4LMsdGenerator::SampleT ( const G4HadProjectile aParticle,
G4double  Mx 
)

Definition at line 434 of file G4LMsdGenerator.cc.

436 {
437  G4double t=0., b=0., rTkin = 50.*CLHEP::GeV, eTkin = aParticle->GetKineticEnergy();
438  G4int i;
439 
440  for( i = 0; i < 23; ++i)
441  {
442  if( Mx <= fMxBdata[i][0] ) break;
443  }
444  if( i <= 0 ) b = fMxBdata[0][1];
445  else if( i >= 22 ) b = fMxBdata[22][1];
446  else b = fMxBdata[i][1];
447 
448  if( eTkin > rTkin ) b *= 1. + G4Log(eTkin/rTkin);
449 
450  G4double rand = G4UniformRand();
451 
452  t = -G4Log(rand)/b;
453 
454  t *= (CLHEP::GeV*CLHEP::GeV); // in G4 internal units
455 
456  return t;
457 }
int G4int
Definition: G4Types.hh:78
static const double GeV
#define G4UniformRand()
Definition: Randomize.hh:97
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
static const G4double fMxBdata[23][2]
Here is the call graph for this function:

Member Data Documentation

◆ fMxBdata

const G4double G4LMsdGenerator::fMxBdata
staticprivate
Initial value:
=
{
{1.09014, 17.8620},
{1.12590, 19.2831},
{1.18549, 17.6907},
{1.21693, 16.4760},
{1.25194, 15.3867},
{1.26932, 14.4236},
{1.29019, 13.2931},
{1.30755, 12.2882},
{1.31790, 11.4509},
{1.33888, 10.6969},
{1.34911, 9.44130},
{1.37711, 8.56148},
{1.39101, 7.76593},
{1.42608, 6.88582},
{1.48593, 6.13019},
{1.53179, 5.87723},
{1.58111, 5.37308},
{1.64105, 4.95217},
{1.69037, 4.44803},
{1.81742, 3.89879},
{1.88096, 3.68693},
{1.95509, 3.43278},
{2.02219, 3.30445}
}

Definition at line 95 of file G4LMsdGenerator.hh.

◆ fPDGencoding

G4int G4LMsdGenerator::fPDGencoding
private

Definition at line 93 of file G4LMsdGenerator.hh.

◆ fProbMx

const G4double G4LMsdGenerator::fProbMx
staticprivate

Definition at line 96 of file G4LMsdGenerator.hh.


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