Geant4  10.03.p03
 All Classes 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 105120 2017-07-13 13:24:43Z 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"
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  explicit 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 
164  // main method to compute cross section per atom
166  G4double kinEnergy,
167  G4double Z,
168  G4double A = 0., /* amu */
169  G4double cutEnergy = 0.0,
170  G4double maxEnergy = DBL_MAX);
171 
172  // main method to compute cross section per atomic shell
174  G4int Z, G4int shellIdx,
175  G4double kinEnergy,
176  G4double cutEnergy = 0.0,
177  G4double maxEnergy = DBL_MAX);
178 
179  // Compute effective ion charge square
180  virtual G4double ChargeSquareRatio(const G4Track&);
181 
182  // Compute effective ion charge square
184  const G4Material*,
185  G4double kineticEnergy);
186 
187  // Compute ion charge
189  const G4Material*,
190  G4double kineticEnergy);
191 
192  // Initialisation for a new track
193  virtual void StartTracking(G4Track*);
194 
195  // add correction to energy loss and compute non-ionizing energy loss
196  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
197  const G4DynamicParticle*,
198  G4double& eloss,
199  G4double& niel,
200  G4double length);
201 
202  // value which may be tabulated (by default cross section)
203  virtual G4double Value(const G4MaterialCutsCouple*,
204  const G4ParticleDefinition*,
205  G4double kineticEnergy);
206 
207  // threshold for zero value
208  virtual G4double MinPrimaryEnergy(const G4Material*,
209  const G4ParticleDefinition*,
210  G4double cut = 0.0);
211 
212  // model can define low-energy limit for the cut
214  const G4MaterialCutsCouple*);
215 
216  // initilisation at run time for a given material
217  virtual void SetupForMaterial(const G4ParticleDefinition*,
218  const G4Material*,
219  G4double kineticEnergy);
220 
221  // add a region for the model
222  virtual void DefineForRegion(const G4Region*);
223 
224  // for automatic documentation
225  virtual void ModelDescription(std::ostream& outFile) const;
226 
227 protected:
228 
229  // initialisation of the ParticleChange for the model
231 
232  // initialisation of the ParticleChange for the model
234 
235  // kinematically allowed max kinetic energy of a secondary
237  G4double kineticEnergy);
238 
239 public:
240 
241  //------------------------------------------------------------------------
242  // Generic methods common to all models
243  //------------------------------------------------------------------------
244 
245  // should be called at initialisation to build element selectors
247  const G4DataVector&);
248 
249  // should be called at initialisation to access element selectors
250  inline std::vector<G4EmElementSelector*>* GetElementSelectors();
251 
252  // should be called at initialisation to set element selectors
253  inline void SetElementSelectors(std::vector<G4EmElementSelector*>*);
254 
255  // dEdx per unit length
256  virtual inline G4double ComputeDEDX(const G4MaterialCutsCouple*,
257  const G4ParticleDefinition*,
258  G4double kineticEnergy,
259  G4double cutEnergy = DBL_MAX);
260 
261  // cross section per volume
263  const G4ParticleDefinition*,
264  G4double kineticEnergy,
265  G4double cutEnergy = 0.0,
266  G4double maxEnergy = DBL_MAX);
267 
268  // compute mean free path via cross section per volume
270  G4double kineticEnergy,
271  const G4Material*,
272  G4double cutEnergy = 0.0,
273  G4double maxEnergy = DBL_MAX);
274 
275  // generic cross section per element
277  const G4Element*,
278  G4double kinEnergy,
279  G4double cutEnergy = 0.0,
280  G4double maxEnergy = DBL_MAX);
281 
282  // select isotope in order to have precise mass of the nucleus
283  inline G4int SelectIsotopeNumber(const G4Element*);
284 
285  // atom can be selected effitiantly if element selectors are initialised
286  inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*,
287  const G4ParticleDefinition*,
288  G4double kineticEnergy,
289  G4double cutEnergy = 0.0,
290  G4double maxEnergy = DBL_MAX);
291 
292  // to select atom cross section per volume is recomputed for each element
293  const G4Element* SelectRandomAtom(const G4Material*,
294  const G4ParticleDefinition*,
295  G4double kineticEnergy,
296  G4double cutEnergy = 0.0,
297  G4double maxEnergy = DBL_MAX);
298 
299  // to select atom if cross section is proportional number of electrons
300  inline G4int SelectRandomAtomNumber(const G4Material*);
301 
302  //------------------------------------------------------------------------
303  // Get/Set methods
304  //------------------------------------------------------------------------
305 
307 
309 
310  inline G4ElementData* GetElementData();
311 
313 
315 
317 
319 
320  inline G4double HighEnergyLimit() const;
321 
322  inline G4double LowEnergyLimit() const;
323 
324  inline G4double HighEnergyActivationLimit() const;
325 
326  inline G4double LowEnergyActivationLimit() const;
327 
328  inline G4double PolarAngleLimit() const;
329 
330  inline G4double SecondaryThreshold() const;
331 
332  inline G4bool LPMFlag() const;
333 
334  inline G4bool DeexcitationFlag() const;
335 
336  inline G4bool ForceBuildTableFlag() const;
337 
338  inline G4bool UseAngularGeneratorFlag() const;
339 
340  inline void SetAngularGeneratorFlag(G4bool);
341 
342  inline void SetHighEnergyLimit(G4double);
343 
344  inline void SetLowEnergyLimit(G4double);
345 
347 
349 
350  inline G4bool IsActive(G4double kinEnergy);
351 
352  inline void SetPolarAngleLimit(G4double);
353 
354  inline void SetSecondaryThreshold(G4double);
355 
356  inline void SetLPMFlag(G4bool val);
357 
358  inline void SetDeexcitationFlag(G4bool val);
359 
360  inline void SetForceBuildTable(G4bool val);
361 
362  inline void SetFluctuationFlag(G4bool val);
363 
364  inline void SetMasterThread(G4bool val);
365 
366  inline G4bool IsMaster() const;
367 
368  inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
369 
370  inline const G4String& GetName() const;
371 
372  inline void SetCurrentCouple(const G4MaterialCutsCouple*);
373 
374  inline const G4Element* GetCurrentElement() const;
375 
376  inline const G4Isotope* GetCurrentIsotope() const;
377 
378  inline G4bool IsLocked() const;
379 
380  inline void SetLocked(G4bool);
381 
382 protected:
383 
384  inline const G4MaterialCutsCouple* CurrentCouple() const;
385 
386  inline void SetCurrentElement(const G4Element*);
387 
388 private:
389 
390  // hide assignment operator
391  G4VEmModel & operator=(const G4VEmModel &right) = delete;
392  G4VEmModel(const G4VEmModel&) = delete;
393 
394  // ======== Parameters of the class fixed at construction =========
395 
396  G4VEmFluctuationModel* flucModel;
397  G4VEmAngularDistribution* anglModel;
398  const G4String name;
399 
400  // ======== Parameters of the class fixed at initialisation =======
401 
402  G4double lowLimit;
403  G4double highLimit;
404  G4double eMinActive;
405  G4double eMaxActive;
406  G4double polarAngleLimit;
407  G4double secondaryThreshold;
408  G4bool theLPMflag;
409  G4bool flagDeexcitation;
410  G4bool flagForceBuildTable;
411  G4bool isMaster;
412 
413  G4bool localTable;
414  G4bool localElmSelectors;
415  G4bool useAngularGenerator;
416  G4bool isLocked;
417  G4int nSelectors;
418  std::vector<G4EmElementSelector*>* elmSelectors;
419  G4LossTableManager* fEmManager;
420 
421 protected:
422 
426  const std::vector<G4double>* theDensityFactor;
427  const std::vector<G4int>* theDensityIdx;
428  size_t idxTable;
430  const static G4double inveplus;
431 
432  // ======== Cached values - may be state dependent ================
433 
434 private:
435 
436  const G4MaterialCutsCouple* fCurrentCouple;
437  const G4Element* fCurrentElement;
438  const G4Isotope* fCurrentIsotope;
439 
440  G4int nsec;
441  std::vector<G4double> xsec;
442 
443 };
444 
445 // ======== Run time inline methods ================
446 
448 {
449  fCurrentCouple = p;
450 }
451 
452 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
453 
455 {
456  return fCurrentCouple;
457 }
458 
459 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
460 
462 {
463  fCurrentElement = elm;
464  fCurrentIsotope = nullptr;
465 }
466 
467 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
468 
470 {
471  return fCurrentElement;
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
475 
477 {
478  return fCurrentIsotope;
479 }
480 
481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
482 
483 inline
485 {
486  return MaxSecondaryEnergy(dynPart->GetParticleDefinition(),
487  dynPart->GetKineticEnergy());
488 }
489 
490 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
491 
493  const G4ParticleDefinition* part,
494  G4double kinEnergy,
495  G4double cutEnergy)
496 {
497  SetCurrentCouple(couple);
498  return ComputeDEDXPerVolume(couple->GetMaterial(),part,kinEnergy,cutEnergy);
499 }
500 
501 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
502 
504  const G4ParticleDefinition* part,
505  G4double kinEnergy,
506  G4double cutEnergy,
507  G4double maxEnergy)
508 {
509  SetCurrentCouple(couple);
510  return CrossSectionPerVolume(couple->GetMaterial(),part,kinEnergy,
511  cutEnergy,maxEnergy);
512 }
513 
514 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
515 
516 inline
518  G4double ekin,
519  const G4Material* material,
520  G4double emin,
521  G4double emax)
522 {
523  G4double cross = CrossSectionPerVolume(material,part,ekin,emin,emax);
524  return cross > 0.0 ? 1./cross : DBL_MAX;
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 = nullptr;
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 = (*elmv)[0]->GetZasInt();
571  if(1 < nn) {
572  const G4double* at = mat->GetVecNbOfAtomsPerVolume();
574  for( size_t i=0; i<nn; ++i) {
575  tot -= at[i];
576  if(tot <= 0.0) {
577  Z = (*elmv)[i]->GetZasInt();
578  break;
579  }
580  }
581  }
582  return Z;
583 }
584 
585 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
586 
588 {
589  SetCurrentElement(elm);
590  G4int N = G4lrint(elm->GetN());
591  G4int ni = elm->GetNumberOfIsotopes();
592  fCurrentIsotope = nullptr;
593  if(ni > 0) {
594  G4int idx = 0;
595  if(ni > 1) {
598  for(; idx<ni; ++idx) {
599  x -= ab[idx];
600  if (x <= 0.0) { break; }
601  }
602  if(idx >= ni) { idx = ni - 1; }
603  }
604  fCurrentIsotope = elm->GetIsotope(idx);
605  N = fCurrentIsotope->GetN();
606  }
607  return N;
608 }
609 
610 // ======== Get/Set inline methods used at initialisation ================
611 
613 {
614  return flucModel;
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618 
620 {
621  return anglModel;
622 }
623 
624 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
625 
627 {
628  if(p != anglModel) {
629  delete anglModel;
630  anglModel = p;
631  }
632 }
633 
634 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635 
637 {
638  return highLimit;
639 }
640 
641 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
642 
644 {
645  return lowLimit;
646 }
647 
648 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
649 
651 {
652  return eMaxActive;
653 }
654 
655 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656 
658 {
659  return eMinActive;
660 }
661 
662 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
663 
665 {
666  return polarAngleLimit;
667 }
668 
669 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
670 
672 {
673  return secondaryThreshold;
674 }
675 
676 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
677 
679 {
680  return theLPMflag;
681 }
682 
683 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
684 
686 {
687  return flagDeexcitation;
688 }
689 
690 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691 
693 {
694  return flagForceBuildTable;
695 }
696 
697 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698 
700 {
701  return useAngularGenerator;
702 }
703 
704 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705 
707 {
708  useAngularGenerator = val;
709 }
710 
712 {
713  lossFlucFlag = val;
714 }
715 
716 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
717 
719 {
720  isMaster = val;
721 }
722 
723 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
724 
726 {
727  return isMaster;
728 }
729 
730 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
731 
733 {
734  highLimit = val;
735 }
736 
737 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
738 
740 {
741  lowLimit = val;
742 }
743 
744 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
745 
747 {
748  eMaxActive = val;
749 }
750 
751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
752 
754 {
755  eMinActive = val;
756 }
757 
758 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
759 
761 {
762  return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
763 }
764 
765 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
766 
768 {
769  if(!isLocked) { polarAngleLimit = val; }
770 }
771 
772 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
773 
775 {
776  secondaryThreshold = val;
777 }
778 
779 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
780 
782 {
783  theLPMflag = val;
784 }
785 
786 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
787 
789 {
790  flagDeexcitation = val;
791 }
792 
793 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794 
796 {
797  flagForceBuildTable = val;
798 }
799 
800 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
801 
802 inline const G4String& G4VEmModel::GetName() const
803 {
804  return name;
805 }
806 
807 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
808 
809 inline std::vector<G4EmElementSelector*>* G4VEmModel::GetElementSelectors()
810 {
811  return elmSelectors;
812 }
813 
814 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
815 
816 inline void
817 G4VEmModel::SetElementSelectors(std::vector<G4EmElementSelector*>* p)
818 {
819  elmSelectors = p;
820  if(elmSelectors) { nSelectors = elmSelectors->size(); }
821  localElmSelectors = false;
822 }
823 
824 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
825 
827 {
828  return fElementData;
829 }
830 
831 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
832 
834 {
835  return xSectionTable;
836 }
837 
838 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
839 
841 {
842  return isLocked;
843 }
844 
845 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
846 
848 {
849  isLocked = val;
850 }
851 
852 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
853 
854 #endif
855 
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:657
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:650
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:484
size_t idxTable
Definition: G4VEmModel.hh:428
const XML_Char * name
Definition: expat.h:151
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:257
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:243
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:292
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:671
std::vector< G4Element * > G4ElementVector
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:692
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:370
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:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
G4double GetN() const
Definition: G4Element.hh:135
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:418
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:340
const char * p
Definition: xmltok.h:285
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:386
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:218
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:411
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:619
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:774
G4bool LPMFlag() const
Definition: G4VEmModel.hh:678
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:833
tuple x
Definition: test.py:50
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:612
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:395
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:826
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
string material
Definition: eplot.py:19
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:706
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:427
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:454
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:224
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:440
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:699
G4bool lossFlucFlag
Definition: G4VEmModel.hh:429
void SetLocked(G4bool)
Definition: G4VEmModel.hh:847
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:503
static const G4double inveplus
Definition: G4VEmModel.hh:430
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:92
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:426
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:587
const G4ParticleDefinition * GetParticleDefinition() const
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:283
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:321
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:345
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:426
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:760
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
static const G4double emax
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:718
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:353
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:362
void SetFluctuationFlag(G4bool val)
Definition: G4VEmModel.hh:711
G4bool DeexcitationFlag() const
Definition: G4VEmModel.hh:685
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:248
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817
static const G4double ab
virtual G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:492
int G4lrint(double ad)
Definition: templates.hh:163
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:781
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:447
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:626
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:664
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:331
void SetForceBuildTable(G4bool val)
Definition: G4VEmModel.hh:795
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:424
**D E S C R I P T I O N
Definition: HEPEvtcom.cc:77
const G4String & GetName() const
Definition: G4VEmModel.hh:802
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:403
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4bool IsLocked() const
Definition: G4VEmModel.hh:840
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:425
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788
G4double ComputeMeanFreePath(const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:517
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:461
#define DBL_MAX
Definition: templates.hh:83
G4ElementData * fElementData
Definition: G4VEmModel.hh:423
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:767
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377
const G4Material * GetMaterial() const
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.hh:564
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:476
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:469