Geant4 11.2.2
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 const G4double& length,
197 G4double& eloss);
198
199 // value which may be tabulated (by default cross section)
200 virtual G4double Value(const G4MaterialCutsCouple*,
202 G4double kineticEnergy);
203
204 // threshold for zero value
205 virtual G4double MinPrimaryEnergy(const G4Material*,
207 G4double cut = 0.0);
208
209 // model can define low-energy limit for the cut
211 const G4MaterialCutsCouple*);
212
213 // initialisation at run time for a given material
214 virtual void SetupForMaterial(const G4ParticleDefinition*,
215 const G4Material*,
216 G4double kineticEnergy);
217
218 // add a region for the model
219 virtual void DefineForRegion(const G4Region*);
220
221 // fill number of different type of secondaries after SampleSecondaries(...)
222 virtual void FillNumberOfSecondaries(G4int& numberOfTriplets,
223 G4int& numberOfRecoil);
224
225 // for automatic documentation
226 virtual void ModelDescription(std::ostream& outFile) const;
227
228protected:
229
230 // initialisation of the ParticleChange for the model
232
233 // initialisation of the ParticleChange for the model
235
236 // kinematically allowed max kinetic energy of a secondary
238 G4double kineticEnergy);
239
240public:
241
242 //------------------------------------------------------------------------
243 // Generic methods common to all models
244 //------------------------------------------------------------------------
245
246 // should be called at initialisation to build element selectors
248 const G4DataVector&);
249
250 // should be called at initialisation to access element selectors
251 inline std::vector<G4EmElementSelector*>* GetElementSelectors();
252
253 // should be called at initialisation to set element selectors
254 inline void SetElementSelectors(std::vector<G4EmElementSelector*>*);
255
256 // dEdx per unit length, base material approach may be used
259 G4double kineticEnergy,
260 G4double cutEnergy = DBL_MAX);
261
262 // cross section per volume, base material approach may be used
265 G4double kineticEnergy,
266 G4double cutEnergy = 0.0,
267 G4double maxEnergy = DBL_MAX);
268
269 // compute mean free path via cross section per volume
271 G4double kineticEnergy,
272 const G4Material*,
273 G4double cutEnergy = 0.0,
274 G4double maxEnergy = DBL_MAX);
275
276 // generic cross section per element
278 const G4Element*,
279 G4double kinEnergy,
280 G4double cutEnergy = 0.0,
281 G4double maxEnergy = DBL_MAX);
282
283 // atom can be selected effitiantly if element selectors are initialised
286 G4double kineticEnergy,
287 G4double cutEnergy = 0.0,
288 G4double maxEnergy = DBL_MAX);
289 // same as SelectRandomAtom above but more efficient since log-ekin is known
292 G4double kineticEnergy,
293 G4double logKineticEnergy,
294 G4double cutEnergy = 0.0,
295 G4double maxEnergy = DBL_MAX);
296
297 // to select atom cross section per volume is recomputed for each element
300 G4double kineticEnergy,
301 G4double cutEnergy = 0.0,
302 G4double maxEnergy = DBL_MAX);
303
304 // to select atom if cross section is proportional number of electrons
305 const G4Element* GetCurrentElement(const G4Material* mat = nullptr) const;
307
308 // select isotope in order to have precise mass of the nucleus
309 const G4Isotope* GetCurrentIsotope(const G4Element* elm = nullptr) const;
310 G4int SelectIsotopeNumber(const G4Element*) const;
311
312 //------------------------------------------------------------------------
313 // Get/Set methods
314 //------------------------------------------------------------------------
315
317
319
321
323
325
327
328 inline G4VEmModel* GetTripletModel();
329
330 inline void SetTripletModel(G4VEmModel*);
331
333
334 inline G4double HighEnergyLimit() const;
335
336 inline G4double LowEnergyLimit() const;
337
338 inline G4double HighEnergyActivationLimit() const;
339
340 inline G4double LowEnergyActivationLimit() const;
341
342 inline G4double PolarAngleLimit() const;
343
344 inline G4double SecondaryThreshold() const;
345
346 inline G4bool DeexcitationFlag() const;
347
348 inline G4bool ForceBuildTableFlag() const;
349
350 inline G4bool UseAngularGeneratorFlag() const;
351
352 inline void SetAngularGeneratorFlag(G4bool);
353
354 inline void SetHighEnergyLimit(G4double);
355
356 inline void SetLowEnergyLimit(G4double);
357
359
361
362 inline G4bool IsActive(G4double kinEnergy) const;
363
364 inline void SetPolarAngleLimit(G4double);
365
366 inline void SetSecondaryThreshold(G4double);
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 G4bool IsLocked() const;
389
390 inline void SetLocked(G4bool);
391
392 // obsolete method
393 [[deprecated("Use G4EmParameters::Instance()->SetLPM instead")]]
394 void SetLPMFlag(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 = nullptr;
411 G4VEmAngularDistribution* anglModel = nullptr;
412 G4VEmModel* fTripletModel = nullptr;
413 const G4MaterialCutsCouple* fCurrentCouple = nullptr;
414 const G4Element* fCurrentElement = nullptr;
415 std::vector<G4EmElementSelector*>* elmSelectors = nullptr;
416 G4LossTableManager* fEmManager;
417
418protected:
419
423 const G4Material* pBaseMaterial = nullptr;
424 const std::vector<G4double>* theDensityFactor = nullptr;
425 const std::vector<G4int>* theDensityIdx = nullptr;
426
429
430private:
431
432 G4double lowLimit;
433 G4double highLimit;
434 G4double eMinActive = 0.0;
435 G4double eMaxActive = DBL_MAX;
436 G4double secondaryThreshold = DBL_MAX;
437 G4double polarAngleLimit;
438
439 G4int nSelectors = 0;
440 G4int nsec = 5;
441
442protected:
443
444 std::size_t currentCoupleIndex = 0;
445 std::size_t basedCoupleIndex = 0;
447
448private:
449
450 G4bool flagDeexcitation = false;
451 G4bool flagForceBuildTable = false;
452 G4bool isMaster = true;
453
454 G4bool localTable = true;
455 G4bool localElmSelectors = true;
456 G4bool useAngularGenerator = false;
457 G4bool useBaseMaterials = false;
458 G4bool isLocked = false;
459
460 const G4String name;
461 std::vector<G4double> xsec;
462
463};
464
465// ======== Run time inline methods ================
466
468{
469 if(fCurrentCouple != ptr) {
470 fCurrentCouple = ptr;
472 pBaseMaterial = ptr->GetMaterial();
473 pFactor = 1.0;
474 if(useBaseMaterials) {
475 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
476 if(nullptr != pBaseMaterial->GetBaseMaterial())
478 pFactor = (*theDensityFactor)[currentCoupleIndex];
479 }
480 }
481}
482
483//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
484
486{
487 return fCurrentCouple;
488}
489
490//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
491
493{
494 fCurrentElement = elm;
495}
496
497//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
498
499inline
505
506//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
507
509 const G4ParticleDefinition* part,
510 G4double kinEnergy,
511 G4double cutEnergy)
512{
513 SetCurrentCouple(couple);
514 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
515}
516
517//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
518
520 const G4ParticleDefinition* part,
521 G4double kinEnergy,
522 G4double cutEnergy,
523 G4double maxEnergy)
524{
525 SetCurrentCouple(couple);
526 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
527 cutEnergy,maxEnergy);
528}
529
530//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
531
532inline
534 G4double ekin,
535 const G4Material* material,
536 G4double emin,
537 G4double emax)
538{
539 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
540 return (cross > 0.0) ? 1./cross : DBL_MAX;
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
544
545inline G4double
547 const G4Element* elm,
548 G4double kinEnergy,
549 G4double cutEnergy,
550 G4double maxEnergy)
551{
552 fCurrentElement = elm;
553 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
554 cutEnergy,maxEnergy);
555}
556
557//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
558
559inline const G4Element*
561 const G4ParticleDefinition* part,
562 G4double kinEnergy,
563 G4double cutEnergy,
564 G4double maxEnergy)
565{
566 SetCurrentCouple(couple);
567 fCurrentElement = (nSelectors > 0) ?
568 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
569 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
570 return fCurrentElement;
571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
574
575inline const G4Element*
577 const G4ParticleDefinition* part,
578 G4double kinEnergy,
579 G4double logKinE,
580 G4double cutEnergy,
581 G4double maxEnergy)
582{
583 SetCurrentCouple(couple);
584 fCurrentElement = (nSelectors > 0)
585 ? ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy,logKinE)
586 : SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
587 return fCurrentElement;
588}
589
590// ======== Get/Set inline methods used at initialisation ================
591
593{
594 return flucModel;
595}
596
597//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
598
600{
601 return anglModel;
602}
603
604//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
605
607{
608 if(p != anglModel) {
609 delete anglModel;
610 anglModel = p;
611 }
612}
613
614//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
615
617{
618 return fTripletModel;
619}
620
621//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
622
624{
625 if(p != fTripletModel) {
626 delete fTripletModel;
627 fTripletModel = p;
628 }
629}
630
631//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
632
634{
635 return highLimit;
636}
637
638//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
639
641{
642 return lowLimit;
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646
648{
649 return eMaxActive;
650}
651
652//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
653
655{
656 return eMinActive;
657}
658
659//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
660
662{
663 return polarAngleLimit;
664}
665
666//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
667
669{
670 return secondaryThreshold;
671}
672
673//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
674
676{
677 return flagDeexcitation;
678}
679
680//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
681
683{
684 return flagForceBuildTable;
685}
686
687//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
688
690{
691 return useAngularGenerator;
692}
693
694//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
695
697{
698 useAngularGenerator = val;
699}
700
701//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
702
704{
705 lossFlucFlag = val;
706}
707
708//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
709
711{
712 isMaster = val;
713}
714
715//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
716
718{
719 return isMaster;
720}
721
722//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
723
725{
726 useBaseMaterials = val;
727}
728
729//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
730
732{
733 return useBaseMaterials;
734}
735
736//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
737
739{
740 highLimit = val;
741}
742
743//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
744
746{
747 lowLimit = val;
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
751
753{
754 eMaxActive = val;
755}
756
757//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
758
760{
761 eMinActive = val;
762}
763
764//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
765
766inline G4bool G4VEmModel::IsActive(G4double kinEnergy) const
767{
768 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
769}
770
771//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
772
774{
775 if(!isLocked) { polarAngleLimit = val; }
776}
777
778//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
779
781{
782 secondaryThreshold = val;
783}
784
785//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
786
788{
789 flagDeexcitation = val;
790}
791
792//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
793
795{
796 flagForceBuildTable = val;
797}
798
799//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
800
801inline const G4String& G4VEmModel::GetName() const
802{
803 return name;
804}
805
806//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
807
808inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
809{
810 return elmSelectors;
811}
812
813//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
814
815inline void
816G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
817{
818 if(p != elmSelectors) {
819 elmSelectors = p;
820 nSelectors = (nullptr != elmSelectors) ? G4int(elmSelectors->size()) : 0;
821 localElmSelectors = false;
822 }
823}
824
825//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
826
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
833
838
839//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
840
842{
843 return isLocked;
844}
845
846//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
847
849{
850 isLocked = val;
851}
852
853//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
854
855#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
const G4Material * GetBaseMaterial() const
void SetLPMFlag(G4bool)
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
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
void SetMasterThread(G4bool val)
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
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 CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
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:85
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