Geant4  10.02
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 93264 2015-10-14 09:30:04Z 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 "G4Isotope.hh"
88 #include "G4DataVector.hh"
89 #include "G4VEmFluctuationModel.hh"
91 #include "G4EmElementSelector.hh"
92 #include <CLHEP/Random/RandomEngine.h>
93 #include <vector>
94 
95 class G4ElementData;
96 class G4PhysicsTable;
97 class G4Region;
98 class G4VParticleChange;
101 class G4Track;
102 class G4LossTableManager;
103 
105 {
106 
107 public:
108 
109  G4VEmModel(const G4String& nam);
110 
111  virtual ~G4VEmModel();
112 
113  //------------------------------------------------------------------------
114  // Virtual methods to be implemented for any concrete model
115  //------------------------------------------------------------------------
116 
117  virtual void Initialise(const G4ParticleDefinition*, 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  // method to get partial cross section
160  G4int /*level*/,
161  const G4ParticleDefinition*,
162  G4double /*kineticEnergy*/)
163  { return 0;}
164 
165  // main method to compute cross section per atom
167  G4double kinEnergy,
168  G4double Z,
169  G4double A = 0., /* amu */
170  G4double cutEnergy = 0.0,
171  G4double maxEnergy = DBL_MAX);
172 
173  // main method to compute cross section per atomic shell
175  G4int Z, G4int shellIdx,
176  G4double kinEnergy,
177  G4double cutEnergy = 0.0,
178  G4double maxEnergy = DBL_MAX);
179 
180  // Compute effective ion charge square
181  virtual G4double ChargeSquareRatio(const G4Track&);
182 
183  // Compute effective ion charge square
185  const G4Material*,
186  G4double kineticEnergy);
187 
188  // Compute ion charge
190  const G4Material*,
191  G4double kineticEnergy);
192 
193  // Initialisation for a new track
194  virtual void StartTracking(G4Track*);
195 
196  // add correction to energy loss and compute non-ionizing energy loss
197  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
198  const G4DynamicParticle*,
199  G4double& eloss,
200  G4double& niel,
201  G4double length);
202 
203  // value which may be tabulated (by default cross section)
204  virtual G4double Value(const G4MaterialCutsCouple*,
205  const G4ParticleDefinition*,
206  G4double kineticEnergy);
207 
208  // threshold for zero value
209  virtual G4double MinPrimaryEnergy(const G4Material*,
210  const G4ParticleDefinition*,
211  G4double cut = 0.0);
212 
213  // model can define low-energy limit for the cut
215  const G4MaterialCutsCouple*);
216 
217  // initilisation at run time for a given material
218  virtual void SetupForMaterial(const G4ParticleDefinition*,
219  const G4Material*,
220  G4double kineticEnergy);
221 
222  // add a region for the model
223  virtual void DefineForRegion(const G4Region*);
224 
225  // for automatic documentation
226  virtual void ModelDescription(std::ostream& outFile) const; //=0;
227 
228 protected:
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 
240 public:
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
258  const G4ParticleDefinition*,
259  G4double kineticEnergy,
260  G4double cutEnergy = DBL_MAX);
261 
262  // cross section per volume
264  const G4ParticleDefinition*,
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  // select isotope in order to have precise mass of the nucleus
284  inline G4int SelectIsotopeNumber(const G4Element*);
285 
286  // atom can be selected effitiantly if element selectors are initialised
287  inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*,
288  const G4ParticleDefinition*,
289  G4double kineticEnergy,
290  G4double cutEnergy = 0.0,
291  G4double maxEnergy = DBL_MAX);
292 
293  // to select atom cross section per volume is recomputed for each element
294  const G4Element* SelectRandomAtom(const G4Material*,
295  const G4ParticleDefinition*,
296  G4double kineticEnergy,
297  G4double cutEnergy = 0.0,
298  G4double maxEnergy = DBL_MAX);
299 
300  // to select atom if cross section is proportional number of electrons
301  inline G4int SelectRandomAtomNumber(const G4Material*);
302 
303  //------------------------------------------------------------------------
304  // Get/Set methods
305  //------------------------------------------------------------------------
306 
308 
310 
311  inline G4ElementData* GetElementData();
312 
314 
316 
318 
320 
321  inline G4double HighEnergyLimit() const;
322 
323  inline G4double LowEnergyLimit() const;
324 
325  inline G4double HighEnergyActivationLimit() const;
326 
327  inline G4double LowEnergyActivationLimit() const;
328 
329  inline G4double PolarAngleLimit() const;
330 
331  inline G4double SecondaryThreshold() const;
332 
333  inline G4bool LPMFlag() const;
334 
335  inline G4bool DeexcitationFlag() const;
336 
337  inline G4bool ForceBuildTableFlag() const;
338 
339  inline G4bool UseAngularGeneratorFlag() const;
340 
341  inline void SetAngularGeneratorFlag(G4bool);
342 
343  inline void SetHighEnergyLimit(G4double);
344 
345  inline void SetLowEnergyLimit(G4double);
346 
348 
350 
351  inline G4bool IsActive(G4double kinEnergy);
352 
353  inline void SetPolarAngleLimit(G4double);
354 
355  inline void SetSecondaryThreshold(G4double);
356 
357  inline void SetLPMFlag(G4bool val);
358 
359  inline void SetDeexcitationFlag(G4bool val);
360 
361  inline void SetForceBuildTable(G4bool val);
362 
363  inline void SetMasterThread(G4bool val);
364 
365  inline G4bool IsMaster() const;
366 
367  inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
368 
369  inline const G4String& GetName() const;
370 
371  inline void SetCurrentCouple(const G4MaterialCutsCouple*);
372 
373  inline const G4Element* GetCurrentElement() const;
374 
375  inline const G4Isotope* GetCurrentIsotope() const;
376 
377  inline G4bool IsLocked() const;
378 
379  inline void SetLocked(G4bool);
380 
381 protected:
382 
383  inline const G4MaterialCutsCouple* CurrentCouple() const;
384 
385  inline void SetCurrentElement(const G4Element*);
386 
387 private:
388 
389  // hide assignment operator
391  G4VEmModel(const G4VEmModel&);
392 
393  // ======== Parameters of the class fixed at construction =========
394 
397  const G4String name;
398 
399  // ======== Parameters of the class fixed at initialisation =======
400 
411 
417  std::vector<G4EmElementSelector*>* elmSelectors;
418 
419 protected:
420 
424  const std::vector<G4double>* theDensityFactor;
425  const std::vector<G4int>* theDensityIdx;
426  size_t idxTable;
427  const static G4double inveplus;
428 
429  // ======== Cashed values - may be state dependent ================
430 
431 private:
432 
437 
439  std::vector<G4double> xsec;
440 
441 };
442 
443 // ======== Run time inline methods ================
444 
446 {
447  fCurrentCouple = p;
448 }
449 
450 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
451 
453 {
454  return fCurrentCouple;
455 }
456 
457 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458 
460 {
461  fCurrentElement = elm;
462  fCurrentIsotope = 0;
463 }
464 
465 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466 
468 {
469  return fCurrentElement;
470 }
471 
472 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473 
475 {
476  return fCurrentIsotope;
477 }
478 
479 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480 
481 inline
483 {
484  return MaxSecondaryEnergy(dynPart->GetParticleDefinition(),
485  dynPart->GetKineticEnergy());
486 }
487 
488 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
489 
491  const G4ParticleDefinition* part,
492  G4double kinEnergy,
493  G4double cutEnergy)
494 {
495  SetCurrentCouple(couple);
496  return ComputeDEDXPerVolume(couple->GetMaterial(),part,kinEnergy,cutEnergy);
497 }
498 
499 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
500 
502  const G4ParticleDefinition* part,
503  G4double kinEnergy,
504  G4double cutEnergy,
505  G4double maxEnergy)
506 {
507  SetCurrentCouple(couple);
508  return CrossSectionPerVolume(couple->GetMaterial(),part,kinEnergy,
509  cutEnergy,maxEnergy);
510 }
511 
512 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
513 
514 inline
516  G4double ekin,
517  const G4Material* material,
518  G4double emin,
519  G4double emax)
520 {
521  G4double mfp = DBL_MAX;
522  G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
523  if (cross > DBL_MIN) { mfp = 1./cross; }
524  return mfp;
525 }
526 
527 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528 
529 inline G4double
531  const G4Element* elm,
532  G4double kinEnergy,
533  G4double cutEnergy,
534  G4double maxEnergy)
535 {
536  SetCurrentElement(elm);
537  return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
538  cutEnergy,maxEnergy);
539 }
540 
541 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
542 
543 inline const G4Element*
545  const G4ParticleDefinition* part,
546  G4double kinEnergy,
547  G4double cutEnergy,
548  G4double maxEnergy)
549 {
550  fCurrentCouple = couple;
551  if(nSelectors > 0) {
552  fCurrentElement =
553  ((*elmSelectors)[couple->GetIndex()])->SelectRandomAtom(kinEnergy);
554  } else {
555  fCurrentElement = SelectRandomAtom(couple->GetMaterial(),part,kinEnergy,
556  cutEnergy,maxEnergy);
557  }
558  fCurrentIsotope = 0;
559  return fCurrentElement;
560 }
561 
562 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563 
565 {
566  // this algorith assumes that cross section is proportional to
567  // number electrons multiplied by number of atoms
568  size_t nn = mat->GetNumberOfElements();
569  const G4ElementVector* elmv = mat->GetElementVector();
570  G4int Z = G4lrint((*elmv)[0]->GetZ());
571  if(1 < nn) {
572  const G4double* at = mat->GetVecNbOfAtomsPerVolume();
574  for( size_t i=0; i<nn; ++i) {
575  Z = G4lrint((*elmv)[0]->GetZ());
576  tot -= Z*at[i];
577  if(tot <= 0.0) { break; }
578  }
579  }
580  return Z;
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584 
586 {
587  SetCurrentElement(elm);
588  G4int N = G4lrint(elm->GetN());
589  G4int ni = elm->GetNumberOfIsotopes();
590  fCurrentIsotope = 0;
591  if(ni > 0) {
592  G4int idx = 0;
593  if(ni > 1) {
596  for(; idx<ni; ++idx) {
597  x -= ab[idx];
598  if (x <= 0.0) { break; }
599  }
600  if(idx >= ni) { idx = ni - 1; }
601  }
602  fCurrentIsotope = elm->GetIsotope(idx);
603  N = fCurrentIsotope->GetN();
604  }
605  return N;
606 }
607 
608 // ======== Get/Set inline methods used at initialisation ================
609 
611 {
612  return flucModel;
613 }
614 
615 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
616 
618 {
619  return anglModel;
620 }
621 
622 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
623 
625 {
626  if(p != anglModel) {
627  delete anglModel;
628  anglModel = 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  isMaster = val;
714 }
715 
716 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
717 
719 {
720  return isMaster;
721 }
722 
723 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724 
726 {
727  highLimit = val;
728 }
729 
730 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731 
733 {
734  lowLimit = val;
735 }
736 
737 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738 
740 {
741  eMaxActive = val;
742 }
743 
744 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745 
747 {
748  eMinActive = val;
749 }
750 
751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752 
754 {
755  return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
756 }
757 
758 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759 
761 {
762  if(!isLocked) { polarAngleLimit = val; }
763 }
764 
765 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766 
768 {
769  secondaryThreshold = val;
770 }
771 
772 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773 
775 {
776  theLPMflag = val;
777 }
778 
779 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780 
782 {
783  flagDeexcitation = val;
784 }
785 
786 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787 
789 {
790  flagForceBuildTable = val;
791 }
792 
793 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794 
795 inline const G4String& G4VEmModel::GetName() const
796 {
797  return name;
798 }
799 
800 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801 
802 inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
803 {
804  return elmSelectors;
805 }
806 
807 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808 
809 inline void
810 G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
811 {
812  elmSelectors = p;
813  if(elmSelectors) { nSelectors = elmSelectors->size(); }
814  localElmSelectors = false;
815 }
816 
817 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
818 
820 {
821  return fElementData;
822 }
823 
824 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
825 
827 {
828  return xSectionTable;
829 }
830 
831 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
832 
834 {
835  return isLocked;
836 }
837 
838 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
839 
841 {
842  isLocked = val;
843 }
844 
845 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
846 
847 #endif
848 
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:655
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:648
G4bool localTable
Definition: G4VEmModel.hh:412
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:482
size_t idxTable
Definition: G4VEmModel.hh:426
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
G4bool flagForceBuildTable
Definition: G4VEmModel.hh:409
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:258
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:244
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:284
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
std::vector< G4Element * > G4ElementVector
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:690
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:362
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:148
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4bool flagDeexcitation
Definition: G4VEmModel.hh:408
G4double GetN() const
Definition: G4Element.hh:134
G4bool theLPMflag
Definition: G4VEmModel.hh:407
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:332
G4bool isMaster
Definition: G4VEmModel.hh:410
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:378
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:219
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:767
G4bool LPMFlag() const
Definition: G4VEmModel.hh:676
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:826
const G4String name
Definition: G4VEmModel.hh:397
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:610
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:387
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:819
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
const G4Element * fCurrentElement
Definition: G4VEmModel.hh:435
G4VEmFluctuationModel * flucModel
Definition: G4VEmModel.hh:395
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4int nSelectors
Definition: G4VEmModel.hh:416
G4bool isLocked
Definition: G4VEmModel.hh:415
G4LossTableManager * fManager
Definition: G4VEmModel.hh:433
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:704
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:425
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4double eMaxActive
Definition: G4VEmModel.hh:404
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:225
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
G4double secondaryThreshold
Definition: G4VEmModel.hh:406
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:410
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:432
const G4MaterialCutsCouple * fCurrentCouple
Definition: G4VEmModel.hh:434
G4double polarAngleLimit
Definition: G4VEmModel.hh:405
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
void SetLocked(G4bool)
Definition: G4VEmModel.hh:840
std::vector< G4EmElementSelector * > * elmSelectors
Definition: G4VEmModel.hh:417
bool G4bool
Definition: G4Types.hh:79
G4VEmModel & operator=(const G4VEmModel &right)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:501
static const G4double inveplus
Definition: G4VEmModel.hh:427
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:93
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:396
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:418
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:585
const G4ParticleDefinition * GetParticleDefinition() const
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
G4int nsec
Definition: G4VEmModel.hh:438
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:313
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:337
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:424
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:753
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const G4double emax
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:711
G4double eMinActive
Definition: G4VEmModel.hh:403
std::vector< G4double > xsec
Definition: G4VEmModel.hh:439
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:345
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:354
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:683
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:249
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:810
static const G4double ab
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:490
int G4lrint(double ad)
Definition: templates.hh:163
G4bool localElmSelectors
Definition: G4VEmModel.hh:413
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:774
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:445
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
const G4double x[NPOINTSGL]
G4double highLimit
Definition: G4VEmModel.hh:402
G4bool useAngularGenerator
Definition: G4VEmModel.hh:414
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:662
#define DBL_MIN
Definition: templates.hh:75
G4double lowLimit
Definition: G4VEmModel.hh:401
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:323
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:788
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:422
const G4String & GetName() const
Definition: G4VEmModel.hh:795
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:395
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:423
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
const G4Isotope * fCurrentIsotope
Definition: G4VEmModel.hh:436
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:515
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:459
#define DBL_MAX
Definition: templates.hh:83
virtual G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double)
Definition: G4VEmModel.hh:159
G4ElementData * fElementData
Definition: G4VEmModel.hh:421
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:760
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:369
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:564
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:474
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:467