Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEmProcess.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// GEANT4 Class header file
29//
30//
31// File name: G4VEmProcess
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 01.10.2003
36//
37// Modifications: Vladimir Ivanchenko
38//
39// Class Description:
40//
41// It is the base class - EM discrete and rest/discrete process
42
43// -------------------------------------------------------------------
44//
45
46#ifndef G4VEmProcess_h
47#define G4VEmProcess_h 1
48
50
51#include "G4VDiscreteProcess.hh"
52#include "globals.hh"
53#include "G4Material.hh"
55#include "G4Track.hh"
56#include "G4EmModelManager.hh"
57#include "G4UnitsTable.hh"
60#include "G4EmDataHandler.hh"
61#include "G4EmParameters.hh"
62
63class G4Step;
64class G4VEmModel;
65class G4DataVector;
67class G4PhysicsTable;
68class G4PhysicsVector;
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75{
76public:
77
79
80 virtual ~G4VEmProcess();
81
82 //------------------------------------------------------------------------
83 // Virtual methods to be implemented in concrete processes
84 //------------------------------------------------------------------------
85
86 virtual G4bool IsApplicable(const G4ParticleDefinition& p) override = 0;
87
88 // obsolete
89 virtual void PrintInfo() {};
90
91 virtual void ProcessDescription(std::ostream& outFile) const override;
92
93protected:
94
95 virtual void StreamProcessInfo(std::ostream&) const {};
96
97 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
98
99 //------------------------------------------------------------------------
100 // Method with standard implementation; may be overwritten if needed
101 //------------------------------------------------------------------------
102
104 const G4Material*);
105
106 //------------------------------------------------------------------------
107 // Implementation of virtual methods common to all Discrete processes
108 //------------------------------------------------------------------------
109
110public:
111
112 // Initialise for build of tables
113 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
114
115 // Build physics table during initialisation
116 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
117
118 // Called before tracking of each new G4Track
119 virtual void StartTracking(G4Track*) override;
120
121 // implementation of virtual method, specific for G4VEmProcess
123 const G4Track& track,
124 G4double previousStepSize,
125 G4ForceCondition* condition) override;
126
127 // implementation of virtual method, specific for G4VEmProcess
128 virtual G4VParticleChange* PostStepDoIt(const G4Track&,
129 const G4Step&) override;
130
131 // Store PhysicsTable in a file.
132 // Return false in case of failure at I/O
134 const G4String& directory,
135 G4bool ascii = false) override;
136
137 // Retrieve Physics from a file.
138 // (return true if the Physics Table can be build by using file)
139 // (return false if the process has no functionality or in case of failure)
140 // File name should is constructed as processName+particleName and the
141 // should be placed under the directory specified by the argument.
143 const G4String& directory,
144 G4bool ascii) override;
145
146 // allowing check process name
147 virtual G4VEmProcess* GetEmProcess(const G4String& name);
148
149 //------------------------------------------------------------------------
150 // Specific methods for Discrete EM post step simulation
151 //------------------------------------------------------------------------
152
153 // It returns the cross section per volume for energy/ material
155 const G4MaterialCutsCouple* couple,
156 G4double logKinEnergy = DBL_MAX);
157
158 // It returns the cross section of the process per atom
160 G4double Z, G4double A=0.,
161 G4double cut=0.0);
162
163 G4double MeanFreePath(const G4Track& track);
164
165 // Obsolete method to access cross section per volume
166 G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple* couple);
167
168 // The main method to access cross section per volume
169 inline G4double GetLambda(G4double kinEnergy,
170 const G4MaterialCutsCouple* couple,
171 G4double logKinEnergy);
172
173 //------------------------------------------------------------------------
174 // Specific methods to build and access Physics Tables
175 //------------------------------------------------------------------------
176
177 // Binning for lambda table
178 void SetLambdaBinning(G4int nbins);
179
180 // Min kinetic energy for tables
182
183 // Min kinetic energy for high energy table
185
186 // Max kinetic energy for tables
188
189 // Cross section table pointers
190 inline G4PhysicsTable* LambdaTable() const;
191 inline G4PhysicsTable* LambdaTablePrim() const;
192
193 //------------------------------------------------------------------------
194 // Define and access particle type
195 //------------------------------------------------------------------------
196
197 inline const G4ParticleDefinition* Particle() const;
198 inline const G4ParticleDefinition* SecondaryParticle() const;
199
200 //------------------------------------------------------------------------
201 // Specific methods to set, access, modify models and basic parameters
202 //------------------------------------------------------------------------
203
204protected:
205 // Select model in run time
206 inline G4VEmModel* SelectModel(G4double kinEnergy, size_t index);
207
208public:
209 // Select model by energy and region index
210 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
211 size_t idxRegion) const;
212
213 // Add model for region, smaller value of order defines which
214 // model will be selected for a given energy interval
215 void AddEmModel(G4int, G4VEmModel*, const G4Region* region = nullptr);
216
217 // Assign a model to a process local list, to enable the list in run time
218 // the derived process should execute AddEmModel(..) for all such models
219 void SetEmModel(G4VEmModel*, G4int index = 0);
220
221 // return a model from the local list
222 G4VEmModel* EmModel(size_t index = 0) const;
223
224 // Define new energy range for the model identified by the name
226
227 // Access to models
228 G4int GetNumberOfModels() const;
229 G4int GetNumberOfRegionModels(size_t couple_index) const;
230 G4VEmModel* GetRegionModel(G4int idx, size_t couple_index) const;
231 G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
232
233 // Access to active model
234 inline const G4VEmModel* GetCurrentModel() const;
235
236 // Access to the current G4Element
237 const G4Element* GetCurrentElement() const;
238
239 // Biasing parameters
240 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
241 inline G4double CrossSectionBiasingFactor() const;
242
243 // Activate forced interaction
244 void ActivateForcedInteraction(G4double length = 0.0,
245 const G4String& r = "",
246 G4bool flag = true);
247
248 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
249 G4double energyLimit);
250
251 inline void SetEmMasterProcess(const G4VEmProcess*);
252
253 inline void SetIntegral(G4bool val);
254
255 inline void SetBuildTableFlag(G4bool val);
256
257 inline void CurrentSetup(const G4MaterialCutsCouple*, G4double energy);
258
259 //------------------------------------------------------------------------
260 // Other generic methods
261 //------------------------------------------------------------------------
262
263protected:
264
265 virtual G4double GetMeanFreePath(const G4Track& track,
266 G4double previousStepSize,
267 G4ForceCondition* condition) override;
268
270
271 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
272
273 inline G4int LambdaBinning() const;
274
275 inline G4double MinKinEnergy() const;
276
277 inline G4double MaxKinEnergy() const;
278
279 // Single scattering parameters
280 inline G4double PolarAngleLimit() const;
281
282 inline G4bool IsIntegral() const;
283
284 inline G4double RecalculateLambda(G4double kinEnergy,
285 const G4MaterialCutsCouple* couple);
286
288
289 inline void SetParticle(const G4ParticleDefinition* p);
290
291 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
292
293 inline size_t CurrentMaterialCutsCoupleIndex() const;
294
295 inline const G4MaterialCutsCouple* MaterialCutsCouple() const;
296
297 inline G4bool ApplyCuts() const;
298
300
302
303 inline void SetStartFromNullFlag(G4bool val);
304
305 inline void SetSplineFlag(G4bool val);
306
307 inline const G4Element* GetTargetElement() const;
308
309 inline const G4Isotope* GetTargetIsotope() const;
310
311private:
312
313 void Clear();
314
315 void BuildLambdaTable();
316
317 void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&,
318 G4bool rst=false) const;
319
320 void FindLambdaMax();
321
322 void PrintWarning(G4String tit, G4double val);
323
324 void ComputeIntegralLambda(G4double kinEnergy, G4double logKinEnergy);
325
326 // inline G4double GetLambdaFromTable(G4double kinEnergy);
327 inline G4double GetLambdaFromTable(G4double kinEnergy, G4double logKinEnergy);
328
329 // inline G4double GetLambdaFromTablePrim(G4double kinEnergy);
330 inline G4double GetLambdaFromTablePrim(G4double kinEnergy,
331 G4double logKinEnergy);
332
333 // inline G4double GetCurrentLambda(G4double kinEnergy);
334 inline G4double GetCurrentLambda(G4double kinEnergy, G4double logKinEnergy);
335
336 inline G4double ComputeCurrentLambda(G4double kinEnergy);
337
338 // hide copy constructor and assignment operator
339 G4VEmProcess(G4VEmProcess &) = delete;
340 G4VEmProcess & operator=(const G4VEmProcess &right) = delete;
341
342 // ======== Parameters of the class fixed at construction =========
343
344 G4EmModelManager* modelManager;
345 const G4ParticleDefinition* thePositron;
346 const G4ParticleDefinition* secondaryParticle;
347
348 G4bool buildLambdaTable;
349
350 // ======== Parameters of the class fixed at initialisation =======
351
352 std::vector<G4VEmModel*> emModels;
353 G4int numberOfModels;
354
355 // tables and vectors
356 G4PhysicsTable* theLambdaTable;
357 G4PhysicsTable* theLambdaTablePrim;
358 std::vector<G4double> theEnergyOfCrossSectionMax;
359 std::vector<G4double> theCrossSectionMax;
360
361 const std::vector<G4double>* theCuts;
362 const std::vector<G4double>* theCutsGamma;
363 const std::vector<G4double>* theCutsElectron;
364 const std::vector<G4double>* theCutsPositron;
365
366 G4int nLambdaBins;
367
368 G4double minKinEnergy;
369 G4double minKinEnergyPrim;
370 G4double maxKinEnergy;
371 G4double lambdaFactor;
372 G4double logLambdaFactor;
373 G4double biasFactor;
374 G4double massRatio;
375
376 G4bool integral;
377 G4bool applyCuts;
378 G4bool startFromNull;
379 G4bool splineFlag;
380 G4bool actMinKinEnergy;
381 G4bool actMaxKinEnergy;
382 G4bool actBinning;
383 G4bool actSpline;
384 G4bool isIon;
385
386 // ======== Cashed values - may be state dependent ================
387
388protected:
389
392
397 std::vector<G4DynamicParticle*> secParticles;
400 const std::vector<G4double>* theDensityFactor;
401 const std::vector<G4int>* theDensityIdx;
402
405
411
413
418
419private:
420
421 const G4VEmProcess* masterProc;
422 G4EmDataHandler* theData;
423 G4VEmModel* currentModel;
424
425 const G4ParticleDefinition* particle;
426
427 // cache
428 const G4ParticleDefinition* currentParticle;
429 const G4Material* baseMaterial;
430
431 G4double fFactor;
432 G4bool biasFlag;
433 G4bool weightFlag;
434};
435
436// ======== Run time inline methods ================
437
439{
440 return applyCuts;
441}
442
443//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
444
446{
447 return currentCoupleIndex;
448}
449
450//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
451
453{
454 return currentCouple;
455}
456
457//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458
460{
461 return (*theCutsGamma)[currentCoupleIndex];
462}
463
464//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
465
467{
468 return (*theCutsElectron)[currentCoupleIndex];
469}
470
471//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
472
474{
475 if(couple != currentCouple) {
476 currentCouple = couple;
477 currentMaterial = couple->GetMaterial();
478 baseMaterial = (currentMaterial->GetBaseMaterial())
480 currentCoupleIndex = couple->GetIndex();
481 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
482 fFactor = biasFactor*(*theDensityFactor)[currentCoupleIndex];
484 }
485}
486
487//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
488
489inline
491{
492 if(1 < numberOfModels) {
493 currentModel = modelManager->SelectModel(kinEnergy, index);
494 }
495 currentModel->SetCurrentCouple(currentCouple);
496 return currentModel;
497}
498
499//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
500
501inline
503 size_t idxRegion) const
504{
505 return modelManager->SelectModel(kinEnergy, idxRegion);
506}
507
508//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
509
510inline G4double G4VEmProcess::GetLambdaFromTable(G4double e, G4double loge)
511{
512 return ((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
513}
514
515//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
516
517inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e, G4double loge)
518{
519 return ((*theLambdaTablePrim)[basedCoupleIndex])->LogVectorValue(e, loge)/e;
520}
521
522//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
523
524inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
525{
526 return currentModel->CrossSectionPerVolume(baseMaterial, currentParticle, e);
527}
528
529//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
530
531inline G4double G4VEmProcess::GetCurrentLambda(G4double e, G4double loge)
532{
533 G4double x(0.0);
534 if(e >= minKinEnergyPrim) { x = GetLambdaFromTablePrim(e, loge); }
535 else if(theLambdaTable) { x = GetLambdaFromTable(e, loge); }
536 else if(currentModel) { x = ComputeCurrentLambda(e); }
537 return fFactor*x;
538}
539
540//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
541
542inline void
544{
545 DefineMaterial(couple);
546 SelectModel(energy*massRatio, currentCoupleIndex);
547}
548
549inline G4double
551 G4double logKinEnergy)
552{
553 CurrentSetup(couple, kinEnergy);
554 return GetCurrentLambda(kinEnergy, logKinEnergy);
555}
556
557//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
558
559inline G4double
561{
562 CurrentSetup(couple, e);
563 return fFactor*ComputeCurrentLambda(e);
564}
565
566// ======== Get/Set inline methods used at initialisation ================
567
568//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
569
571{
572 return nLambdaBins;
573}
574
575//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
576
578{
579 return minKinEnergy;
580}
581
582//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
583
585{
586 return maxKinEnergy;
587}
588
589//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
590
592{
594}
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
597
599{
600 return biasFactor;
601}
602
603//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
604
606{
607 return theLambdaTable;
608}
609
610//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
611
613{
614 return theLambdaTablePrim;
615}
616
617//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618
620{
621 return particle;
622}
623
624//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
625
627{
628 return secondaryParticle;
629}
630
631//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
632
634{
635 integral = val;
636}
637
638//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
639
641{
642 return integral;
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646
648{
649 buildLambdaTable = val;
650}
651
652//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
653
655{
656 return &fParticleChange;
657}
658
659//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
660
662{
663 particle = p;
664 currentParticle = p;
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668
670{
671 secondaryParticle = p;
672}
673
674//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
675
677{
678 startFromNull = val;
679}
680
681//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682
684{
685 splineFlag = val;
686 actSpline = true;
687}
688
689//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
690
692{
693 return currentModel->GetCurrentElement();
694}
695
696//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
697
699{
700 return currentModel->GetCurrentIsotope();
701}
702
703//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
704
706{
707 return currentModel;
708}
709
710//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
711
713{
714 masterProc = ptr;
715}
716
717//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
718
719#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4double MscThetaLimit() const
const G4Material * GetMaterial() const
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:231
Definition: G4Step.hh:62
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:465
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:495
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:502
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4double RecalculateLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4bool ApplyCuts() const
G4double GetGammaEnergyCut()
void SetIntegral(G4bool val)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double MeanFreePath(const G4Track &track)
G4int mainSecondaries
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
const G4ParticleDefinition * Particle() const
virtual void StreamProcessInfo(std::ostream &) const
Definition: G4VEmProcess.hh:95
G4LossTableManager * lManager
G4EmBiasingManager * biasManager
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double preStepLambda
G4VEmModel * EmModel(size_t index=0) const
void SetBuildTableFlag(G4bool val)
G4double preStepLogKinEnergy
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4PhysicsTable * LambdaTable() const
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
void SetMinKinEnergy(G4double e)
size_t basedCoupleIndex
G4double GetElectronEnergyCut()
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t idxRegion) const
virtual void PrintInfo()
Definition: G4VEmProcess.hh:89
virtual void StartTracking(G4Track *) override
G4double CrossSectionBiasingFactor() const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4double mfpKinEnergy
void SetEmMasterProcess(const G4VEmProcess *)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
void SetLambdaBinning(G4int nbins)
const std::vector< G4double > * theDensityFactor
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void SetSecondaryParticle(const G4ParticleDefinition *p)
const G4VEmModel * GetCurrentModel() const
void SetSplineFlag(G4bool val)
virtual void ProcessDescription(std::ostream &outFile) const override
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void UpdateEmModel(const G4String &, G4double, G4double)
const G4Element * GetTargetElement() const
virtual G4VEmProcess * GetEmProcess(const G4String &name)
G4double MaxKinEnergy() const
G4ParticleChangeForGamma * GetParticleChange()
const G4Isotope * GetTargetIsotope() const
G4VEmModel * SelectModel(G4double kinEnergy, size_t index)
G4bool isTheMaster
G4double MinKinEnergy() const
std::vector< G4DynamicParticle * > secParticles
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
const G4ParticleDefinition * theElectron
G4EmParameters * theParameters
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * theGamma
void SetStartFromNullFlag(G4bool val)
G4PhysicsTable * LambdaTablePrim() const
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
const G4MaterialCutsCouple * MaterialCutsCouple() const
G4double preStepKinEnergy
G4double PolarAngleLimit() const
size_t currentCoupleIndex
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4ParticleChangeForGamma fParticleChange
const std::vector< G4int > * theDensityIdx
G4bool IsIntegral() const
G4int LambdaBinning() const
void SetParticle(const G4ParticleDefinition *p)
G4VEmModel * GetRegionModel(G4int idx, size_t couple_index) const
void SetMinKinEnergyPrim(G4double e)
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMaxKinEnergy(G4double e)
G4int GetNumberOfRegionModels(size_t couple_index) const
const G4Material * currentMaterial
size_t CurrentMaterialCutsCoupleIndex() const
const G4Element * GetCurrentElement() const
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4int GetNumberOfModels() const
virtual ~G4VEmProcess()
const G4ParticleDefinition * SecondaryParticle() const
#define DBL_MAX
Definition: templates.hh:62