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

#include <G4UrbanFluctuation.hh>

+ Inheritance diagram for G4UrbanFluctuation:

Public Member Functions

 G4UrbanFluctuation (const G4String &nam="UrbanFluc")
 
 ~G4UrbanFluctuation () override
 
G4UrbanFluctuationoperator= (const G4UrbanFluctuation &right)=delete
 
 G4UrbanFluctuation (const G4UrbanFluctuation &)=delete
 
- Public Member Functions inherited from G4UniversalFluctuation
 G4UniversalFluctuation (const G4String &nam="UniFluc")
 
 ~G4UniversalFluctuation () override
 
G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double, const G4double, const G4double, const G4double) override
 
G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double, const G4double, const G4double) override
 
void InitialiseMe (const G4ParticleDefinition *) override
 
void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2) override
 
G4UniversalFluctuationoperator= (const G4UniversalFluctuation &right)=delete
 
 G4UniversalFluctuation (const G4UniversalFluctuation &)=delete
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length)=0
 
virtual void InitialiseMe (const G4ParticleDefinition *)
 
virtual void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2)
 
const G4StringGetName () const
 
G4VEmFluctuationModeloperator= (const G4VEmFluctuationModel &right)=delete
 
 G4VEmFluctuationModel (const G4VEmFluctuationModel &)=delete
 

Protected Member Functions

G4double SampleGlandz (CLHEP::HepRandomEngine *rndm, const G4Material *, const G4double tcut) override
 
- Protected Member Functions inherited from G4UniversalFluctuation
virtual G4double SampleGlandz (CLHEP::HepRandomEngine *rndm, const G4Material *, const G4double tcut)
 
void AddExcitation (CLHEP::HepRandomEngine *rndm, const G4double ax, const G4double ex, G4double &eav, G4double &eloss, G4double &esig2)
 
void SampleGauss (CLHEP::HepRandomEngine *rndm, const G4double eav, const G4double esig2, G4double &eloss)
 

Additional Inherited Members

- Protected Attributes inherited from G4UniversalFluctuation
G4double particleMass = 0.0
 
G4double m_Inv_particleMass = DBL_MAX
 
G4double m_massrate = DBL_MAX
 
G4double chargeSquare = 1.0
 
G4double ipotFluct = 0.0
 
G4double ipotLogFluct = 0.0
 
G4double e0 = 0.0
 
G4double minNumberInteractionsBohr = 10.0
 
G4double minLoss
 
G4double nmaxCont = 8.0
 
G4double rate = 0.56
 
G4double fw = 4.0
 
G4double a0 = 42.0
 
G4double w2 = 0.0
 
G4double meanLoss = 0.0
 
const G4ParticleDefinitionparticle = nullptr
 
G4doublerndmarray = nullptr
 
G4int sizearray = 30
 

Detailed Description

Definition at line 54 of file G4UrbanFluctuation.hh.

Constructor & Destructor Documentation

◆ G4UrbanFluctuation() [1/2]

G4UrbanFluctuation::G4UrbanFluctuation ( const G4String nam = "UrbanFluc")
explicit

◆ ~G4UrbanFluctuation()

G4UrbanFluctuation::~G4UrbanFluctuation ( )
overridedefault

◆ G4UrbanFluctuation() [2/2]

G4UrbanFluctuation::G4UrbanFluctuation ( const G4UrbanFluctuation )
delete

Member Function Documentation

◆ operator=()

G4UrbanFluctuation & G4UrbanFluctuation::operator= ( const G4UrbanFluctuation right)
delete

◆ SampleGlandz()

G4double G4UrbanFluctuation::SampleGlandz ( CLHEP::HepRandomEngine rndm,
const G4Material material,
const G4double  tcut 
)
overrideprotectedvirtual

Reimplemented from G4UniversalFluctuation.

Definition at line 65 of file G4UrbanFluctuation.cc.

