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