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

#include <G4LivermorePolarizedGammaConversionModel.hh>

+ Inheritance diagram for G4LivermorePolarizedGammaConversionModel:

Public Member Functions

 G4LivermorePolarizedGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedGammaConversion")
 
virtual ~G4LivermorePolarizedGammaConversionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 45 of file G4LivermorePolarizedGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedGammaConversionModel()

G4LivermorePolarizedGammaConversionModel::G4LivermorePolarizedGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermorePolarizedGammaConversion" 
)

Definition at line 43 of file G4LivermorePolarizedGammaConversionModel.cc.

46 isInitialised(false),meanFreePathTable(0),crossSectionHandler(0)
47{
48 lowEnergyLimit = 2*electron_mass_c2;
49 highEnergyLimit = 100 * GeV;
50 SetLowEnergyLimit(lowEnergyLimit);
51 SetHighEnergyLimit(highEnergyLimit);
52 smallEnergy = 2.*MeV;
53
54 Phi=0.;
55 Psi=0.;
56
57 verboseLevel= 0;
58 // Verbosity scale:
59 // 0 = nothing
60 // 1 = warning for energy non-conservation
61 // 2 = details of energy budget
62 // 3 = calculation of cross sections, file openings, samping of atoms
63 // 4 = entering in methods
64
65 if(verboseLevel > 0) {
66 G4cout << "Livermore Polarized GammaConversion is constructed " << G4endl
67 << "Energy range: "
68 << lowEnergyLimit / keV << " keV - "
69 << highEnergyLimit / GeV << " GeV"
70 << G4endl;
71 }
72
73 crossSectionHandler = new G4CrossSectionHandler();
74 crossSectionHandler->Initialise(0,lowEnergyLimit,highEnergyLimit,400);
75}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4LivermorePolarizedGammaConversionModel()

G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel ( )
virtual

Definition at line 79 of file G4LivermorePolarizedGammaConversionModel.cc.

80{
81 delete crossSectionHandler;
82}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 149 of file G4LivermorePolarizedGammaConversionModel.cc.

154{
155 if (verboseLevel > 3) {
156 G4cout << "G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom()"
157 << G4endl;
158 }
159 if(Z < 0.9 || GammaEnergy <= lowEnergyLimit) { return 0.0; }
160 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
161 return cs;
162}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double FindValue(G4int Z, G4double e) const

◆ GetMeanFreePath()

G4double G4LivermorePolarizedGammaConversionModel::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protected

◆ Initialise()

void G4LivermorePolarizedGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 88 of file G4LivermorePolarizedGammaConversionModel.cc.

90{
91 if (verboseLevel > 3)
92 G4cout << "Calling G4LivermorePolarizedGammaConversionModel::Initialise()" << G4endl;
93
94 if (crossSectionHandler)
95 {
96 crossSectionHandler->Clear();
97 delete crossSectionHandler;
98 }
99
100 // Energy limits
101 /*
102 // V.Ivanchenko: this was meanless check
103 if (LowEnergyLimit() < lowEnergyLimit)
104 {
105 G4cout << "G4LivermorePolarizedGammaConversionModel: low energy limit increased from " <<
106 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
107 // SetLowEnergyLimit(lowEnergyLimit);
108 }
109 */
110 if (HighEnergyLimit() > highEnergyLimit)
111 {
112 G4cout << "G4LivermorePolarizedGammaConversionModel: high energy limit decreased from " <<
113 HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
114 // V.Ivanchenko: this is forbidden
115 // SetHighEnergyLimit(highEnergyLimit);
116 }
117
118 // Reading of data files - all materials are read
119
120 crossSectionHandler = new G4CrossSectionHandler;
121 crossSectionHandler->Clear();
122 G4String crossSectionFile = "pair/pp-cs-";
123 crossSectionHandler->LoadData(crossSectionFile);
124
125 //
126 if (verboseLevel > 2) {
127 G4cout << "Loaded cross section files for Livermore Polarized GammaConversion model"
128 << G4endl;
129 }
130 InitialiseElementSelectors(particle,cuts);
131
132 if(verboseLevel > 0) {
133 G4cout << "Livermore Polarized GammaConversion model is initialized " << G4endl
134 << "Energy range: "
135 << LowEnergyLimit() / keV << " keV - "
136 << HighEnergyLimit() / GeV << " GeV"
137 << G4endl;
138 }
139
140 //
141 if(!isInitialised) {
142 isInitialised = true;
144 }
145}
void LoadData(const G4String &dataFile)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ SampleSecondaries()

