Geant4 10.7.0
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 45 of file G4LMsdGenerator.cc.

47
48{
49 fPDGencoding = 0;
50
51 // theParticleChange = new G4HadFinalState;
52}

◆ ~G4LMsdGenerator()

G4LMsdGenerator::~G4LMsdGenerator ( )

Definition at line 54 of file G4LMsdGenerator.cc.

55{
56 // delete theParticleChange;
57}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 112 of file G4LMsdGenerator.cc.

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

78{
79 G4bool applied = false;
80
81 if( ( aTrack.GetDefinition() == G4Proton::Proton() ||
82 aTrack.GetDefinition() == G4Neutron::Neutron() ) &&
83 targetNucleus.GetA_asInt() >= 1 &&
84 // aTrack.GetKineticEnergy() > 1800*CLHEP::MeV
85 aTrack.GetKineticEnergy() > 300*CLHEP::MeV
86 ) // 750*CLHEP::MeV )
87 {
88 applied = true;
89 }
90 else if( ( aTrack.GetDefinition() == G4PionPlus::PionPlus() ||
91 aTrack.GetDefinition() == G4PionMinus::PionMinus() ) &&
92 targetNucleus.GetA_asInt() >= 1 &&
93 aTrack.GetKineticEnergy() > 2340*CLHEP::MeV )
94 {
95 applied = true;
96 }
97 else if( ( aTrack.GetDefinition() == G4KaonPlus::KaonPlus() ||
98 aTrack.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
99 targetNucleus.GetA_asInt() >= 1 &&
100 aTrack.GetKineticEnergy() > 1980*CLHEP::MeV )
101 {
102 applied = true;
103 }
104 return applied;
105}
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 59 of file G4LMsdGenerator.cc.

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

◆ SampleMx()

G4double G4LMsdGenerator::SampleMx ( const G4HadProjectile aParticle)

Definition at line 290 of file G4LMsdGenerator.cc.

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

Referenced by ApplyYourself().

◆ SampleT()

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

Definition at line 440 of file G4LMsdGenerator.cc.

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

Referenced by ApplyYourself().


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