Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4VEmModel.hh 76333 2013-11-08 14:31:50Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4VEmModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 23-12-02 V.Ivanchenko change interface before move to cut per region
42 // 24-01-03 Cut per region (V.Ivanchenko)
43 // 13-02-03 Add name (V.Ivanchenko)
44 // 25-02-03 Add sample theta and displacement (V.Ivanchenko)
45 // 23-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection
46 // calculation (V.Ivanchenko)
47 // 01-03-04 L.Urban signature changed in SampleCosineTheta
48 // 23-04-04 L.urban signature of SampleCosineTheta changed back
49 // 17-11-04 Add method CrossSectionPerAtom (V.Ivanchenko)
50 // 14-03-05 Reduce number of pure virtual methods and make inline part
51 // separate (V.Ivanchenko)
52 // 24-03-05 Remove IsInCharge and add G4VParticleChange in the constructor (VI)
53 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
54 // 15-04-05 optimize internal interface for msc (V.Ivanchenko)
55 // 08-05-05 A -> N (V.Ivanchenko)
56 // 25-07-05 Move constructor and destructor to the body (V.Ivanchenko)
57 // 02-02-06 ComputeCrossSectionPerAtom: default value A=0. (mma)
58 // 06-02-06 add method ComputeMeanFreePath() (mma)
59 // 07-03-06 Optimize msc methods (V.Ivanchenko)
60 // 29-06-06 Add member currentElement and Get/Set methods (V.Ivanchenko)
61 // 29-10-07 Added SampleScattering (V.Ivanchenko)
62 // 15-07-08 Reorder class members and improve comments (VI)
63 // 21-07-08 Added vector of G4ElementSelector and methods to use it (VI)
64 // 12-09-08 Added methods GetParticleCharge, GetChargeSquareRatio,
65 // CorrectionsAlongStep, ActivateNuclearStopping (VI)
66 // 16-02-09 Moved implementations of virtual methods to source (VI)
67 // 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
68 // 13-10-10 Added G4VEmAngularDistribution (VI)
69 //
70 // Class Description:
71 //
72 // Abstract interface to energy loss models
73 
74 // -------------------------------------------------------------------
75 //
76 
77 #ifndef G4VEmModel_h
78 #define G4VEmModel_h 1
79 
80 #include "globals.hh"
81 #include "G4DynamicParticle.hh"
82 #include "G4ParticleDefinition.hh"
83 #include "G4MaterialCutsCouple.hh"
84 #include "G4Material.hh"
85 #include "G4Element.hh"
86 #include "G4ElementVector.hh"
87 #include "G4DataVector.hh"
88 #include "G4VEmFluctuationModel.hh"
90 #include "G4EmElementSelector.hh"
91 #include "Randomize.hh"
92 #include <vector>
93 
94 class G4ElementData;
95 class G4PhysicsTable;
96 class G4Region;
97 class G4VParticleChange;
100 class G4Track;
101 class G4LossTableManager;
102 
104 {
105 
106 public:
107 
108  G4VEmModel(const G4String& nam);
109 
110  virtual ~G4VEmModel();
111 
112  //------------------------------------------------------------------------
113  // Virtual methods to be implemented for any concrete model
114  //------------------------------------------------------------------------
115 
116  virtual void Initialise(const G4ParticleDefinition*,
117  const G4DataVector&) = 0;
118 
119  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
120  const G4MaterialCutsCouple*,
121  const G4DynamicParticle*,
122  G4double tmin = 0.0,
123  G4double tmax = DBL_MAX) = 0;
124 
125  //------------------------------------------------------------------------
126  // Methods for initialisation of MT; may be overwritten if needed
127  //------------------------------------------------------------------------
128 
129  // initilisation in local thread
130  virtual void InitialiseLocal(const G4ParticleDefinition*,
131  G4VEmModel* masterModel);
132 
133  // initilisation of a new material at run time
134  virtual void InitialiseForMaterial(const G4ParticleDefinition*,
135  const G4Material*);
136 
137  // initilisation of a new element at run time
138  virtual void InitialiseForElement(const G4ParticleDefinition*,
139  G4int Z);
140 
141  //------------------------------------------------------------------------
142  // Methods with standard implementation; may be overwritten if needed
143  //------------------------------------------------------------------------
144 
145  // main method to compute dEdx
147  const G4ParticleDefinition*,
148  G4double kineticEnergy,
149  G4double cutEnergy = DBL_MAX);
150 
151  // main method to compute cross section per Volume
153  const G4ParticleDefinition*,
154  G4double kineticEnergy,
155  G4double cutEnergy = 0.0,
156  G4double maxEnergy = DBL_MAX);
157 
158  // main method to compute cross section per atom
160  G4double kinEnergy,
161  G4double Z,
162  G4double A = 0., /* amu */
163  G4double cutEnergy = 0.0,
164  G4double maxEnergy = DBL_MAX);
165 
166  // Compute effective ion charge square
167  virtual G4double ChargeSquareRatio(const G4Track&);
168 
169  // Compute effective ion charge square
171  const G4Material*,
172  G4double kineticEnergy);
173 
174  // Compute ion charge
176  const G4Material*,
177  G4double kineticEnergy);
178 
179  // Initialisation for a new track
180  virtual void StartTracking(G4Track*);
181 
182  // add correction to energy loss and compute non-ionizing energy loss
183  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
184  const G4DynamicParticle*,
185  G4double& eloss,
186  G4double& niel,
187  G4double length);
188 
189  // value which may be tabulated (by default cross section)
190  virtual G4double Value(const G4MaterialCutsCouple*,
191  const G4ParticleDefinition*,
192  G4double kineticEnergy);
193 
194  // threshold for zero value
195  virtual G4double MinPrimaryEnergy(const G4Material*,
196  const G4ParticleDefinition*,
197  G4double cut = 0.0);
198 
199  // model can define low-energy limit for the cut
201  const G4MaterialCutsCouple*);
202 
203  // initilisation at run time for a given material
204  virtual void SetupForMaterial(const G4ParticleDefinition*,
205  const G4Material*,
206  G4double kineticEnergy);
207 
208  // add a region for the model
209  virtual void DefineForRegion(const G4Region*);
210 
211 protected:
212 
213  // initialisation of the ParticleChange for the model
215 
216  // initialisation of the ParticleChange for the model
218 
219  // kinematically allowed max kinetic energy of a secondary
221  G4double kineticEnergy);
222 
223 public:
224 
225  //------------------------------------------------------------------------
226  // Generic methods common to all models
227  //------------------------------------------------------------------------
228 
229  // should be called at initialisation to build element selectors
231  const G4DataVector&);
232 
233  // should be called at initialisation to access element selectors
234  inline std::vector<G4EmElementSelector*>* GetElementSelectors();
235 
236  // should be called at initialisation to set element selectors
237  inline void SetElementSelectors(std::vector<G4EmElementSelector*>*);
238 
239  // dEdx per unit length
241  const G4ParticleDefinition*,
242  G4double kineticEnergy,
243  G4double cutEnergy = DBL_MAX);
244 
245  // cross section per volume
247  const G4ParticleDefinition*,
248  G4double kineticEnergy,
249  G4double cutEnergy = 0.0,
250  G4double maxEnergy = DBL_MAX);
251 
252  // compute mean free path via cross section per volume
254  G4double kineticEnergy,
255  const G4Material*,
256  G4double cutEnergy = 0.0,
257  G4double maxEnergy = DBL_MAX);
258 
259  // generic cross section per element
261  const G4Element*,
262  G4double kinEnergy,
263  G4double cutEnergy = 0.0,
264  G4double maxEnergy = DBL_MAX);
265 
266  // select isotope in order to have precise mass of the nucleus
267  inline G4int SelectIsotopeNumber(const G4Element*);
268 
269  // atom can be selected effitiantly if element selectors are initialised
270  inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*,
271  const G4ParticleDefinition*,
272  G4double kineticEnergy,
273  G4double cutEnergy = 0.0,
274  G4double maxEnergy = DBL_MAX);
275 
276  // to select atom cross section per volume is recomputed for each element
277  const G4Element* SelectRandomAtom(const G4Material*,
278  const G4ParticleDefinition*,
279  G4double kineticEnergy,
280  G4double cutEnergy = 0.0,
281  G4double maxEnergy = DBL_MAX);
282 
283  // to select atom if cross section is proportional number of electrons
284  inline G4int SelectRandomAtomNumber(const G4Material*);
285 
286  //------------------------------------------------------------------------
287  // Get/Set methods
288  //------------------------------------------------------------------------
289 
291 
293 
294  inline G4ElementData* GetElementData();
295 
297 
299 
301 
303 
304  inline G4double HighEnergyLimit() const;
305 
306  inline G4double LowEnergyLimit() const;
307 
308  inline G4double HighEnergyActivationLimit() const;
309 
310  inline G4double LowEnergyActivationLimit() const;
311 
312  inline G4double PolarAngleLimit() const;
313 
314  inline G4double SecondaryThreshold() const;
315 
316  inline G4bool LPMFlag() const;
317 
318  inline G4bool DeexcitationFlag() const;
319 
320  inline G4bool ForceBuildTableFlag() const;
321 
322  inline G4bool UseAngularGeneratorFlag() const;
323 
324  inline void SetAngularGeneratorFlag(G4bool);
325 
326  inline void SetHighEnergyLimit(G4double);
327 
328  inline void SetLowEnergyLimit(G4double);
329 
331 
333 
334  inline G4bool IsActive(G4double kinEnergy);
335 
336  inline void SetPolarAngleLimit(G4double);
337 
338  inline void SetSecondaryThreshold(G4double);
339 
340  inline void SetLPMFlag(G4bool val);
341 
342  inline void SetDeexcitationFlag(G4bool val);
343 
344  inline void SetForceBuildTable(G4bool val);
345 
346  inline void SetMasterThread(G4bool val);
347 
348  inline G4bool IsMaster() const;
349 
350  inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
351 
352  inline const G4String& GetName() const;
353 
354  inline void SetCurrentCouple(const G4MaterialCutsCouple*);
355 
356  inline const G4Element* GetCurrentElement() const;
357 
358 protected:
359 
360  inline const G4MaterialCutsCouple* CurrentCouple() const;
361 
362  inline void SetCurrentElement(const G4Element*);
363 
364 private:
365 
366  // hide assignment operator
367  G4VEmModel & operator=(const G4VEmModel &right);
368  G4VEmModel(const G4VEmModel&);
369 
370  // ======== Parameters of the class fixed at construction =========
371 
372  G4VEmFluctuationModel* flucModel;
373  G4VEmAngularDistribution* anglModel;
374  const G4String name;
375 
376  // ======== Parameters of the class fixed at initialisation =======
377 
378  G4double lowLimit;
379  G4double highLimit;
380  G4double eMinActive;
381  G4double eMaxActive;
382  G4double polarAngleLimit;
383  G4double secondaryThreshold;
384  G4bool theLPMflag;
385  G4bool flagDeexcitation;
386  G4bool flagForceBuildTable;
387  G4bool isMaster;
388 
389  G4bool localTable;
390  G4bool localElmSelectors;
391  G4bool useAngularGenerator;
392  G4int nSelectors;
393  std::vector<G4EmElementSelector*>* elmSelectors;
394 
395 protected:
396 
400  const std::vector<G4double>* theDensityFactor;
401  const std::vector<G4int>* theDensityIdx;
402  size_t idxTable;
403 
404  // ======== Cashed values - may be state dependent ================
405 
406 private:
407 
408  G4LossTableManager* fManager;
409  const G4MaterialCutsCouple* fCurrentCouple;
410  const G4Element* fCurrentElement;
411 
412  G4int nsec;
413  std::vector<G4double> xsec;
414 
415 };
416 
417 // ======== Run time inline methods ================
418 
420 {
421  fCurrentCouple = p;
422 }
423 
424 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
425 
427 {
428  return fCurrentCouple;
429 }
430 
431 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
432 
434 {
435  fCurrentElement = elm;
436 }
437 
438 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
439 
441 {
442  return fCurrentElement;
443 }
444 
445 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
446 
447 inline
449 {
450  return MaxSecondaryEnergy(dynPart->GetParticleDefinition(),
451  dynPart->GetKineticEnergy());
452 }
453 
454 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455 
457  const G4ParticleDefinition* part,
458  G4double kinEnergy,
459  G4double cutEnergy)
460 {
461  SetCurrentCouple(couple);
462  return ComputeDEDXPerVolume(couple->GetMaterial(),part,kinEnergy,cutEnergy);
463 }
464 
465 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466 
468  const G4ParticleDefinition* part,
469  G4double kinEnergy,
470  G4double cutEnergy,
471  G4double maxEnergy)
472 {
473  SetCurrentCouple(couple);
474  return CrossSectionPerVolume(couple->GetMaterial(),part,kinEnergy,
475  cutEnergy,maxEnergy);
476 }
477 
478 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479 
480 inline
482  G4double ekin,
483  const G4Material* material,
484  G4double emin,
485  G4double emax)
486 {
487  G4double mfp = DBL_MAX;
488  G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
489  if (cross > DBL_MIN) { mfp = 1./cross; }
490  return mfp;
491 }
492 
493 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
494 
496  const G4ParticleDefinition* part,
497  const G4Element* elm,
498  G4double kinEnergy,
499  G4double cutEnergy,
500  G4double maxEnergy)
501 {
502  SetCurrentElement(elm);
503  return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
504  cutEnergy,maxEnergy);
505 }
506 
507 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
508 
509 inline const G4Element*
511  const G4ParticleDefinition* part,
512  G4double kinEnergy,
513  G4double cutEnergy,
514  G4double maxEnergy)
515 {
516  fCurrentCouple = couple;
517  if(nSelectors > 0) {
518  fCurrentElement =
519  ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy);
520  } else {
521  fCurrentElement = SelectRandomAtom(couple->GetMaterial(),part,kinEnergy,
522  cutEnergy,maxEnergy);
523  }
524  return fCurrentElement;
525 }
526 
527 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528 
530 {
531  size_t nn = mat->GetNumberOfElements();
532  const G4ElementVector* elmv = mat->GetElementVector();
533  G4int Z = G4int((*elmv)[0]->GetZ());
534  if(1 < nn) {
535  const G4double* at = mat->GetVecNbOfAtomsPerVolume();
537  for( size_t i=0; i<nn; ++i) {
538  Z = G4int((*elmv)[0]->GetZ());
539  tot -= Z*at[i];
540  if(tot <= 0.0) { break; }
541  }
542  }
543  return Z;
544 }
545 
546 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
547 
549 {
550  SetCurrentElement(elm);
551  G4int N = G4int(elm->GetN() + 0.5);
552  G4int ni = elm->GetNumberOfIsotopes();
553  if(ni > 0) {
554  G4int idx = 0;
555  if(ni > 1) {
558  for(; idx<ni; ++idx) {
559  x -= ab[idx];
560  if (x <= 0.0) { break; }
561  }
562  if(idx >= ni) { idx = ni - 1; }
563  }
564  N = elm->GetIsotope(idx)->GetN();
565  }
566  return N;
567 }
568 
569 // ======== Get/Set inline methods used at initialisation ================
570 
572 {
573  return flucModel;
574 }
575 
576 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
577 
579 {
580  return anglModel;
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584 
586 {
587  anglModel = p;
588 }
589 
590 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
591 
593 {
594  return highLimit;
595 }
596 
597 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
598 
600 {
601  return lowLimit;
602 }
603 
604 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
605 
607 {
608  return eMaxActive;
609 }
610 
611 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
612 
614 {
615  return eMinActive;
616 }
617 
618 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
619 
621 {
622  return polarAngleLimit;
623 }
624 
625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
626 
628 {
629  return secondaryThreshold;
630 }
631 
632 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
633 
635 {
636  return theLPMflag;
637 }
638 
639 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
640 
642 {
643  return flagDeexcitation;
644 }
645 
646 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
647 
649 {
650  return flagForceBuildTable;
651 }
652 
653 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
654 
656 {
657  return useAngularGenerator;
658 }
659 
660 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
661 
663 {
664  useAngularGenerator = val;
665 }
666 
667 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
668 
670 {
671  isMaster = val;
672 }
673 
674 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
675 
677 {
678  return isMaster;
679 }
680 
681 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
682 
684 {
685  highLimit = val;
686 }
687 
688 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
689 
691 {
692  lowLimit = val;
693 }
694 
695 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
696 
698 {
699  eMaxActive = val;
700 }
701 
702 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
703 
705 {
706  eMinActive = val;
707 }
708 
709 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
710 
712 {
713  return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
714 }
715 
716 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
717 
719 {
720  polarAngleLimit = val;
721 }
722 
723 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724 
726 {
727  secondaryThreshold = val;
728 }
729 
730 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731 
733 {
734  theLPMflag = val;
735 }
736 
737 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738 
740 {
741  flagDeexcitation = val;
742 }
743 
744 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745 
747 {
748  flagForceBuildTable = val;
749 }
750 
751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752 
753 inline const G4String& G4VEmModel::GetName() const
754 {
755  return name;
756 }
757 
758 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759 
760 inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
761 {
762  return elmSelectors;
763 }
764 
765 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766 
767 inline void
768 G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
769 {
770  elmSelectors = p;
771  if(elmSelectors) { nSelectors = elmSelectors->size(); }
772  localElmSelectors = false;
773 }
774 
775 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
776 
778 {
779  return fElementData;
780 }
781 
782 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
783 
785 {
786  return xSectionTable;
787 }
788 
789 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
790 
791 #endif
792 
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:613
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:606
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
size_t idxTable
Definition: G4VEmModel.hh:402
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:697
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:231
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:271
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:627
std::vector< G4Element * > G4ElementVector
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:648
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:339
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4double GetN() const
Definition: G4Element.hh:134
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:309
const char * p
Definition: xmltok.h:285
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:355
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:206
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:380
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:725
G4bool LPMFlag() const
Definition: G4VEmModel.hh:634
const XML_Char * name
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:784
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:571
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:364
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:777
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
string material
Definition: eplot.py:19
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:662
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:401
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:212
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:87
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:655
bool G4bool
Definition: G4Types.hh:79
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:467
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:88
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:395
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
const G4ParticleDefinition * GetParticleDefinition() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:314
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:400
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:669
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:322
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:331
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:641
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:236
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:456
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:732
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
#define DBL_MIN
Definition: templates.hh:75
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:746
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:398
**D E S C R I P T I O N
const G4String & GetName() const
Definition: G4VEmModel.hh:753
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:372
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:399
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:481
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:433
#define DBL_MAX
Definition: templates.hh:83
G4ElementData * fElementData
Definition: G4VEmModel.hh:397
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:529
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440