Geant4-11
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
298 // to select atom cross section per volume is recomputed for each element
301 G4double kineticEnergy,
302 G4double cutEnergy = 0.0,
303 G4double maxEnergy = DBL_MAX);
304
305 // to select atom if cross section is proportional number of electrons
307
308 // select isotope in order to have precise mass of the nucleus
310
311 //------------------------------------------------------------------------
312 // Get/Set methods
313 //------------------------------------------------------------------------
314
316
318
320
322
324
326
327 inline G4VEmModel* GetTripletModel();
328
329 inline void SetTripletModel(G4VEmModel*);
330
332
333 inline G4double HighEnergyLimit() const;
334
335 inline G4double LowEnergyLimit() const;
336
337 inline G4double HighEnergyActivationLimit() const;
338
339 inline G4double LowEnergyActivationLimit() const;
340
341 inline G4double PolarAngleLimit() const;
342
343 inline G4double SecondaryThreshold() const;
344
345 inline G4bool LPMFlag() const;
346
347 inline G4bool DeexcitationFlag() const;
348
349 inline G4bool ForceBuildTableFlag() const;
350
351 inline G4bool UseAngularGeneratorFlag() const;
352
353 inline void SetAngularGeneratorFlag(G4bool);
354
355 inline void SetHighEnergyLimit(G4double);
356
357 inline void SetLowEnergyLimit(G4double);
358
360
362
363 inline G4bool IsActive(G4double kinEnergy) const;
364
365 inline void SetPolarAngleLimit(G4double);
366
367 inline void SetSecondaryThreshold(G4double);
368
369 inline void SetLPMFlag(G4bool val);
370
371 inline void SetDeexcitationFlag(G4bool val);
372
373 inline void SetForceBuildTable(G4bool val);
374
375 inline void SetFluctuationFlag(G4bool val);
376
377 inline void SetMasterThread(G4bool val);
378
379 inline G4bool IsMaster() const;
380
381 inline void SetUseBaseMaterials(G4bool val);
382
383 inline G4bool UseBaseMaterials() const;
384
385 inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
386
387 inline const G4String& GetName() const;
388
389 inline void SetCurrentCouple(const G4MaterialCutsCouple*);
390
391 inline const G4Element* GetCurrentElement() const;
392
393 inline const G4Isotope* GetCurrentIsotope() const;
394
395 inline G4bool IsLocked() const;
396
397 inline void SetLocked(G4bool);
398
399 // hide assignment operator
400 G4VEmModel & operator=(const G4VEmModel &right) = delete;
401 G4VEmModel(const G4VEmModel&) = delete;
402
403protected:
404
405 inline const G4MaterialCutsCouple* CurrentCouple() const;
406
407 inline void SetCurrentElement(const G4Element*);
408
409private:
410
411 // ======== Parameters of the class fixed at construction =========
412
417 const G4Element* fCurrentElement = nullptr;
418 const G4Isotope* fCurrentIsotope = nullptr;
419 std::vector<G4EmElementSelector*>* elmSelectors = nullptr;
421
422protected:
423
427 const G4Material* pBaseMaterial = nullptr;
428 const std::vector<G4double>* theDensityFactor = nullptr;
429 const std::vector<G4int>* theDensityIdx = nullptr;
430
433
434private:
435
442
445
446protected:
447
451
452private:
453
458
464
466 std::vector<G4double> xsec;
467
468};
469
470// ======== Run time inline methods ================
471
473{
474 if(fCurrentCouple != ptr) {
475 fCurrentCouple = ptr;
477 pBaseMaterial = ptr->GetMaterial();
478 pFactor = 1.0;
479 if(useBaseMaterials) {
480 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
481 if(nullptr != pBaseMaterial->GetBaseMaterial())
483 pFactor = (*theDensityFactor)[currentCoupleIndex];
484 }
485 }
486}
487
488//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
489
491{
492 return fCurrentCouple;
493}
494
495//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
496
498{
499 fCurrentElement = elm;
500 fCurrentIsotope = nullptr;
501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
506{
507 return fCurrentElement;
508}
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
511
513{
514 return fCurrentIsotope;
515}
516
517//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
518
519inline
521{
523 dynPart->GetKineticEnergy());
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
527
529 const G4ParticleDefinition* part,
530 G4double kinEnergy,
531 G4double cutEnergy)
532{
533 SetCurrentCouple(couple);
534 return pFactor*ComputeDEDXPerVolume(pBaseMaterial,part,kinEnergy,cutEnergy);
535}
536
537//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
538
540 const G4ParticleDefinition* part,
541 G4double kinEnergy,
542 G4double cutEnergy,
543 G4double maxEnergy)
544{
545 SetCurrentCouple(couple);
546 return pFactor*CrossSectionPerVolume(pBaseMaterial,part,kinEnergy,
547 cutEnergy,maxEnergy);
548}
549
550//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
551
552inline
554 G4double ekin,
555 const G4Material* material,
556 G4double emin,
558{
559 G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
560 return (cross > 0.0) ? 1./cross : DBL_MAX;
561}
562
563//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
564
565inline G4double
567 const G4Element* elm,
568 G4double kinEnergy,
569 G4double cutEnergy,
570 G4double maxEnergy)
571{
573 return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
574 cutEnergy,maxEnergy);
575}
576
577//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
578
579inline const G4Element*
581 const G4ParticleDefinition* part,
582 G4double kinEnergy,
583 G4double cutEnergy,
584 G4double maxEnergy)
585{
586 SetCurrentCouple(couple);
588 ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy) :
589 SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
590 fCurrentIsotope = nullptr;
591 return fCurrentElement;
592}
593
594//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
595
596inline const G4Element*
598 const G4ParticleDefinition* part,
599 G4double kinEnergy,
600 G4double logKinE,
601 G4double cutEnergy,
602 G4double maxEnergy)
603{
604 SetCurrentCouple(couple);
606 ? ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy,logKinE)
607 : SelectRandomAtom(pBaseMaterial,part,kinEnergy,cutEnergy,maxEnergy);
608 fCurrentIsotope = nullptr;
609 return fCurrentElement;
610}
611
612// ======== Get/Set inline methods used at initialisation ================
613
615{
616 return flucModel;
617}
618
619//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
620
622{
623 return anglModel;
624}
625
626//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
627
629{
630 if(p != anglModel) {
631 delete anglModel;
632 anglModel = p;
633 }
634}
635
636//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
637
639{
640 return fTripletModel;
641}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
644
646{
647 if(p != fTripletModel) {
648 delete fTripletModel;
649 fTripletModel = p;
650 }
651}
652
653//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
654
656{
657 return highLimit;
658}
659
660//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
661
663{
664 return lowLimit;
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668
670{
671 return eMaxActive;
672}
673
674//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
675
677{
678 return eMinActive;
679}
680
681//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682
684{
685 return polarAngleLimit;
686}
687
688//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
689
691{
692 return secondaryThreshold;
693}
694
695//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
696
698{
699 return theLPMflag;
700}
701
702//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
703
705{
706 return flagDeexcitation;
707}
708
709//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710
712{
713 return flagForceBuildTable;
714}
715
716//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
717
719{
720 return useAngularGenerator;
721}
722
723//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724
726{
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731
733{
734 lossFlucFlag = val;
735}
736
737//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738
740{
741 isMaster = val;
742}
743
744//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745
747{
748 return isMaster;
749}
750
751//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752
754{
755 useBaseMaterials = val;
756}
757
758//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759
761{
762 return useBaseMaterials;
763}
764
765//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766
768{
769 highLimit = val;
770}
771
772//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773
775{
776 lowLimit = val;
777}
778
779//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780
782{
783 eMaxActive = val;
784}
785
786//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787
789{
790 eMinActive = val;
791}
792
793//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794
795inline G4bool G4VEmModel::IsActive(G4double kinEnergy) const
796{
797 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
798}
799
800//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801
803{
804 if(!isLocked) { polarAngleLimit = val; }
805}
806
807//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808
810{
811 secondaryThreshold = val;
812}
813
814//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
815
817{
818 theLPMflag = val;
819}
820
821//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
822
824{
825 flagDeexcitation = val;
826}
827
828//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
829
831{
833}
834
835//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
836
837inline const G4String& G4VEmModel::GetName() const
838{
839 return name;
840}
841
842//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
843
844inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
845{
846 return elmSelectors;
847}
848
849//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
850
851inline void
852G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
853{
854 if(p != elmSelectors) {
855 elmSelectors = p;
856 nSelectors = (nullptr != elmSelectors) ? G4int(elmSelectors->size()) : 0;
857 localElmSelectors = false;
858 }
859}
860
861//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
862
864{
865 return fElementData;
866}
867
868//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
869
871{
872 return xSectionTable;
873}
874
875//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
876
878{
879 return isLocked;
880}
881
882//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
883
885{
886 isLocked = val;
887}
888
889//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
890
891#endif
static const G4double emax
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:229
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:365
G4double highLimit
Definition: G4VEmModel.hh:437
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:455
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:802
G4VEmModel(const G4VEmModel &)=delete
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:415
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:426
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:428
G4double polarAngleLimit
Definition: G4VEmModel.hh:441
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:739
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:614
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:852
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:683
std::vector< G4double > xsec
Definition: G4VEmModel.hh:466
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:621
G4bool theLPMflag
Definition: G4VEmModel.hh:454
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:209
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:830
G4double inveplus
Definition: G4VEmModel.hh:431
G4bool localTable
Definition: G4VEmModel.hh:459
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:844
G4bool isMaster
Definition: G4VEmModel.hh:457
G4double lowLimit
Definition: G4VEmModel.hh:436
const G4Element * fCurrentElement
Definition: G4VEmModel.hh:417
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4bool useAngularGenerator
Definition: G4VEmModel.hh:461
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:809
G4VEmModel * GetTripletModel()
Definition: G4VEmModel.hh:638
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
G4bool isLocked
Definition: G4VEmModel.hh:463
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:553
G4int nSelectors
Definition: G4VEmModel.hh:443
G4bool lossFlucFlag
Definition: G4VEmModel.hh:450
G4bool localElmSelectors
Definition: G4VEmModel.hh:460
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:505
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:725
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
void SetFluctuationFlag(G4bool val)
Definition: G4VEmModel.hh:732
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:816
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:223
G4int nsec
Definition: G4VEmModel.hh:444
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:297
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:406
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
G4double eMinActive
Definition: G4VEmModel.hh:438
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:414
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:528
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:497
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:425
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:580
const G4String name
Definition: G4VEmModel.hh:465
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:427
G4bool LPMFlag() const
Definition: G4VEmModel.hh:697
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:360
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:469
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:319
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:440
G4bool useBaseMaterials
Definition: G4VEmModel.hh:462
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:781
G4double pFactor
Definition: G4VEmModel.hh:432
void SetTripletModel(G4VEmModel *)
Definition: G4VEmModel.hh:645
std::vector< G4EmElementSelector * > * elmSelectors
Definition: G4VEmModel.hh:419
G4ElementData * fElementData
Definition: G4VEmModel.hh:424
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:669
void SetLocked(G4bool)
Definition: G4VEmModel.hh:884
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:823
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:539
const G4Isotope * fCurrentIsotope
Definition: G4VEmModel.hh:418
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:66
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:795
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:351
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:490
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:399
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:429
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:704
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:382
G4VEmFluctuationModel * flucModel
Definition: G4VEmModel.hh:413
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:718
const G4String & GetName() const
Definition: G4VEmModel.hh:837
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:597
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:261
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:870
G4bool flagForceBuildTable
Definition: G4VEmModel.hh:456
G4double eMaxActive
Definition: G4VEmModel.hh:439
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:753
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:512
size_t currentCoupleIndex
Definition: G4VEmModel.hh:448
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:138
G4LossTableManager * fEmManager
Definition: G4VEmModel.hh:420
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4VEmModel * fTripletModel
Definition: G4VEmModel.hh:415
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:84
G4bool UseBaseMaterials() const
Definition: G4VEmModel.hh:760
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:676
G4bool IsLocked() const
Definition: G4VEmModel.hh:877
const G4MaterialCutsCouple * fCurrentCouple
Definition: G4VEmModel.hh:416
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:270
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:204
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:424
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:711
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:520
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:690
G4VEmModel & operator=(const G4VEmModel &right)=delete
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:432
size_t basedCoupleIndex
Definition: G4VEmModel.hh:449
G4bool flagDeexcitation
Definition: G4VEmModel.hh:455
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:374
G4double secondaryThreshold
Definition: G4VEmModel.hh:440
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:863
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:108
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62