void G4LivermorePolarizedGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 167 of file G4LivermorePolarizedGammaConversionModel.cc.

172{
173
174 // Fluorescence generated according to:
175 // J. Stepanek ,"A program to determine the radiation spectra due to a single atomic
176 // subshell ionisation by a particle or due to deexcitation or decay of radionuclides",
177 // Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
178
179 if (verboseLevel > 3)
180 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedGammaConversionModel" << G4endl;
181
182 G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
183 // Within energy limit?
184
185 if(photonEnergy <= lowEnergyLimit)
186 {
189 return;
190 }
191
192
193 G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
194 G4ThreeVector gammaDirection0 = aDynamicGamma->GetMomentumDirection();
195
196 // Make sure that the polarization vector is perpendicular to the
197 // gamma direction. If not
198
199 if(!(gammaPolarization0.isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.mag()==0))
200 { // only for testing now
201 gammaPolarization0 = GetRandomPolarization(gammaDirection0);
202 }
203 else
204 {
205 if ( gammaPolarization0.howOrthogonal(gammaDirection0) != 0)
206 {
207 gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
208 }
209 }
210
211 // End of Protection
212
213
214 G4double epsilon ;
215 G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
216
217 // Do it fast if photon energy < 2. MeV
218
219 if (photonEnergy < smallEnergy )
220 {
221 epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
222 }
223 else
224 {
225
226 // Select randomly one element in the current material
227
228 // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy);
229 //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
230
231 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
232 const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
233
234 /*
235 if (element == 0)
236 {
237 G4cout << "G4LivermorePolarizedGammaConversionModel::PostStepDoIt - element = 0" << G4endl;
238 }
239 */
240
241 G4IonisParamElm* ionisation = element->GetIonisation();
242
243 /*
244 if (ionisation == 0)
245 {
246 G4cout << "G4LivermorePolarizedGammaConversionModel::PostStepDoIt - ionisation = 0" << G4endl;
247 }
248 */
249
250
251 // Extract Coulomb factor for this Element
252
253 G4double fZ = 8. * (ionisation->GetlogZ3());
254 if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
255
256 // Limits of the screening variable
257 G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
258 G4double screenMax = exp ((42.24 - fZ)/8.368) - 0.952 ;
259 G4double screenMin = std::min(4.*screenFactor,screenMax) ;
260
261 // Limits of the energy sampling
262 G4double epsilon1 = 0.5 - 0.5 * sqrt(1. - screenMin / screenMax) ;
263 G4double epsilonMin = std::max(epsilon0Local,epsilon1);
264 G4double epsilonRange = 0.5 - epsilonMin ;
265
266 // Sample the energy rate of the created electron (or positron)
267 G4double screen;
268 G4double gReject ;
269
270 G4double f10 = ScreenFunction1(screenMin) - fZ;
271 G4double f20 = ScreenFunction2(screenMin) - fZ;
272 G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
273 G4double normF2 = std::max(1.5 * f20,0.);
274
275 do {
276 if (normF1 / (normF1 + normF2) > G4UniformRand() )
277 {
278 epsilon = 0.5 - epsilonRange * pow(G4UniformRand(), 0.3333) ;
279 screen = screenFactor / (epsilon * (1. - epsilon));
280 gReject = (ScreenFunction1(screen) - fZ) / f10 ;
281 }
282 else
283 {
284 epsilon = epsilonMin + epsilonRange * G4UniformRand();
285 screen = screenFactor / (epsilon * (1 - epsilon));
286 gReject = (ScreenFunction2(screen) - fZ) / f20 ;
287
288
289 }
290 } while ( gReject < G4UniformRand() );
291
292 } // End of epsilon sampling
293
294 // Fix charges randomly
295
296 G4double electronTotEnergy;
297 G4double positronTotEnergy;
298
299
300 if (G4int(2*G4UniformRand()))
301 {
302 electronTotEnergy = (1. - epsilon) * photonEnergy;
303 positronTotEnergy = epsilon * photonEnergy;
304 }
305 else
306 {
307 positronTotEnergy = (1. - epsilon) * photonEnergy;
308 electronTotEnergy = epsilon * photonEnergy;
309 }
310
311 // Scattered electron (positron) angles. ( Z - axis along the parent photon)
312 // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
313 // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
314
315/*
316 G4double u;
317 const G4double a1 = 0.625;
318 G4double a2 = 3. * a1;
319
320 if (0.25 > G4UniformRand())
321 {
322 u = - log(G4UniformRand() * G4UniformRand()) / a1 ;
323 }
324 else
325 {
326 u = - log(G4UniformRand() * G4UniformRand()) / a2 ;
327 }
328*/
329
330 G4double Ene = electronTotEnergy/electron_mass_c2; // Normalized energy
331
332 G4double cosTheta = 0.;
333 G4double sinTheta = 0.;
334
335 SetTheta(&cosTheta,&sinTheta,Ene);
336
337 // G4double theta = u * electron_mass_c2 / photonEnergy ;
338 // G4double phi = twopi * G4UniformRand() ;
339
340 G4double phi,psi=0.;
341
342 //corrected e+ e- angular angular distribution //preliminary!
343
344 // if(photonEnergy>50*MeV)
345 // {
346 phi = SetPhi(photonEnergy);
347 psi = SetPsi(photonEnergy,phi);
348 // }
349 //else
350 // {
351 //psi = G4UniformRand()*2.*pi;
352 //phi = pi; // coplanar
353 // }
354
355 Psi = psi;
356 Phi = phi;
357 //G4cout << "PHI " << phi << G4endl;
358 //G4cout << "PSI " << psi << G4endl;
359
360 G4double phie = psi; //azimuthal angle for the electron
361
362 G4double dirX = sinTheta*cos(phie);
363 G4double dirY = sinTheta*sin(phie);
364 G4double dirZ = cosTheta;
365 G4ThreeVector electronDirection(dirX,dirY,dirZ);
366 // Kinematics of the created pair:
367 // the electron and positron are assumed to have a symetric angular
368 // distribution with respect to the Z axis along the parent photon
369
370 //G4double localEnergyDeposit = 0. ;
371
372 G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
373
374 SystemOfRefChange(gammaDirection0,electronDirection,
375 gammaPolarization0);
376
378 electronDirection,
379 electronKineEnergy);
380
381 // The e+ is always created (even with kinetic energy = 0) for further annihilation
382
383 Ene = positronTotEnergy/electron_mass_c2; // Normalized energy
384
385 cosTheta = 0.;
386 sinTheta = 0.;
387
388 SetTheta(&cosTheta,&sinTheta,Ene);
389 G4double phip = phie+phi; //azimuthal angle for the positron
390
391 dirX = sinTheta*cos(phip);
392 dirY = sinTheta*sin(phip);
393 dirZ = cosTheta;
394 G4ThreeVector positronDirection(dirX,dirY,dirZ);
395
396 G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
397 SystemOfRefChange(gammaDirection0,positronDirection,
398 gammaPolarization0);
399
400 // Create G4DynamicParticle object for the particle2
402 positronDirection, positronKineEnergy);
403
404
405 fvect->push_back(particle1);
406 fvect->push_back(particle2);
407
408
409
410 // Kill the incident photon
411
412
413
414 // Create lists of pointers to DynamicParticles (photons and electrons)
415 // (Is the electron vector necessary? To be checked)
416 // std::vector<G4DynamicParticle*>* photonVector = 0;
417 //std::vector<G4DynamicParticle*> electronVector;
418
422
423}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
double mag() const
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:209
G4double GetlogZ3() const
G4double GetZ3() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LivermorePolarizedGammaConversionModel::fParticleChange
protected

Definition at line 73 of file G4LivermorePolarizedGammaConversionModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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