Geant4 11.2.2
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 "G4UnitsTable.hh"
59#include "G4EmParameters.hh"
60#include "G4EmDataHandler.hh"
61#include "G4EmTableType.hh"
62#include "G4EmModelManager.hh"
64
65class G4Step;
66class G4VEmModel;
67class G4DataVector;
69class G4PhysicsTable;
70class G4PhysicsVector;
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77{
78public:
79
81
82 ~G4VEmProcess() override;
83
84 //------------------------------------------------------------------------
85 // Virtual methods to be implemented in concrete processes
86 //------------------------------------------------------------------------
87
88 virtual G4bool IsApplicable(const G4ParticleDefinition& p) override = 0;
89
90 void ProcessDescription(std::ostream& outFile) const override;
91
92protected:
93
94 virtual void StreamProcessInfo(std::ostream&) const {};
95
96 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
97
98 //------------------------------------------------------------------------
99 // Implementation of virtual methods common to all Discrete processes
100 //------------------------------------------------------------------------
101
102public:
103
104 // Initialise for build of tables
105 void PreparePhysicsTable(const G4ParticleDefinition&) override;
106
107 // Build physics table during initialisation
108 void BuildPhysicsTable(const G4ParticleDefinition&) override;
109
110 // Called before tracking of each new G4Track
111 void StartTracking(G4Track*) override;
112
113 // implementation of virtual method, specific for G4VEmProcess
115 const G4Track& track,
116 G4double previousStepSize,
117 G4ForceCondition* condition) override;
118
119 // implementation of virtual method, specific for G4VEmProcess
120 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
121
122 // Store PhysicsTable in a file.
123 // Return false in case of failure at I/O
125 const G4String& directory,
126 G4bool ascii = false) override;
127
128 // Retrieve Physics from a file.
129 // (return true if the Physics Table can be build by using file)
130 // (return false if the process has no functionality or in case of failure)
131 // File name should is constructed as processName+particleName and the
132 // should be placed under the directory specified by the argument.
134 const G4String& directory,
135 G4bool ascii) override;
136
137 // allowing check process name
138 virtual G4VEmProcess* GetEmProcess(const G4String& name);
139
140 //------------------------------------------------------------------------
141 // Specific methods for Discrete EM post step simulation
142 //------------------------------------------------------------------------
143
144 // The main method to access cross section per volume
145 inline G4double GetLambda(G4double kinEnergy,
146 const G4MaterialCutsCouple* couple,
147 G4double logKinEnergy);
148
149 // It returns the cross section per volume for energy/material
150 G4double GetCrossSection(const G4double kinEnergy,
151 const G4MaterialCutsCouple* couple) override;
152
153 // It returns the cross section of the process per atom
155 G4double Z, G4double A=0.,
156 G4double cut=0.0);
157
158 inline G4double MeanFreePath(const G4Track& track);
159
160 //------------------------------------------------------------------------
161 // Specific methods to build and access Physics Tables
162 //------------------------------------------------------------------------
163
164 // Binning for lambda table
165 void SetLambdaBinning(G4int nbins);
166
167 // Min kinetic energy for tables
169
170 // Min kinetic energy for high energy table
172
173 // Max kinetic energy for tables
175
176 // Cross section table pointers
177 inline G4PhysicsTable* LambdaTable() const;
178 inline G4PhysicsTable* LambdaTablePrim() const;
179 inline void SetLambdaTable(G4PhysicsTable*);
181
182 // Integral method type and peak positions
183 inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
184 inline void SetEnergyOfCrossSectionMax(std::vector<G4double>*);
187
188 //------------------------------------------------------------------------
189 // Define and access particle type
190 //------------------------------------------------------------------------
191
192 inline const G4ParticleDefinition* Particle() const;
193 inline const G4ParticleDefinition* SecondaryParticle() const;
194
195protected:
196
197 //------------------------------------------------------------------------
198 // Specific methods to set, access, modify models and basic parameters
199 //------------------------------------------------------------------------
200
201 // Select model in run time
202 inline G4VEmModel* SelectModel(G4double kinEnergy, size_t);
203
204public:
205
206 // Select model by energy and couple index
207 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
208 size_t idxCouple) const;
209
210 // Add model for region, smaller value of order defines which
211 // model will be selected for a given energy interval
212 void AddEmModel(G4int, G4VEmModel*, const G4Region* region = nullptr);
213
214 // Assign a model to a process local list, to enable the list in run time
215 // the derived process should execute AddEmModel(..) for all such models
216 void SetEmModel(G4VEmModel*, G4int index = 0);
217
218 inline G4int NumberOfModels() const;
219
220 // return a model from the local list
221 inline G4VEmModel* EmModel(size_t index = 0) const;
222
223 // Access to active model
224 inline const G4VEmModel* GetCurrentModel() const;
225
226 // Access to models
227 inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
228
229 // Access to the current G4Element
230 const G4Element* GetCurrentElement() const;
231
232 // Biasing parameters
233 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
234 inline G4double CrossSectionBiasingFactor() const;
235
236 // Activate forced interaction
237 void ActivateForcedInteraction(G4double length = 0.0,
238 const G4String& r = "",
239 G4bool flag = true);
240
241 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
242 G4double energyLimit);
243
244 inline void SetEmMasterProcess(const G4VEmProcess*);
245
246 inline void SetBuildTableFlag(G4bool val);
247
248 inline void CurrentSetup(const G4MaterialCutsCouple*, G4double energy);
249
250 inline G4bool UseBaseMaterial() const;
251
252 void BuildLambdaTable();
253
254 void StreamInfo(std::ostream& outFile, const G4ParticleDefinition&,
255 G4bool rst=false) const;
256
257 // hide copy constructor and assignment operator
259 G4VEmProcess & operator=(const G4VEmProcess &right) = delete;
260
261 //------------------------------------------------------------------------
262 // Other generic methods
263 //------------------------------------------------------------------------
264
265protected:
266
267 G4double GetMeanFreePath(const G4Track& track,
268 G4double previousStepSize,
269 G4ForceCondition* condition) override;
270
272
273 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
274
275 inline G4int LambdaBinning() const;
276
277 inline G4double MinKinEnergy() const;
278
279 inline G4double MaxKinEnergy() const;
280
281 // Single scattering parameters
282 inline G4double PolarAngleLimit() const;
283
285
286 inline void SetParticle(const G4ParticleDefinition* p);
287
288 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
289
290 inline size_t CurrentMaterialCutsCoupleIndex() const;
291
292 inline const G4MaterialCutsCouple* MaterialCutsCouple() const;
293
294 inline G4bool ApplyCuts() const;
295
297
299
300 inline void SetStartFromNullFlag(G4bool val);
301
302 inline void SetSplineFlag(G4bool val);
303
304 const G4Element* GetTargetElement() const;
305
306 const G4Isotope* GetTargetIsotope() const;
307
308 // these two methods assume that vectors are initilized
309 // and idx is within vector length
310 inline G4int DensityIndex(G4int idx) const;
311 inline G4double DensityFactor(G4int idx) const;
312
313private:
314
315 void PrintWarning(G4String tit, G4double val);
316
317 void ComputeIntegralLambda(G4double kinEnergy, const G4Track&);
318
319 inline G4double LogEkin(const G4Track&);
320
321 inline G4double GetLambdaFromTable(G4double kinEnergy);
322
323 inline G4double GetLambdaFromTable(G4double kinEnergy, G4double logKinEnergy);
324
325 inline G4double GetLambdaFromTablePrim(G4double kinEnergy);
326
327 inline G4double GetLambdaFromTablePrim(G4double kinEnergy, G4double logKinEnergy);
328
329 inline G4double GetCurrentLambda(G4double kinEnergy);
330
331 inline G4double GetCurrentLambda(G4double kinEnergy, G4double logKinEnergy);
332
333 inline G4double ComputeCurrentLambda(G4double kinEnergy);
334
335 // ======== pointers =========
336
337 G4EmModelManager* modelManager = nullptr;
338 const G4ParticleDefinition* particle = nullptr;
339 const G4ParticleDefinition* currentParticle = nullptr;
340 const G4ParticleDefinition* theGamma = nullptr;
341 const G4ParticleDefinition* theElectron = nullptr;
342 const G4ParticleDefinition* thePositron = nullptr;
343 const G4ParticleDefinition* secondaryParticle = nullptr;
344 const G4VEmProcess* masterProc = nullptr;
345 G4EmDataHandler* theData = nullptr;
346 G4VEmModel* currentModel = nullptr;
347 G4LossTableManager* lManager = nullptr;
348 G4EmParameters* theParameters = nullptr;
349 const G4Material* baseMaterial = nullptr;
350
351 // ======== tables and vectors ========
352 G4PhysicsTable* theLambdaTable = nullptr;
353 G4PhysicsTable* theLambdaTablePrim = nullptr;
354
355 const std::vector<G4double>* theCuts = nullptr;
356 const std::vector<G4double>* theCutsGamma = nullptr;
357 const std::vector<G4double>* theCutsElectron = nullptr;
358 const std::vector<G4double>* theCutsPositron = nullptr;
359
360protected:
361
362 // ======== pointers =========
363
365 const G4Material* currentMaterial = nullptr;
367 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
368
369private:
370
371 const std::vector<G4double>* theDensityFactor = nullptr;
372 const std::vector<G4int>* theDensityIdx = nullptr;
373
374 // ======== parameters =========
375 G4double minKinEnergy;
376 G4double maxKinEnergy;
377 G4double minKinEnergyPrim = DBL_MAX;
378 G4double lambdaFactor = 0.8;
379 G4double invLambdaFactor;
380 G4double biasFactor = 1.0;
381 G4double massRatio = 1.0;
382 G4double fFactor = 1.0;
383 G4double fLambda = 0.0;
384 G4double fLambdaEnergy = 0.0;
385
386protected:
387
391
392private:
393
395
396 G4int numberOfModels = 0;
397 G4int nLambdaBins = 84;
398
399protected:
400
409 size_t coupleIdxLambda = 0;
410 size_t idxLambda = 0;
411
414
415private:
416
417 G4bool buildLambdaTable = true;
418 G4bool applyCuts = false;
419 G4bool startFromNull = false;
420 G4bool splineFlag = true;
421 G4bool actMinKinEnergy = false;
422 G4bool actMaxKinEnergy = false;
423 G4bool actBinning = false;
424 G4bool isIon = false;
425 G4bool biasFlag = false;
426 G4bool weightFlag = false;
427
428protected:
429
430 // ======== particle change =========
431 std::vector<G4DynamicParticle*> secParticles;
433
434private:
435
436 // ======== local vectors =========
437 std::vector<G4VEmModel*> emModels;
438
439};
440
441// ======== Run time inline methods ================
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 baseMaterial = currentMaterial = couple->GetMaterial();
479 fFactor = biasFactor;
481 if(baseMat) {
482 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
483 if(nullptr != currentMaterial->GetBaseMaterial())
484 baseMaterial = currentMaterial->GetBaseMaterial();
485 fFactor *= (*theDensityFactor)[currentCoupleIndex];
486 }
487 }
488}
489
490//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
491
492inline
494{
495 if(1 < numberOfModels) {
496 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
497 }
498 currentModel->SetCurrentCouple(currentCouple);
499 return currentModel;
500}
501
502//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
503
504inline
506 size_t idxCouple) const
507{
508 return modelManager->SelectModel(kinEnergy, idxCouple);
509}
510
511//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
512
513inline G4double G4VEmProcess::GetLambdaFromTable(G4double e)
514{
515 return ((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
516}
517
518//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
519
520inline G4double G4VEmProcess::LogEkin(const G4Track& track)
521{
522 return track.GetDynamicParticle()->GetLogKineticEnergy();
523}
524
525//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
526
527inline G4double G4VEmProcess::GetLambdaFromTable(G4double e, G4double loge)
528{
529 return ((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
530}
531
532//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
533
534inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e)
535{
536 return ((*theLambdaTablePrim)[basedCoupleIndex])->Value(e, idxLambda)/e;
537}
538
539//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
540
541inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e, G4double loge)
542{
543 return ((*theLambdaTablePrim)[basedCoupleIndex])->LogVectorValue(e, loge)/e;
544}
545
546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547
548inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
549{
550 return currentModel->CrossSectionPerVolume(baseMaterial, currentParticle, e);
551}
552
553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
554
555inline G4double G4VEmProcess::GetCurrentLambda(G4double e)
556{
557 if(currentCoupleIndex != coupleIdxLambda || fLambdaEnergy != e) {
559 fLambdaEnergy = e;
560 if(e >= minKinEnergyPrim) { fLambda = GetLambdaFromTablePrim(e); }
561 else if(nullptr != theLambdaTable) { fLambda = GetLambdaFromTable(e); }
562 else { fLambda = ComputeCurrentLambda(e); }
563 fLambda *= fFactor;
564 }
565 return fLambda;
566}
567
568//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
569
570inline G4double G4VEmProcess::GetCurrentLambda(G4double e, G4double loge)
571{
572 if(currentCoupleIndex != coupleIdxLambda || fLambdaEnergy != e) {
574 fLambdaEnergy = e;
575 if(e >= minKinEnergyPrim) { fLambda = GetLambdaFromTablePrim(e, loge); }
576 else if(nullptr != theLambdaTable) { fLambda = GetLambdaFromTable(e, loge); }
577 else { fLambda = ComputeCurrentLambda(e); }
578 fLambda *= fFactor;
579 }
580 return fLambda;
581}
582
583//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584
585inline void
587{
588 DefineMaterial(couple);
589 SelectModel(energy*massRatio, currentCoupleIndex);
590}
591
592//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
593
594inline G4double
596 G4double logKinEnergy)
597{
598 CurrentSetup(couple, kinEnergy);
599 return GetCurrentLambda(kinEnergy, logKinEnergy);
600}
601
602//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
603
605{
606 const G4double kinEnergy = track.GetKineticEnergy();
607 CurrentSetup(track.GetMaterialCutsCouple(), kinEnergy);
608 const G4double xs = GetCurrentLambda(kinEnergy,
610 return (0.0 < xs) ? 1.0/xs : DBL_MAX;
611}
612
613// ======== Get/Set inline methods used at initialisation ================
614
616{
617 return applyCuts;
618}
619
620//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
621
623{
624 return nLambdaBins;
625}
626
627//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
628
630{
631 return minKinEnergy;
632}
633
634//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635
637{
638 return maxKinEnergy;
639}
640
641//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
642
644{
645 return biasFactor;
646}
647
648//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
649
651{
652 return theLambdaTable;
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
658{
659 return theLambdaTablePrim;
660}
661
662//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
663
665{
666 theLambdaTable = ptr;
667}
668
669//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
670
672{
673 theLambdaTablePrim = ptr;
674}
675
676//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677
678inline std::vector<G4double>* G4VEmProcess::EnergyOfCrossSectionMax() const
679{
681}
682
683//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
684
685inline void
687{
689}
690
691//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
692
694{
695 return particle;
696}
697
698//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
699
701{
702 return secondaryParticle;
703}
704
705//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
706
708{
709 fXSType = val;
710}
711
712//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
713
715{
716 return fXSType;
717}
718
719//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
720
722{
723 buildLambdaTable = val;
724}
725
726//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
727
732
733//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
734
736{
737 particle = p;
738 currentParticle = p;
739}
740
741//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742
744{
745 secondaryParticle = p;
746}
747
748//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
749
751{
752 startFromNull = val;
753}
754
755//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
756
758{
759 splineFlag = val;
760}
761
762//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
763
765{
766 return (*theDensityIdx)[idx];
767}
768
769//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
770
772{
773 return (*theDensityFactor)[idx];
774}
775
776//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
777
779{
780 return baseMat;
781}
782
783//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
784
786{
787 return currentModel;
788}
789
790//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
791
793{
794 masterProc = ptr;
795}
796
797//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
798
800{
801 return numberOfModels;
802}
803
804//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
805
806inline G4VEmModel* G4VEmProcess::EmModel(size_t index) const
807{
808 return (index < emModels.size()) ? emModels[index] : nullptr;
809}
810
811//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
812
814{
815 return modelManager->GetModel(idx, ver);
816}
817
818//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
819
820#endif
G4CrossSectionType
@ fEmNoIntegral
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
const G4double A[17]
G4double GetLogKineticEnergy() const
G4VEmModel * SelectModel(G4double energy, std::size_t index)
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
const G4Material * GetMaterial() const
const G4Material * GetBaseMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4CrossSectionType CrossSectionType() const
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4bool ApplyCuts() const
G4double GetGammaEnergyCut()
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)
const G4ParticleDefinition * Particle() const
virtual void StreamProcessInfo(std::ostream &) const
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
G4EmBiasingManager * biasManager
void BuildPhysicsTable(const G4ParticleDefinition &) override
void SetLambdaTablePrim(G4PhysicsTable *)
G4double preStepLambda
G4VEmModel * EmModel(size_t index=0) const
void SetBuildTableFlag(G4bool val)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4VEmProcess(G4VEmProcess &)=delete
void SetEnergyOfCrossSectionMax(std::vector< G4double > *)
G4PhysicsTable * LambdaTable() const
std::vector< G4double > * theEnergyOfCrossSectionMax
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
void SetMinKinEnergy(G4double e)
size_t basedCoupleIndex
G4double GetElectronEnergyCut()
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4int DensityIndex(G4int idx) const
void StartTracking(G4Track *) override
G4double GetCrossSection(const G4double kinEnergy, const G4MaterialCutsCouple *couple) 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 StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
size_t coupleIdxLambda
~G4VEmProcess() override
void SetLambdaBinning(G4int nbins)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void SetSecondaryParticle(const G4ParticleDefinition *p)
const G4VEmModel * GetCurrentModel() const
void SetSplineFlag(G4bool val)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t idxCouple) const
std::vector< G4double > * EnergyOfCrossSectionMax() const
G4bool UseBaseMaterial() const
G4int NumberOfModels() const
void SetCrossSectionType(G4CrossSectionType val)
void ProcessDescription(std::ostream &outFile) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VEmProcess & operator=(const G4VEmProcess &right)=delete
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
const G4Element * GetTargetElement() const
G4double DensityFactor(G4int idx) const
virtual G4VEmProcess * GetEmProcess(const G4String &name)
G4double MaxKinEnergy() const
G4ParticleChangeForGamma * GetParticleChange()
const G4Isotope * GetTargetIsotope() const
G4double MinKinEnergy() const
std::vector< G4DynamicParticle * > secParticles
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
const G4MaterialCutsCouple * currentCouple
void SetStartFromNullFlag(G4bool val)
G4PhysicsTable * LambdaTablePrim() const
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
G4ParticleChangeForGamma fParticleChange
void SetLambdaTable(G4PhysicsTable *)
G4int LambdaBinning() const
void SetParticle(const G4ParticleDefinition *p)
void SetMinKinEnergyPrim(G4double e)
void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMaxKinEnergy(G4double e)
const G4Material * currentMaterial
size_t CurrentMaterialCutsCoupleIndex() const
const G4Element * GetCurrentElement() const
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void BuildLambdaTable()
const G4ParticleDefinition * SecondaryParticle() const
#define DBL_MAX
Definition templates.hh:62