Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VEmProcess.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: G4VEmProcess.hh 95657 2016-02-17 13:03:36Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4VEmProcess
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 01.10.2003
38 //
39 // Modifications:
40 // 30-06-04 make destructor virtual (V.Ivanchenko)
41 // 09-08-04 optimise integral option (V.Ivanchenko)
42 // 11-08-04 add protected methods to access cuts (V.Ivanchenko)
43 // 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
44 // 16-09-04 Add flag for LambdaTable and method RecalculateLambda (VI)
45 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
46 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47 // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
48 // 09-05-05 Fix problem in logic when path boundary between materials (VI)
49 // 11-01-06 add A to parameters of ComputeCrossSectionPerAtom (VI)
50 // 01-02-06 put default value A=0. to keep compatibility with v5.2 (mma)
51 // 13-05-06 Add method to access model by index (V.Ivanchenko)
52 // 12-09-06 add SetModel() (mma)
53 // 25-09-07 More accurate handling zero xsect in
54 // PostStepGetPhysicalInteractionLength (V.Ivanchenko)
55 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
56 // 15-07-08 Reorder class members for further multi-thread development (VI)
57 // 17-02-10 Added pointer currentParticle (VI)
58 //
59 // Class Description:
60 //
61 // It is the unified Discrete process
62 
63 // -------------------------------------------------------------------
64 //
65 
66 #ifndef G4VEmProcess_h
67 #define G4VEmProcess_h 1
68 
70 
71 #include "G4VDiscreteProcess.hh"
72 #include "globals.hh"
73 #include "G4Material.hh"
74 #include "G4MaterialCutsCouple.hh"
75 #include "G4Track.hh"
76 #include "G4EmModelManager.hh"
77 #include "G4UnitsTable.hh"
78 #include "G4ParticleDefinition.hh"
80 #include "G4EmParameters.hh"
81 
82 class G4Step;
83 class G4VEmModel;
84 class G4DataVector;
85 class G4VParticleChange;
86 class G4PhysicsTable;
87 class G4PhysicsVector;
88 class G4EmBiasingManager;
89 class G4LossTableManager;
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94 {
95 public:
96 
97  G4VEmProcess(const G4String& name,
99 
100  virtual ~G4VEmProcess();
101 
102  //------------------------------------------------------------------------
103  // Virtual methods to be implemented in concrete processes
104  //------------------------------------------------------------------------
105 
106  virtual G4bool IsApplicable(const G4ParticleDefinition& p) override = 0;
107 
108  virtual void PrintInfo() = 0;
109 
110  virtual void ProcessDescription(std::ostream& outFile) const; // = 0;
111 
112 protected:
113 
114  virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
115 
116  //------------------------------------------------------------------------
117  // Method with standard implementation; may be overwritten if needed
118  //------------------------------------------------------------------------
119 
121  const G4Material*);
122 
123  //------------------------------------------------------------------------
124  // Implementation of virtual methods common to all Discrete processes
125  //------------------------------------------------------------------------
126 
127 public:
128 
129  // Initialise for build of tables
130  virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
131 
132  // Build physics table during initialisation
133  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
134 
135  // Called before tracking of each new G4Track
136  virtual void StartTracking(G4Track*) override;
137 
138  // implementation of virtual method, specific for G4VEmProcess
140  const G4Track& track,
141  G4double previousStepSize,
142  G4ForceCondition* condition) override;
143 
144  // implementation of virtual method, specific for G4VEmProcess
145  virtual G4VParticleChange* PostStepDoIt(const G4Track&,
146  const G4Step&) override;
147 
148  // Store PhysicsTable in a file.
149  // Return false in case of failure at I/O
151  const G4String& directory,
152  G4bool ascii = false) override;
153 
154  // Retrieve Physics from a file.
155  // (return true if the Physics Table can be build by using file)
156  // (return false if the process has no functionality or in case of failure)
157  // File name should is constructed as processName+particleName and the
158  // should be placed under the directory specifed by the argument.
160  const G4String& directory,
161  G4bool ascii) override;
162 
163  //------------------------------------------------------------------------
164  // Specific methods for Discrete EM post step simulation
165  //------------------------------------------------------------------------
166 
167  // It returns the cross section per volume for energy/ material
169  const G4MaterialCutsCouple* couple);
170 
171  // It returns the cross section of the process per atom
173  G4double Z, G4double A=0.,
174  G4double cut=0.0);
175 
176  G4double MeanFreePath(const G4Track& track);
177 
178  // It returns cross section per volume
179  inline G4double GetLambda(G4double& kinEnergy,
180  const G4MaterialCutsCouple* couple);
181 
182  //------------------------------------------------------------------------
183  // Specific methods to build and access Physics Tables
184  //------------------------------------------------------------------------
185 
186  // Binning for lambda table
187  void SetLambdaBinning(G4int nbins);
188 
189  // Min kinetic energy for tables
190  void SetMinKinEnergy(G4double e);
191 
192  // Min kinetic energy for high energy table
194 
195  // Max kinetic energy for tables
196  void SetMaxKinEnergy(G4double e);
197 
198  // Cross section table pointers
199  inline G4PhysicsTable* LambdaTable() const;
200  inline G4PhysicsTable* LambdaTablePrim() const;
201 
202  //------------------------------------------------------------------------
203  // Define and access particle type
204  //------------------------------------------------------------------------
205 
206  inline const G4ParticleDefinition* Particle() const;
207  inline const G4ParticleDefinition* SecondaryParticle() const;
208 
209  //------------------------------------------------------------------------
210  // Specific methods to set, access, modify models and basic parameters
211  //------------------------------------------------------------------------
212 
213 protected:
214  // Select model in run time
215  inline G4VEmModel* SelectModel(G4double& kinEnergy, size_t index);
216 
217 public:
218  // Select model by energy and region index
219  inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
220  size_t& idxRegion) const;
221 
222  // Add model for region, smaller value of order defines which
223  // model will be selected for a given energy interval
224  void AddEmModel(G4int, G4VEmModel*, const G4Region* region = nullptr);
225 
226  // return the assigned model
227  G4VEmModel* EmModel(G4int index = 1) const;
228 
229  // Assign a model to a process
230  void SetEmModel(G4VEmModel*, G4int index = 1);
231 
232  // Define new energy range for the model identified by the name
233  void UpdateEmModel(const G4String&, G4double, G4double);
234 
235  // Access to models
236  G4int GetNumberOfModels() const;
237  G4int GetNumberOfRegionModels(size_t couple_index) const;
238  G4VEmModel* GetRegionModel(G4int idx, size_t couple_index) const;
239  G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false) const;
240 
241  // Access to active model
242  inline const G4VEmModel* GetCurrentModel() const;
243 
244  // Access to the current G4Element
245  const G4Element* GetCurrentElement() const;
246 
247  // Biasing parameters
248  void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
249  inline G4double CrossSectionBiasingFactor() const;
250 
251  // Activate forced interaction
252  void ActivateForcedInteraction(G4double length = 0.0,
253  const G4String& r = "",
254  G4bool flag = true);
255 
256  void ActivateSecondaryBiasing(const G4String& region, G4double factor,
257  G4double energyLimit);
258 
259  inline void SetIntegral(G4bool val);
260 
261  inline void SetBuildTableFlag(G4bool val);
262 
263  //------------------------------------------------------------------------
264  // Other generic methods
265  //------------------------------------------------------------------------
266 
267 protected:
268 
269  virtual G4double GetMeanFreePath(const G4Track& track,
270  G4double previousStepSize,
271  G4ForceCondition* condition) override;
272 
274 
275  inline G4int LambdaBinning() const;
276 
277  inline G4double MinKinEnergy() const;
278 
279  inline G4double MaxKinEnergy() const;
280 
281  // Single scattering parameters
282  inline G4double PolarAngleLimit() const;
283 
284  inline G4bool IsIntegral() const;
285 
286  inline G4double RecalculateLambda(G4double kinEnergy,
287  const G4MaterialCutsCouple* couple);
288 
290 
291  inline void SetParticle(const G4ParticleDefinition* p);
292 
293  inline void SetSecondaryParticle(const G4ParticleDefinition* p);
294 
295  inline size_t CurrentMaterialCutsCoupleIndex() const;
296 
297  inline G4double GetGammaEnergyCut();
298 
300 
301  inline void SetStartFromNullFlag(G4bool val);
302 
303  inline void SetSplineFlag(G4bool val);
304 
305  inline const G4Element* GetTargetElement() const;
306 
307  inline const G4Isotope* GetTargetIsotope() const;
308 
309 private:
310 
311  void Clear();
312 
313  void BuildLambdaTable();
314 
315  void PrintInfoProcess(const G4ParticleDefinition&);
316 
317  void FindLambdaMax();
318 
319  void PrintWarning(G4String tit, G4double val);
320 
321  inline void DefineMaterial(const G4MaterialCutsCouple* couple);
322 
323  inline void ComputeIntegralLambda(G4double kinEnergy);
324 
325  inline G4double GetLambdaFromTable(G4double kinEnergy);
326 
327  inline G4double GetLambdaFromTablePrim(G4double kinEnergy);
328 
329  inline G4double GetCurrentLambda(G4double kinEnergy);
330 
331  inline G4double ComputeCurrentLambda(G4double kinEnergy);
332 
333  // hide copy constructor and assignment operator
334  G4VEmProcess(G4VEmProcess &) = delete;
335  G4VEmProcess & operator=(const G4VEmProcess &right) = delete;
336 
337  // ======== Parameters of the class fixed at construction =========
338 
339  G4LossTableManager* lManager;
340  G4EmParameters* theParameters;
341  G4EmModelManager* modelManager;
342  G4EmBiasingManager* biasManager;
343  const G4ParticleDefinition* theGamma;
344  const G4ParticleDefinition* theElectron;
345  const G4ParticleDefinition* thePositron;
346  const G4ParticleDefinition* secondaryParticle;
347 
348  G4bool buildLambdaTable;
349 
350  // ======== Parameters of the class fixed at initialisation =======
351 
352  std::vector<G4VEmModel*> emModels;
353  G4int numberOfModels;
354 
355  // tables and vectors
356  G4PhysicsTable* theLambdaTable;
357  G4PhysicsTable* theLambdaTablePrim;
358  std::vector<G4double> theEnergyOfCrossSectionMax;
359  std::vector<G4double> theCrossSectionMax;
360 
361  size_t idxLambda;
362  size_t idxLambdaPrim;
363 
364  const std::vector<G4double>* theCuts;
365  const std::vector<G4double>* theCutsGamma;
366  const std::vector<G4double>* theCutsElectron;
367  const std::vector<G4double>* theCutsPositron;
368  const std::vector<G4double>* theDensityFactor;
369  const std::vector<G4int>* theDensityIdx;
370 
371  G4int nLambdaBins;
372 
373  G4double minKinEnergy;
374  G4double minKinEnergyPrim;
375  G4double maxKinEnergy;
376  G4double lambdaFactor;
377  G4double biasFactor;
378 
379  G4bool integral;
380  G4bool applyCuts;
381  G4bool startFromNull;
382  G4bool splineFlag;
383  G4bool actMinKinEnergy;
384  G4bool actMaxKinEnergy;
385  G4bool actBinning;
386  G4bool actSpline;
387 
388  // ======== Cashed values - may be state dependent ================
389 
390 protected:
391 
393 
394 private:
395 
396  std::vector<G4DynamicParticle*> secParticles;
397 
398  G4VEmModel* currentModel;
399 
400  const G4ParticleDefinition* particle;
401  const G4ParticleDefinition* currentParticle;
402 
403  // cache
404  const G4Material* baseMaterial;
405  const G4Material* currentMaterial;
406  const G4MaterialCutsCouple* currentCouple;
407  size_t currentCoupleIndex;
408  size_t basedCoupleIndex;
409 
410  G4double mfpKinEnergy;
411  G4double preStepKinEnergy;
412  G4double preStepLambda;
413  G4double fFactor;
414  G4bool biasFlag;
415  G4bool weightFlag;
416 
417  G4int mainSecondaries;
418  G4int secID;
419  G4int fluoID;
420  G4int augerID;
421  G4int biasID;
422 };
423 
424 // ======== Run time inline methods ================
425 
427 {
428  return currentCoupleIndex;
429 }
430 
431 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
432 
434 {
435  return (*theCutsGamma)[currentCoupleIndex];
436 }
437 
438 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
439 
441 {
442  return (*theCutsElectron)[currentCoupleIndex];
443 }
444 
445 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
446 
447 inline void G4VEmProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
448 {
449  if(couple != currentCouple) {
450  currentCouple = couple;
451  currentMaterial = couple->GetMaterial();
452  baseMaterial = currentMaterial->GetBaseMaterial();
453  currentCoupleIndex = couple->GetIndex();
454  basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
455  fFactor = biasFactor*(*theDensityFactor)[currentCoupleIndex];
456  if(!baseMaterial) { baseMaterial = currentMaterial; }
457  mfpKinEnergy = DBL_MAX;
458  idxLambda = idxLambdaPrim = 0;
459  }
460 }
461 
462 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
463 
464 inline
466 {
467  if(1 < numberOfModels) {
468  currentModel = modelManager->SelectModel(kinEnergy, index);
469  }
470  currentModel->SetCurrentCouple(currentCouple);
471  return currentModel;
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
475 
476 inline
478  size_t& idxRegion) const
479 {
480  return modelManager->SelectModel(kinEnergy, idxRegion);
481 }
482 
483 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
484 
485 inline G4double G4VEmProcess::GetLambdaFromTable(G4double e)
486 {
487  return ((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
488 }
489 
490 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
491 
492 inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e)
493 {
494  return ((*theLambdaTablePrim)[basedCoupleIndex])->Value(e, idxLambdaPrim)/e;
495 }
496 
497 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
498 
499 inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
500 {
501  return currentModel->CrossSectionPerVolume(
502  baseMaterial,currentParticle, e,(*theCuts)[currentCoupleIndex]);
503 }
504 
505 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
506 
507 inline G4double G4VEmProcess::GetCurrentLambda(G4double e)
508 {
509  G4double x;
510  if(e >= minKinEnergyPrim) { x = GetLambdaFromTablePrim(e); }
511  else if(theLambdaTable) { x = GetLambdaFromTable(e); }
512  else { x = ComputeCurrentLambda(e); }
513  return fFactor*x;
514 }
515 
516 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
517 
518 inline G4double
520  const G4MaterialCutsCouple* couple)
521 {
522  DefineMaterial(couple);
523  SelectModel(kinEnergy, currentCoupleIndex);
524  return GetCurrentLambda(kinEnergy);
525 }
526 
527 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
528 
529 inline G4double
531 {
532  DefineMaterial(couple);
533  SelectModel(e, currentCoupleIndex);
534  return fFactor*ComputeCurrentLambda(e);
535 }
536 
537 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
538 
539 inline void G4VEmProcess::ComputeIntegralLambda(G4double e)
540 {
541  mfpKinEnergy = theEnergyOfCrossSectionMax[currentCoupleIndex];
542  if (e <= mfpKinEnergy) {
543  preStepLambda = GetCurrentLambda(e);
544 
545  } else {
546  G4double e1 = e*lambdaFactor;
547  if(e1 > mfpKinEnergy) {
548  preStepLambda = GetCurrentLambda(e);
549  G4double preStepLambda1 = GetCurrentLambda(e1);
550  if(preStepLambda1 > preStepLambda) {
551  mfpKinEnergy = e1;
552  preStepLambda = preStepLambda1;
553  }
554  } else {
555  preStepLambda = fFactor*theCrossSectionMax[currentCoupleIndex];
556  }
557  }
558 }
559 
560 // ======== Get/Set inline methods used at initialisation ================
561 
562 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563 
565 {
566  return nLambdaBins;
567 }
568 
569 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
570 
572 {
573  return minKinEnergy;
574 }
575 
576 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
577 
579 {
580  return maxKinEnergy;
581 }
582 
583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
584 
586 {
587  return theParameters->MscThetaLimit();
588 }
589 
590 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
591 
593 {
594  return biasFactor;
595 }
596 
597 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
598 
600 {
601  return theLambdaTable;
602 }
603 
604 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
605 
607 {
608  return theLambdaTablePrim;
609 }
610 
611 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
612 
614 {
615  return particle;
616 }
617 
618 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
619 
621 {
622  return secondaryParticle;
623 }
624 
625 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
626 
628 {
629  integral = val;
630 }
631 
632 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
633 
635 {
636  return integral;
637 }
638 
639 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
640 
642 {
643  buildLambdaTable = val;
644 }
645 
646 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
647 
649 {
650  return &fParticleChange;
651 }
652 
653 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
654 
656 {
657  particle = p;
658  currentParticle = p;
659 }
660 
661 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662 
664 {
665  secondaryParticle = p;
666 }
667 
668 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
669 
671 {
672  startFromNull = val;
673 }
674 
675 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
676 
678 {
679  splineFlag = val;
680  actSpline = true;
681 }
682 
683 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
684 
686 {
687  return currentModel->GetCurrentElement();
688 }
689 
690 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
691 
693 {
694  return currentModel->GetCurrentIsotope();
695 }
696 
697 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
698 
700 {
701  return currentModel;
702 }
703 
704 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705 
706 #endif
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="", G4bool flag=true)
G4double condition(const G4ErrorSymMatrix &m)
G4PhysicsTable * LambdaTable() const
const XML_Char * name
Definition: expat.h:151
virtual void PrintInfo()=0
G4double GetLambda(G4double &kinEnergy, const G4MaterialCutsCouple *couple)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:257
virtual void StartTracking(G4Track *) override
G4PhysicsTable * LambdaTablePrim() const
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
G4int GetNumberOfRegionModels(size_t couple_index) const
G4int GetNumberOfModels() const
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:91
void SetBuildTableFlag(G4bool val)
G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4VEmModel * EmModel(G4int index=1) const
const char * p
Definition: xmltok.h:285
G4double MscThetaLimit() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4ParticleChangeForGamma fParticleChange
G4double GetElectronEnergyCut()
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
void SetSplineFlag(G4bool val)
const G4VEmModel * GetCurrentModel() const
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4double MaxKinEnergy() const
virtual ~G4VEmProcess()
const G4ParticleDefinition * SecondaryParticle() const
const G4ParticleDefinition * Particle() const
void SetStartFromNullFlag(G4bool val)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
int G4int
Definition: G4Types.hh:78
void UpdateEmModel(const G4String &, G4double, G4double)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4int LambdaBinning() const
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
void SetMinKinEnergyPrim(G4double e)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
G4double MeanFreePath(const G4Track &track)
void SetLambdaBinning(G4int nbins)
void SetEmModel(G4VEmModel *, G4int index=1)
size_t CurrentMaterialCutsCoupleIndex() const
double A(double temperature)
G4VEmModel * GetRegionModel(G4int idx, size_t couple_index) const
G4double PolarAngleLimit() const
bool G4bool
Definition: G4Types.hh:79
G4double GetGammaEnergyCut()
void SetParticle(const G4ParticleDefinition *p)
G4ParticleChangeForGamma * GetParticleChange()
const G4Element * GetTargetElement() const
Definition: G4Step.hh:76
void SetIntegral(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idxRegion) const
const G4Isotope * GetTargetIsotope() const
void SetMaxKinEnergy(G4double e)
G4double CrossSectionBiasingFactor() const
virtual void ProcessDescription(std::ostream &outFile) const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4bool IsIntegral() const
const G4Element * GetCurrentElement() const
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override=0
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4double MinKinEnergy() const
#define DBL_MAX
Definition: templates.hh:83
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4ProcessType
const G4Isotope * GetCurrentIsotope() const
Definition: G4VEmModel.hh:473
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466
G4double RecalculateLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)