68{
69 if (material != lastMaterial) {
70 auto ioni = material->GetIonisation();
71 f1Fluct = ioni->GetF1fluct();
72 f2Fluct = ioni->GetF2fluct();
73 e1Fluct = ioni->GetEnergy1fluct();
74 e2Fluct = ioni->GetEnergy2fluct();
75 e1LogFluct = ioni->GetLogEnergy1fluct();
76 e2LogFluct = ioni->GetLogEnergy2fluct();
77 esmall = 0.5*std::sqrt(e0*ipotFluct);
78 lastMaterial = material;
79 }
80
81 G4double a1(0.0), a2(0.0), a3(0.0);
82 G4double loss = 0.0;
83 G4double e1 = e1Fluct;
84 G4double e2 = e2Fluct;
85
86 if(tcut > ipotFluct) {
87 if(w2 > ipotLogFluct) {
88 if(w2 > e2LogFluct) {
89 const G4double C = meanLoss*(1.-rate)/(w2-ipotLogFluct);
90 a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
91 a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
92 } else {
93 a1 = meanLoss*(1.-rate)/e1;
94 }
95 if(a1 < a0) {
96 const G4double fwnow = 0.5+(fw-0.5)*std::sqrt(a1/a0);
97 a1 /= fwnow;
98 e1 *= fwnow;
99 } else {
100 a1 /= fw;
101 e1 *= fw;
102 }
103 }
104 }
105
106 const G4double w1 = tcut/e0;
107 a3 = rate*meanLoss*(tcut-e0)/(e0*tcut*G4Log(w1));
108 if(a1+a2 <= 0.) { a3 /= rate; }
109
110 //'nearly' Gaussian fluctuation if a1>nmaxCont&&a2>nmaxCont&&a3>nmaxCont
111 G4double emean = 0.;
112 G4double sig2e = 0.;
113
114 // excitation of type 1
115 if(a1 > 0.0) { AddExcitation(rndmEngineF, a1, e1, emean, loss, sig2e); }
116
117 // excitation of type 2
118 if(a2 > 0.0) { AddExcitation(rndmEngineF, a2, e2, emean, loss, sig2e); }
119
120 if(sig2e > 0.0) { SampleGauss(rndmEngineF, emean, sig2e, loss); }
121
122 // ionisation
123 if(a3 > 0.) {
124 emean = 0.;
125 sig2e = 0.;
126 G4double p3 = a3;
127 G4double alfa = 1.;
128 if(a3 > nmaxCont) {
129 alfa = w1*(nmaxCont+a3)/(w1*nmaxCont+a3);
130 const G4double alfa1 = alfa*G4Log(alfa)/(alfa-1.);
131 const G4double namean = a3*w1*(alfa-1.)/((w1-1.)*alfa);
132 emean += namean*e0*alfa1;
133 sig2e += e0*e0*namean*(alfa-alfa1*alfa1);
134 p3 -= namean;
135 }
136
137 const G4double w3 = alfa*e0;
138 if(tcut > w3) {
139 const G4double w = (tcut-w3)/tcut;
140 const G4int nnb = (G4int)G4Poisson(p3);
141 if(nnb > 0) {
142 if(nnb > sizearray) {
143 sizearray = nnb;
144 delete [] rndmarray;
145 rndmarray = new G4double[nnb];
146 }
147 rndmEngineF->flatArray(nnb, rndmarray);
148 for (G4int k=0; k<nnb; ++k) { loss += w3/(1.-w*rndmarray[k]); }
149 }
150 }
151 if(sig2e > 0.0) { SampleGauss(rndmEngineF, emean, sig2e, loss); }
152 }
153 //G4cout << "### loss=" << loss << G4endl;
154 return loss;
155}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
void AddExcitation(CLHEP::HepRandomEngine *rndm, const G4double ax, const G4double ex, G4double &eav, G4double &eloss, G4double &esig2)
void SampleGauss(CLHEP::HepRandomEngine *rndm, const G4double eav, const G4double esig2, G4double &eloss)

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