Geant4 10.7.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 G4MaterialCutsCouple*,
195 const G4DynamicParticle*,
196 G4double& eloss,
197 G4double& niel,
198 G4double length);
199
200 // value which may be tabulated (by default cross section)
201 virtual G4double Value(const G4MaterialCutsCouple*,
203 G4double kineticEnergy);
204
205 // threshold for zero value
206 virtual G4double MinPrimaryEnergy(const G4Material*,
208 G4double cut = 0.0);
209
210 // model can define low-energy limit for the cut
212 const G4MaterialCutsCouple*);
213
214 // initialisation at run time for a given material
215 virtual void SetupForMaterial(const G4ParticleDefinition*,
216 const G4Material*,
217 G4double kineticEnergy);
218
219 // add a region for the model
220 virtual void DefineForRegion(const G4Region*);
221
222 // for automatic documentation
223 virtual void ModelDescription(std::ostream& outFile) const;
224
225protected:
226
227 // initialisation of the ParticleChange for the model
229
230 // initialisation of the ParticleChange for the model
232
233 // kinematically allowed max kinetic energy of a secondary
235 G4double kineticEnergy);
236
237public:
238
239 //------------------------------------------------------------------------
240 // Generic methods common to all models
241 //------------------------------------------------------------------------
242
243 // should be called at initialisation to build element selectors
245 const G4DataVector&);
246
247 // should be called at initialisation to access element selectors
248 inline std::vector<G4EmElementSelector*>* GetElementSelectors();
249
250 // should be called at initialisation to set element selectors
251 inline void SetElementSelectors(std::vector<G4EmElementSelector*>*);
252
253 // dEdx per unit length, base material approach may be used
254 virtual inline G4double ComputeDEDX(const G4MaterialCutsCouple*,
256 G4double kineticEnergy,
257 G4double cutEnergy = DBL_MAX);
258
259 // cross section per volume, base material approach may be used
262 G4double kineticEnergy,
263 G4double cutEnergy = 0.0,
264 G4double maxEnergy = DBL_MAX);
265
266 // compute mean free path via cross section per volume
268 G4double kineticEnergy,
269 const G4Material*,
270 G4double cutEnergy = 0.0,
271 G4double maxEnergy = DBL_MAX);
272
273 // generic cross section per element
275 const G4Element*,
276 G4double kinEnergy,
277 G4double cutEnergy = 0.0,
278 G4double maxEnergy = DBL_MAX);
279
280 // atom can be selected effitiantly if element selectors are initialised
283 G4double kineticEnergy,
284 G4double cutEnergy = 0.0,
285 G4double maxEnergy = DBL_MAX);
286 // same as SelectRandomAtom above but more efficient since log-ekin is known
289 G4double kineticEnergy,
290 G4double logKineticEnergy,
291 G4double cutEnergy = 0.0,
292 G4double maxEnergy = DBL_MAX);
293
294
295 // to select atom cross section per volume is recomputed for each element
298 G4double kineticEnergy,
299 G4double cutEnergy = 0.0,
300 G4double maxEnergy = DBL_MAX);
301
302 // to select atom if cross section is proportional number of electrons
304
305 // select isotope in order to have precise mass of the nucleus
307
308 //------------------------------------------------------------------------
309 // Get/Set methods
310 //------------------------------------------------------------------------
311
313
315
317
319
321
323
324 inline G4VEmModel* GetTripletModel();
325
326 inline void SetTripletModel(G4VEmModel*);
327
329
330 inline G4double HighEnergyLimit() const;
331
332 inline G4double LowEnergyLimit() const;
333
334 inline G4double HighEnergyActivationLimit() const;
335
336 inline G4double LowEnergyActivationLimit() const;
337
338 inline G4double PolarAngleLimit() const;
339
340 inline G4double SecondaryThreshold() const;
341
342 inline G4bool LPMFlag() const;
343
344 inline G4bool DeexcitationFlag() const;
345
346 inline G4bool ForceBuildTableFlag() const;
347
348 inline G4bool UseAngularGeneratorFlag() const;
349
350 inline void SetAngularGeneratorFlag(G4bool);
351
352 inline void SetHighEnergyLimit(G4double);
353
354 inline void SetLowEnergyLimit(G4double);
355
357
359
360 inline G4bool IsActive(G4double kinEnergy) const;
361
362 inline void SetPolarAngleLimit(G4double);
363
364 inline void SetSecondaryThreshold(G4double);
365
366 inline void SetLPMFlag(G4bool val);
367
368 inline void SetDeexcitationFlag(G4bool val);
369
370 inline void SetForceBuildTable(G4bool val);
371
372 inline void SetFluctuationFlag(G4bool val);
373
374 inline void SetMasterThread(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 const G4Element* GetCurrentElement() const;
389
390 inline const G4Isotope* GetCurrentIsotope() const;
391
392 inline G4bool IsLocked() const;
393
394 inline void SetLocked(G4bool);
395
396 // hide assignment operator
397 G4VEmModel & operator=(const G4VEmModel &right) = delete;
398 G4VEmModel(const G4VEmModel&) = delete;
399
400protected:
401
402 inline const G4MaterialCutsCouple* CurrentCouple() const;
403
404 inline void SetCurrentElement(const G4Element*);
405
406private:
407
408 // ======== Parameters of the class fixed at construction =========
409
410 G4VEmFluctuationModel* flucModel;
411 G4VEmAngularDistribution* anglModel;
412 const G4String name;
413
414 // ======== Parameters of the class fixed at initialisation =======
415
416 G4double lowLimit;
417 G4double highLimit;
418 G4double eMinActive;
419 G4double eMaxActive;
420 G4double polarAngleLimit;
421 G4double secondaryThreshold;
422 G4bool theLPMflag;
423 G4bool flagDeexcitation;
424 G4bool flagForceBuildTable;
425 G4bool isMaster;
426
427 G4bool localTable;
428 G4bool localElmSelectors;
429 G4bool useAngularGenerator;
430 G4bool useBaseMaterials;
431 G4bool isLocked;
432 G4int nSelectors;
433 std::vector<G4EmElementSelector*>* elmSelectors;
434 G4LossTableManager* fEmManager;
435
436protected:
437
442 const std::vector<G4double>* theDensityFactor;
443 const std::vector<G4int>* theDensityIdx;
444 size_t idxTable;
448
449 // ======== Cached values - may be state dependent ================
450
451private:
452
453 const G4MaterialCutsCouple* fCurrentCouple;
454 const G4Element* fCurrentElement;
455 const G4Isotope* fCurrentIsotope;
456 G4VEmModel* fTripletModel;
457
458 G4int nsec;
459 std::vector<G4double> xsec;
460
461};
462
463// ======== Run time inline methods ================
464
466{
467 if(fCurrentCouple != ptr) {
468 fCurrentCouple = ptr;
469 pBaseMaterial = ptr->GetMaterial();
470 pFactor = 1.0;
471 if(useBaseMaterials && pBaseMaterial->GetBaseMaterial()) {
473 pFactor = (*theDensityFactor)[(*theDensityIdx)[ptr->GetIndex()]];
474 }
475 }
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479
481{
482 return fCurrentCouple;
483}
484
485//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
486
488{
489 fCurrentElement = elm;
490 fCurrentIsotope = nullptr;
491}
492
493//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
494
496{
497 return fCurrentElement;
498}
499
500//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
501
503{
504 return fCurrentIsotope;
505}
506
507//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
508
509inline
511{
513 dynPart->GetKineticEnergy());
514}
515
516//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
517
519 const G4ParticleDefinition* part,
520 G4double kinEnergy,
521 G4double cutEnergy)
522{
523 SetCurrentCouple(couple);
524 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
525}
526
527//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528
530 const G4ParticleDefinition* part,
531 G4double kinEnergy,
532 G4double cutEnergy,
533 G4double maxEnergy)
534{
535 SetCurrentCouple(couple);
536 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
537 cutEnergy,maxEnergy);
538}
539
540//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
541
542inline
544 G4double ekin,
545 const G4Material* material,
546 G4double emin,
547 G4double emax)
548{
549 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
550 return (cross > 0.0) ? 1./cross : DBL_MAX;
551}
552
553//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
554
555inline G4double
557 const G4Element* elm,
558 G4double kinEnergy,
559 G4double cutEnergy,
560 G4double maxEnergy)
561{
563 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
564 cutEnergy,maxEnergy);
565}
566
567//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
568
569inline const G4Element*
571 const G4ParticleDefinition* part,
572 G4double kinEnergy,
573 G4double cutEnergy,
574 G4double maxEnergy)
575{
576 SetCurrentCouple(couple);
577 fCurrentElement = (nSelectors > 0) ?
578 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
579 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
580 fCurrentIsotope = nullptr;
581 return fCurrentElement;
582}
583
584//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
585
586inline const G4Element*
588 const G4ParticleDefinition* part,
589 G4double kinEnergy,
590 G4double logKinE,
591 G4double cutEnergy,
592 G4double maxEnergy)
593{
594 SetCurrentCouple(couple);
595 fCurrentElement = (nSelectors > 0)
596 ? ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy,logKinE)
597 : SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
598 fCurrentIsotope = nullptr;
599 return fCurrentElement;
600}
601
602// ======== Get/Set inline methods used at initialisation ================
603
605{
606 return flucModel;
607}
608
609//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
610
612{
613 return anglModel;
614}
615
616//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
617
619{
620 if(p != anglModel) {
621 delete anglModel;
622 anglModel = p;
623 }
624}
625
626//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627
629{
630 return fTripletModel;
631}
632
633//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
634
636{
637 if(p != fTripletModel) {
638 delete fTripletModel;
639 fTripletModel = p;
640 }
641}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
644
646{
647 return highLimit;
648}
649
650//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
651
653{
654 return lowLimit;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
658
660{
661 return eMaxActive;
662}
663
664//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
665
667{
668 return eMinActive;
669}
670
671//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
672
674{
675 return polarAngleLimit;
676}
677
678//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
679
681{
682 return secondaryThreshold;
683}
684
685//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
686
688{
689 return theLPMflag;
690}
691
692//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
693
695{
696 return flagDeexcitation;
697}
698
699//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
700
702{
703 return flagForceBuildTable;
704}
705
706//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
707
709{
710 return useAngularGenerator;
711}
712
713//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
714
716{
717 useAngularGenerator = val;
718}
719
720//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
721
723{
724 lossFlucFlag = val;
725}
726
727//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
728
730{
731 isMaster = val;
732}
733
734//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
735
737{
738 return isMaster;
739}
740
741//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
742
744{
745 useBaseMaterials = val;
746}
747
748//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
749
751{
752 return useBaseMaterials;
753}
754
755//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
756
758{
759 highLimit = val;
760}
761
762//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
763
765{
766 lowLimit = val;
767}
768
769//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
770
772{
773 eMaxActive = val;
774}
775
776//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
777
779{
780 eMinActive = val;
781}
782
783//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
784
785inline G4bool G4VEmModel::IsActive(G4double kinEnergy) const
786{
787 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
788}
789
790//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
791
793{
794 if(!isLocked) { polarAngleLimit = val; }
795}
796
797//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
798
800{
801 secondaryThreshold = val;
802}
803
804//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
805
807{
808 theLPMflag = val;
809}
810
811//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
812
814{
815 flagDeexcitation = val;
816}
817
818//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
819
821{
822 flagForceBuildTable = val;
823}
824
825//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
826
827inline const G4String& G4VEmModel::GetName() const
828{
829 return name;
830}
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
833
834inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
835{
836 return elmSelectors;
837}
838
839//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
840
841inline void
842G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
843{
844 if(p != elmSelectors) {
845 elmSelectors = p;
846 nSelectors = (elmSelectors) ? G4int(elmSelectors->size()) : 0;
847 localElmSelectors = false;
848 }
849}
850
851//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
852
854{
855 return fElementData;
856}
857
858//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
859
861{
862 return xSectionTable;
863}
864
865//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
866
868{
869 return isLocked;
870}
871
872//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
873
875{
876 isLocked = val;
877}
878
879//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
880
881#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:130
G4double GetN() const
Definition: G4Element.hh:134
const G4Material * GetMaterial() const
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:231
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:464
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:792
G4VEmModel(const G4VEmModel &)=delete
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:424
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:440
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:359
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:408
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:442
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:729
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:400
size_t idxTable
Definition: G4VEmModel.hh:444
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:673
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:226
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:820
G4double inveplus
Definition: G4VEmModel.hh:446
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:799
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:628
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
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:543
G4bool lossFlucFlag
Definition: G4VEmModel.hh:445
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:495
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:715
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456
void SetFluctuationFlag(G4bool val)
Definition: G4VEmModel.hh:722
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:806
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:240
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:315
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:415
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
virtual G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:518
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:487
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:439
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:441
G4bool LPMFlag() const
Definition: G4VEmModel.hh:687
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:378
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:478
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:337
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:449
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:771
G4double pFactor
Definition: G4VEmModel.hh:447
void SetTripletModel(G4VEmModel *)
Definition: G4VEmModel.hh:635
G4ElementData * fElementData
Definition: G4VEmModel.hh:438
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:659
void SetLocked(G4bool)
Definition: G4VEmModel.hh:874
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:529
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:369
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:480
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:443
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:694
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
const G4String & GetName() const
Definition: G4VEmModel.hh:827
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:279
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:860
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:743
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:502
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:94
G4bool UseBaseMaterials() const
Definition: G4VEmModel.hh:750
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:666
G4bool IsLocked() const
Definition: G4VEmModel.hh:867
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:288
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:220
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:433
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:701
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
G4VEmModel & operator=(const G4VEmModel &right)=delete
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:441
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:383
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:853
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
#define DBL_MAX
Definition: templates.hh:62