Geant4 11.1.1
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 LPMFlag() 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 SetLPMFlag(G4bool val);
371
372 inline void SetDeexcitationFlag(G4bool val);
373
374 inline void SetForceBuildTable(G4bool val);
375
376 inline void SetFluctuationFlag(G4bool val);
377
378 inline void SetMasterThread(G4bool val);
379
380 inline G4bool IsMaster() const;
381
382 inline void SetUseBaseMaterials(G4bool val);
383
384 inline G4bool UseBaseMaterials() const;
385
386 inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
387
388 inline const G4String& GetName() const;
389
390 inline void SetCurrentCouple(const G4MaterialCutsCouple*);
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 = 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
447
448private:
449
450 G4bool theLPMflag = false;
451 G4bool flagDeexcitation = false;
452 G4bool flagForceBuildTable = false;
453 G4bool isMaster = true;
454
455 G4bool localTable = true;
456 G4bool localElmSelectors = true;
457 G4bool useAngularGenerator = false;
458 G4bool useBaseMaterials = false;
459 G4bool isLocked = false;
460
461 const G4String name;
462 std::vector<G4double> xsec;
463
464};
465
466// ======== Run time inline methods ================
467
469{
470 if(fCurrentCouple != ptr) {
471 fCurrentCouple = ptr;
473 pBaseMaterial = ptr->GetMaterial();
474 pFactor = 1.0;
475 if(useBaseMaterials) {
476 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
477 if(nullptr != pBaseMaterial->GetBaseMaterial())
479 pFactor = (*theDensityFactor)[currentCoupleIndex];
480 }
481 }
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
485
487{
488 return fCurrentCouple;
489}
490
491//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
492
494{
495 fCurrentElement = elm;
496}
497
498//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
499
500inline
502{
504 dynPart->GetKineticEnergy());
505}
506
507//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
508
510 const G4ParticleDefinition* part,
511 G4double kinEnergy,
512 G4double cutEnergy)
513{
514 SetCurrentCouple(couple);
515 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
516}
517
518//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
519
521 const G4ParticleDefinition* part,
522 G4double kinEnergy,
523 G4double cutEnergy,
524 G4double maxEnergy)
525{
526 SetCurrentCouple(couple);
527 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
528 cutEnergy,maxEnergy);
529}
530
531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
532
533inline
535 G4double ekin,
536 const G4Material* material,
537 G4double emin,
538 G4double emax)
539{
540 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
541 return (cross > 0.0) ? 1./cross : DBL_MAX;
542}
543
544//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
545
546inline G4double
548 const G4Element* elm,
549 G4double kinEnergy,
550 G4double cutEnergy,
551 G4double maxEnergy)
552{
553 fCurrentElement = elm;
554 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
555 cutEnergy,maxEnergy);
556}
557
558//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
559
560inline const G4Element*
562 const G4ParticleDefinition* part,
563 G4double kinEnergy,
564 G4double cutEnergy,
565 G4double maxEnergy)
566{
567 SetCurrentCouple(couple);
568 fCurrentElement = (nSelectors > 0) ?
569 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
570 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
571 return fCurrentElement;
572}
573
574//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
575
576inline const G4Element*
578 const G4ParticleDefinition* part,
579 G4double kinEnergy,
580 G4double logKinE,
581 G4double cutEnergy,
582 G4double maxEnergy)
583{
584 SetCurrentCouple(couple);
585 fCurrentElement = (nSelectors > 0)
586 ? ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy,logKinE)
587 : SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
588 return fCurrentElement;
589}
590
591// ======== Get/Set inline methods used at initialisation ================
592
594{
595 return flucModel;
596}
597
598//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
599
601{
602 return anglModel;
603}
604
605//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
606
608{
609 if(p != anglModel) {
610 delete anglModel;
611 anglModel = p;
612 }
613}
614
615//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
616
618{
619 return fTripletModel;
620}
621
622//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
623
625{
626 if(p != fTripletModel) {
627 delete fTripletModel;
628 fTripletModel = p;
629 }
630}
631
632//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
633
635{
636 return highLimit;
637}
638
639//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
640
642{
643 return lowLimit;
644}
645
646//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
647
649{
650 return eMaxActive;
651}
652
653//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
654
656{
657 return eMinActive;
658}
659
660//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
661
663{
664 return polarAngleLimit;
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668
670{
671 return secondaryThreshold;
672}
673
674//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
675
677{
678 return theLPMflag;
679}
680
681//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682
684{
685 return flagDeexcitation;
686}
687
688//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
689
691{
692 return flagForceBuildTable;
693}
694
695//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
696
698{
699 return useAngularGenerator;
700}
701
702//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
703
705{
706 useAngularGenerator = val;
707}
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710
712{
713 lossFlucFlag = val;
714}
715
716//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
717
719{
720 isMaster = val;
721}
722
723//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724
726{
727 return isMaster;
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
733{
734 useBaseMaterials = val;
735}
736
737//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738
740{
741 return useBaseMaterials;
742}
743
744//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745
747{
748 highLimit = val;
749}
750
751//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752
754{
755 lowLimit = val;
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
761{
762 eMaxActive = val;
763}
764
765//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766
768{
769 eMinActive = val;
770}
771
772//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773
774inline G4bool G4VEmModel::IsActive(G4double kinEnergy) const
775{
776 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
777}
778
779//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780
782{
783 if(!isLocked) { polarAngleLimit = val; }
784}
785
786//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787
789{
790 secondaryThreshold = val;
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
796{
797 theLPMflag = val;
798}
799
800//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801
803{
804 flagDeexcitation = val;
805}
806
807//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808
810{
811 flagForceBuildTable = val;
812}
813
814//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
815
816inline const G4String& G4VEmModel::GetName() const
817{
818 return name;
819}
820
821//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
822
823inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
824{
825 return elmSelectors;
826}
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
829
830inline void
831G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
832{
833 if(p != elmSelectors) {
834 elmSelectors = p;
835 nSelectors = (nullptr != elmSelectors) ? G4int(elmSelectors->size()) : 0;
836 localElmSelectors = false;
837 }
838}
839
840//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
841
843{
844 return fElementData;
845}
846
847//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
848
850{
851 return xSectionTable;
852}
853
854//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
855
857{
858 return isLocked;
859}
860
861//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
862
864{
865 isLocked = val;
866}
867
868//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
869
870#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:135
const G4Material * GetMaterial() const
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:228
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:308
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:398
G4int SelectIsotopeNumber(const G4Element *) const
Definition: G4VEmModel.cc:276
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:781
G4VEmModel(const G4VEmModel &)=delete
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:358
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:422
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:284
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:424
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:718
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:593
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:334
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:153
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:809
G4double inveplus
Definition: G4VEmModel.hh:427
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4int SelectRandomAtomNumber(const G4Material *) const
Definition: G4VEmModel.cc:253
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:788
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:617
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:181
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:534
G4bool lossFlucFlag
Definition: G4VEmModel.hh:446
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:704
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:390
void SetFluctuationFlag(G4bool val)
Definition: G4VEmModel.hh:711
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:795
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:167
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:509
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:493
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:421
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:423
G4bool LPMFlag() const
Definition: G4VEmModel.hh:676
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
Definition: G4VEmModel.cc:261
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:303
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:412
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:383
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:760
G4double pFactor
Definition: G4VEmModel.hh:428
void SetTripletModel(G4VEmModel *)
Definition: G4VEmModel.hh:624
G4ElementData * fElementData
Definition: G4VEmModel.hh:420
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:648
void SetLocked(G4bool)
Definition: G4VEmModel.hh:863
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:520
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:774
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:294
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:342
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:425
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:683
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:325
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
const G4String & GetName() const
Definition: G4VEmModel.hh:816
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
Definition: G4VEmModel.cc:242
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:205
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:172
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:849
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:732
size_t currentCoupleIndex
Definition: G4VEmModel.hh:444
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:85
G4bool UseBaseMaterials() const
Definition: G4VEmModel.hh:739
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:655
G4bool IsLocked() const
Definition: G4VEmModel.hh:856
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:214
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:148
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:367
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:690
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:501
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
G4VEmModel & operator=(const G4VEmModel &right)=delete
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:375
size_t basedCoupleIndex
Definition: G4VEmModel.hh:445
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:317
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:842
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
#define DBL_MAX
Definition: templates.hh:62