Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VEmModel.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: G4VEmModel
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 03.01.2002
36//
37// Modifications:
38//
39// 23-12-02 V.Ivanchenko change interface before move to cut per region
40// 24-01-03 Cut per region (V.Ivanchenko)
41// 13-02-03 Add name (V.Ivanchenko)
42// 25-02-03 Add sample theta and displacement (V.Ivanchenko)
43// 23-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection
44// calculation (V.Ivanchenko)
45// 01-03-04 L.Urban signature changed in SampleCosineTheta
46// 23-04-04 L.urban signature of SampleCosineTheta changed back
47// 17-11-04 Add method CrossSectionPerAtom (V.Ivanchenko)
48// 14-03-05 Reduce number of pure virtual methods and make inline part
49// separate (V.Ivanchenko)
50// 24-03-05 Remove IsInCharge and add G4VParticleChange in the constructor (VI)
51// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
52// 15-04-05 optimize internal interface for msc (V.Ivanchenko)
53// 08-05-05 A -> N (V.Ivanchenko)
54// 25-07-05 Move constructor and destructor to the body (V.Ivanchenko)
55// 02-02-06 ComputeCrossSectionPerAtom: default value A=0. (mma)
56// 06-02-06 add method ComputeMeanFreePath() (mma)
57// 07-03-06 Optimize msc methods (V.Ivanchenko)
58// 29-06-06 Add member currentElement and Get/Set methods (V.Ivanchenko)
59// 29-10-07 Added SampleScattering (V.Ivanchenko)
60// 15-07-08 Reorder class members and improve comments (VI)
61// 21-07-08 Added vector of G4ElementSelector and methods to use it (VI)
62// 12-09-08 Added methods GetParticleCharge, GetChargeSquareRatio,
63// CorrectionsAlongStep, ActivateNuclearStopping (VI)
64// 16-02-09 Moved implementations of virtual methods to source (VI)
65// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
66// 13-10-10 Added G4VEmAngularDistribution (VI)
67//
68// Class Description:
69//
70// Abstract interface to energy loss models
71
72// -------------------------------------------------------------------
73//
74
75#ifndef G4VEmModel_h
76#define G4VEmModel_h 1
77
78#include "globals.hh"
79#include "G4DynamicParticle.hh"
82#include "G4Material.hh"
83#include "G4Element.hh"
84#include "G4ElementVector.hh"
85#include "G4Isotope.hh"
86#include "G4DataVector.hh"
91#include <vector>
92
93class G4ElementData;
94class G4PhysicsTable;
95class G4Region;
99class G4Track;
101
103{
104
105public:
106
107 explicit G4VEmModel(const G4String& nam);
108
109 virtual ~G4VEmModel();
110
111 //------------------------------------------------------------------------
112 // Virtual methods to be implemented for any concrete model
113 //------------------------------------------------------------------------
114
115 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) = 0;
116
117 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
119 const G4DynamicParticle*,
120 G4double tmin = 0.0,
121 G4double tmax = DBL_MAX) = 0;
122
123 //------------------------------------------------------------------------
124 // Methods for initialisation of MT; may be overwritten if needed
125 //------------------------------------------------------------------------
126
127 // initialisation in local thread
128 virtual void InitialiseLocal(const G4ParticleDefinition*,
129 G4VEmModel* masterModel);
130
131 // initialisation of a new material at run time
133 const G4Material*);
134
135 // initialisation of a new element at run time
136 virtual void InitialiseForElement(const G4ParticleDefinition*,
137 G4int Z);
138
139 //------------------------------------------------------------------------
140 // Methods with standard implementation; may be overwritten if needed
141 //------------------------------------------------------------------------
142
143 // main method to compute dEdx
146 G4double kineticEnergy,
147 G4double cutEnergy = DBL_MAX);
148
149 // main method to compute cross section per Volume
152 G4double kineticEnergy,
153 G4double cutEnergy = 0.0,
154 G4double maxEnergy = DBL_MAX);
155
156 // method to get partial cross section
158 G4int level,
160 G4double kineticEnergy);
161
162 // main method to compute cross section per atom
164 G4double kinEnergy,
165 G4double Z,
166 G4double A = 0., /* amu */
167 G4double cutEnergy = 0.0,
168 G4double maxEnergy = DBL_MAX);
169
170 // main method to compute cross section per atomic shell
172 G4int Z, G4int shellIdx,
173 G4double kinEnergy,
174 G4double cutEnergy = 0.0,
175 G4double maxEnergy = DBL_MAX);
176
177 // Compute effective ion charge square
178 virtual G4double ChargeSquareRatio(const G4Track&);
179
180 // Compute effective ion charge square
182 const G4Material*,
183 G4double kineticEnergy);
184
185 // Compute ion charge
187 const G4Material*,
188 G4double kineticEnergy);
189
190 // Initialisation for a new track
191 virtual void StartTracking(G4Track*);
192
193 // add correction to energy loss and compute non-ionizing energy loss
194 virtual void CorrectionsAlongStep(const G4Material*,
196 const G4double kinEnergy,
197 const G4double cutEnergy,
198 const G4double& length,
199 G4double& eloss);
200
201 // value which may be tabulated (by default cross section)
202 virtual G4double Value(const G4MaterialCutsCouple*,
204 G4double kineticEnergy);
205
206 // threshold for zero value
207 virtual G4double MinPrimaryEnergy(const G4Material*,
209 G4double cut = 0.0);
210
211 // model can define low-energy limit for the cut
213 const G4MaterialCutsCouple*);
214
215 // initialisation at run time for a given material
216 virtual void SetupForMaterial(const G4ParticleDefinition*,
217 const G4Material*,
218 G4double kineticEnergy);
219
220 // add a region for the model
221 virtual void DefineForRegion(const G4Region*);
222
223 // fill number of different type of secondaries after SampleSecondaries(...)
224 virtual void FillNumberOfSecondaries(G4int& numberOfTriplets,
225 G4int& numberOfRecoil);
226
227 // for automatic documentation
228 virtual void ModelDescription(std::ostream& outFile) const;
229
230protected:
231
232 // initialisation of the ParticleChange for the model
234
235 // initialisation of the ParticleChange for the model
237
238 // kinematically allowed max kinetic energy of a secondary
240 G4double kineticEnergy);
241
242public:
243
244 //------------------------------------------------------------------------
245 // Generic methods common to all models
246 //------------------------------------------------------------------------
247
248 // should be called at initialisation to build element selectors
250 const G4DataVector&);
251
252 // should be called at initialisation to access element selectors
253 inline std::vector<G4EmElementSelector*>* GetElementSelectors();
254
255 // should be called at initialisation to set element selectors
256 inline void SetElementSelectors(std::vector<G4EmElementSelector*>*);
257
258 // dEdx per unit length, base material approach may be used
261 G4double kineticEnergy,
262 G4double cutEnergy = DBL_MAX);
263
264 // cross section per volume, base material approach may be used
267 G4double kineticEnergy,
268 G4double cutEnergy = 0.0,
269 G4double maxEnergy = DBL_MAX);
270
271 // compute mean free path via cross section per volume
273 G4double kineticEnergy,
274 const G4Material*,
275 G4double cutEnergy = 0.0,
276 G4double maxEnergy = DBL_MAX);
277
278 // generic cross section per element
280 const G4Element*,
281 G4double kinEnergy,
282 G4double cutEnergy = 0.0,
283 G4double maxEnergy = DBL_MAX);
284
285 // atom can be selected effitiantly if element selectors are initialised
288 G4double kineticEnergy,
289 G4double cutEnergy = 0.0,
290 G4double maxEnergy = DBL_MAX);
291 // same as SelectRandomAtom above but more efficient since log-ekin is known
294 G4double kineticEnergy,
295 G4double logKineticEnergy,
296 G4double cutEnergy = 0.0,
297 G4double maxEnergy = DBL_MAX);
298
299 // to select atom cross section per volume is recomputed for each element
302 G4double kineticEnergy,
303 G4double cutEnergy = 0.0,
304 G4double maxEnergy = DBL_MAX);
305
306 // to select atom if cross section is proportional number of electrons
307 const G4Element* GetCurrentElement(const G4Material* mat = nullptr) const;
309
310 // select isotope in order to have precise mass of the nucleus
311 const G4Isotope* GetCurrentIsotope(const G4Element* elm = nullptr) const;
312 G4int SelectIsotopeNumber(const G4Element*) const;
313
314 //------------------------------------------------------------------------
315 // Get/Set methods
316 //------------------------------------------------------------------------
317
319
321
323
325
327
329
330 inline G4VEmModel* GetTripletModel();
331
332 inline void SetTripletModel(G4VEmModel*);
333
335
336 inline G4double HighEnergyLimit() const;
337
338 inline G4double LowEnergyLimit() const;
339
340 inline G4double HighEnergyActivationLimit() const;
341
342 inline G4double LowEnergyActivationLimit() const;
343
344 inline G4double PolarAngleLimit() const;
345
346 inline G4double SecondaryThreshold() const;
347
348 inline G4bool DeexcitationFlag() const;
349
350 inline G4bool ForceBuildTableFlag() const;
351
352 inline G4bool UseAngularGeneratorFlag() const;
353
354 inline void SetAngularGeneratorFlag(G4bool);
355
356 inline void SetHighEnergyLimit(G4double);
357
358 inline void SetLowEnergyLimit(G4double);
359
361
363
364 inline G4bool IsActive(G4double kinEnergy) const;
365
366 inline void SetPolarAngleLimit(G4double);
367
368 inline void SetSecondaryThreshold(G4double);
369
370 inline void SetDeexcitationFlag(G4bool val);
371
372 inline void SetForceBuildTable(G4bool val);
373
374 inline void SetFluctuationFlag(G4bool val);
375
376 inline G4bool IsMaster() const;
377
378 inline void SetUseBaseMaterials(G4bool val);
379
380 inline G4bool UseBaseMaterials() const;
381
382 inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
383
384 inline const G4String& GetName() const;
385
386 inline void SetCurrentCouple(const G4MaterialCutsCouple*);
387
388 inline G4bool IsLocked() const;
389
390 inline void SetLocked(G4bool);
391
392 // obsolete methods
393 [[deprecated("Use G4EmParameters::Instance()->SetLPM instead")]]
394 void SetLPMFlag(G4bool);
395
397
398 // hide assignment operator
399 G4VEmModel & operator=(const G4VEmModel &right) = delete;
400 G4VEmModel(const G4VEmModel&) = delete;
401
402protected:
403
404 inline const G4MaterialCutsCouple* CurrentCouple() const;
405
406 inline void SetCurrentElement(const G4Element*);
407
408private:
409
410 // ======== Parameters of the class fixed at construction =========
411
412 G4VEmFluctuationModel* flucModel = nullptr;
413 G4VEmAngularDistribution* anglModel = nullptr;
414 G4VEmModel* fTripletModel = nullptr;
415 const G4MaterialCutsCouple* fCurrentCouple = nullptr;
416 const G4Element* fCurrentElement = nullptr;
417 std::vector<G4EmElementSelector*>* elmSelectors = nullptr;
418 G4LossTableManager* fEmManager;
419
420protected:
421
425 const G4Material* pBaseMaterial = nullptr;
426 const std::vector<G4double>* theDensityFactor = nullptr;
427 const std::vector<G4int>* theDensityIdx = nullptr;
428
431
432private:
433
434 G4double lowLimit;
435 G4double highLimit;
436 G4double eMinActive = 0.0;
437 G4double eMaxActive = DBL_MAX;
438 G4double secondaryThreshold = DBL_MAX;
439 G4double polarAngleLimit;
440
441 G4int nSelectors = 0;
442 G4int nsec = 5;
443
444protected:
445
446 std::size_t currentCoupleIndex = 0;
447 std::size_t basedCoupleIndex = 0;
449
450private:
451
452 G4bool flagDeexcitation = false;
453 G4bool flagForceBuildTable = false;
454 G4bool isMaster = true;
455
456 G4bool localTable = true;
457 G4bool localElmSelectors = true;
458 G4bool useAngularGenerator = false;
459 G4bool useBaseMaterials = false;
460 G4bool isLocked = false;
461 G4bool localChange = false;
462
463 const G4String name;
464 std::vector<G4double> xsec;
465
466};
467
468// ======== Run time inline methods ================
469
471{
472 if(fCurrentCouple != ptr) {
473 fCurrentCouple = ptr;
475 pBaseMaterial = ptr->GetMaterial();
476 pFactor = 1.0;
477 if(useBaseMaterials) {
478 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
479 if(nullptr != pBaseMaterial->GetBaseMaterial())
480 pBaseMaterial = pBaseMaterial->GetBaseMaterial();
481 pFactor = (*theDensityFactor)[currentCoupleIndex];
482 }
483 }
484}
485
486//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
487
489{
490 return fCurrentCouple;
491}
492
493//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
494
496{
497 fCurrentElement = elm;
498}
499
500//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
501
502inline
508
509//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
510
512 const G4ParticleDefinition* part,
513 G4double kinEnergy,
514 G4double cutEnergy)
515{
516 SetCurrentCouple(couple);
517 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
518}
519
520//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
521
523 const G4ParticleDefinition* part,
524 G4double kinEnergy,
525 G4double cutEnergy,
526 G4double maxEnergy)
527{
528 SetCurrentCouple(couple);
529 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
530 cutEnergy,maxEnergy);
531}
532
533//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
534
535inline
537 G4double ekin,
538 const G4Material* material,
539 G4double emin,
540 G4double emax)
541{
542 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
543 return (cross > 0.0) ? 1./cross : DBL_MAX;
544}
545
546//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547
548inline G4double
550 const G4Element* elm,
551 G4double kinEnergy,
552 G4double cutEnergy,
553 G4double maxEnergy)
554{
555 fCurrentElement = elm;
556 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
557 cutEnergy,maxEnergy);
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
562inline const G4Element*
564 const G4ParticleDefinition* part,
565 G4double kinEnergy,
566 G4double cutEnergy,
567 G4double maxEnergy)
568{
569 SetCurrentCouple(couple);
570 fCurrentElement = (nSelectors > 0) ?
571 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
572 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
573 return fCurrentElement;
574}
575
576//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
577
578inline const G4Element*
580 const G4ParticleDefinition* part,
581 G4double kinEnergy,
582 G4double logKinE,
583 G4double cutEnergy,
584 G4double maxEnergy)
585{
586 SetCurrentCouple(couple);
587 fCurrentElement = (nSelectors > 0)
588 ? ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy,logKinE)
589 : SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
590 return fCurrentElement;
591}
592
593// ======== Get/Set inline methods used at initialisation ================
594
596{
597 return flucModel;
598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
603{
604 return anglModel;
605}
606
607//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
608
610{
611 if(p != anglModel) {
612 delete anglModel;
613 anglModel = p;
614 }
615}
616
617//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618
620{
621 return fTripletModel;
622}
623
624//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
625
627{
628 if(p != fTripletModel) {
629 delete fTripletModel;
630 fTripletModel = p;
631 }
632}
633
634//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635
637{
638 return highLimit;
639}
640
641//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
642
644{
645 return lowLimit;
646}
647
648//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
649
651{
652 return eMaxActive;
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
658{
659 return eMinActive;
660}
661
662//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
663
665{
666 return polarAngleLimit;
667}
668
669//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
670
672{
673 return secondaryThreshold;
674}
675
676//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677
679{
680 return flagDeexcitation;
681}
682
683//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
684
686{
687 return flagForceBuildTable;
688}
689
690//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691
693{
694 return useAngularGenerator;
695}
696
697//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698
700{
701 useAngularGenerator = val;
702}
703
704//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705
707{
708 lossFlucFlag = val;
709}
710
711//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
712
714{
715 return isMaster;
716}
717
718//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
719
721{
722 useBaseMaterials = val;
723}
724
725//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
726
728{
729 return useBaseMaterials;
730}
731
732//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
733
735{
736 highLimit = val;
737}
738
739//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
740
742{
743 lowLimit = val;
744}
745
746//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
747
749{
750 eMaxActive = val;
751}
752
753//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
754
756{
757 eMinActive = val;
758}
759
760//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
761
762inline G4bool G4VEmModel::IsActive(G4double kinEnergy) const
763{
764 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
765}
766
767//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
768
770{
771 if(!isLocked) { polarAngleLimit = val; }
772}
773
774//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
775
777{
778 secondaryThreshold = val;
779}
780
781//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
782
784{
785 flagDeexcitation = val;
786}
787
788//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
789
791{
792 flagForceBuildTable = val;
793}
794
795//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
796
797inline const G4String& G4VEmModel::GetName() const
798{
799 return name;
800}
801
802//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
803
804inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
805{
806 return elmSelectors;
807}
808
809//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
810
811inline void
812G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
813{
814 if(p != elmSelectors) {
815 elmSelectors = p;
816 nSelectors = (nullptr != elmSelectors) ? G4int(elmSelectors->size()) : 0;
817 localElmSelectors = false;
818 }
819}
820
821//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
822
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
829
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836
838{
839 return isLocked;
840}
841
842//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
843
845{
846 isLocked = val;
847}
848
849//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
850
851#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition G4Element.hh:119
G4double GetN() const
Definition G4Element.hh:123
const G4Material * GetMaterial() const
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
void SetLPMFlag(G4bool)
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
void SetMasterThread(G4bool)
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
G4int SelectIsotopeNumber(const G4Element *) const
void SetPolarAngleLimit(G4double)
G4VEmModel(const G4VEmModel &)=delete
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
G4PhysicsTable * xSectionTable
void SetHighEnergyLimit(G4double)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const std::vector< G4double > * theDensityFactor
G4VEmFluctuationModel * GetModelOfFluctuations()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
std::size_t currentCoupleIndex
void SetActivationLowEnergyLimit(G4double)
G4double PolarAngleLimit() const
G4VEmAngularDistribution * GetAngularDistribution()
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
void SetForceBuildTable(G4bool val)
G4double inveplus
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4int SelectRandomAtomNumber(const G4Material *) const
G4bool IsMaster() const
void SetSecondaryThreshold(G4double)
G4VEmModel * GetTripletModel()
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4bool lossFlucFlag
void SetAngularGeneratorFlag(G4bool)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
std::size_t basedCoupleIndex
void SetFluctuationFlag(G4bool val)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
G4double HighEnergyLimit() const
virtual void CorrectionsAlongStep(const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss)
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
void SetCurrentElement(const G4Element *)
G4VParticleChange * pParticleChange
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const G4Material * pBaseMaterial
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
virtual void DefineForRegion(const G4Region *)
virtual void ModelDescription(std::ostream &outFile) const
void SetLowEnergyLimit(G4double)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
void SetActivationHighEnergyLimit(G4double)
G4double pFactor
void SetTripletModel(G4VEmModel *)
G4ElementData * fElementData
G4double HighEnergyActivationLimit() const
void SetLocked(G4bool)
void SetDeexcitationFlag(G4bool val)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool IsActive(G4double kinEnergy) const
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const G4MaterialCutsCouple * CurrentCouple() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
const std::vector< G4int > * theDensityIdx
G4bool DeexcitationFlag() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4bool UseAngularGeneratorFlag() const
const G4String & GetName() const
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
G4PhysicsTable * GetCrossSectionTable()
void SetUseBaseMaterials(G4bool val)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
virtual ~G4VEmModel()
Definition G4VEmModel.cc:86
G4bool UseBaseMaterials() const
G4double LowEnergyActivationLimit() const
G4bool IsLocked() const
virtual void StartTracking(G4Track *)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
G4bool ForceBuildTableFlag() const
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
G4double SecondaryThreshold() const
G4VEmModel & operator=(const G4VEmModel &right)=delete
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double ChargeSquareRatio(const G4Track &)
G4ElementData * GetElementData()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
#define DBL_MAX
Definition templates.hh:62