Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LMsdGenerator Class Reference

#include <G4LMsdGenerator.hh>

+ Inheritance 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 G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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)
 
G4int GetVerboseLevel () const
 
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
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

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 59 of file G4LMsdGenerator.hh.

Constructor & Destructor Documentation

◆ G4LMsdGenerator()

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

Definition at line 46 of file G4LMsdGenerator.cc.

47 : G4HadronicInteraction(name), secID(-1)
48
49{
50 fPDGencoding = 0;
51 secID = G4PhysicsModelCatalog::GetModelID( "model_LMsdGenerator" );
52
53 // theParticleChange = new G4HadFinalState;
54}
static G4int GetModelID(const G4int modelIndex)

◆ ~G4LMsdGenerator()

G4LMsdGenerator::~G4LMsdGenerator ( )

Definition at line 56 of file G4LMsdGenerator.cc.

57{
58 // delete theParticleChange;
59}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 114 of file G4LMsdGenerator.cc.

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

◆ IsApplicable()

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

Reimplemented from G4HadronicInteraction.

Definition at line 78 of file G4LMsdGenerator.cc.

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

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 61 of file G4LMsdGenerator.cc.

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

◆ SampleMx()

G4double G4LMsdGenerator::SampleMx ( const G4HadProjectile aParticle)

Definition at line 292 of file G4LMsdGenerator.cc.

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

Referenced by ApplyYourself().

◆ SampleT()

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

Definition at line 442 of file G4LMsdGenerator.cc.

444{
445 G4double t=0., b=0., rTkin = 50.*CLHEP::GeV, eTkin = aParticle->GetKineticEnergy();
446 G4int i;
447
448 for( i = 0; i < 23; ++i)
449 {
450 if( Mx <= fMxBdata[i][0] ) break;
451 }
452 if( i <= 0 ) b = fMxBdata[0][1];
453 else if( i >= 22 ) b = fMxBdata[22][1];
454 else b = fMxBdata[i][1];
455
456 if( eTkin > rTkin ) b *= 1. + G4Log(eTkin/rTkin);
457
458 G4double rand = G4UniformRand();
459
460 t = -G4Log(rand)/b;
461
462 t *= (CLHEP::GeV*CLHEP::GeV); // in G4 internal units
463
464 return t;
465}
G4double G4Log(G4double x)
Definition: G4Log.hh:227

Referenced by ApplyYourself().


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