Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GammaConversionToMuons.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// ------------ G4GammaConversionToMuons physics process ------
28// by H.Burkhardt, S. Kelner and R. Kokoulin, April 2002
29// -----------------------------------------------------------------------------
30//
31// 05-08-04: suppression of .icc file (mma)
32// 13-08-04, public ComputeCrossSectionPerAtom() and ComputeMeanFreePath() (mma)
33//
34// class description
35//
36// gamma ---> mu+ mu-
37// inherit from G4VDiscreteProcess
38//
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
42#ifndef G4GammaConversionToMuons_h
43#define G4GammaConversionToMuons_h 1
44
45#include "G4ios.hh"
46#include "globals.hh"
47#include "Randomize.hh"
48#include "G4VDiscreteProcess.hh"
49#include "G4PhysicsTable.hh"
50#include "G4PhysicsLogVector.hh"
52#include "G4Element.hh"
53#include "G4Step.hh"
54#include <vector>
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
60
62{
63public: // with description
64
66 const G4String& processName ="GammaToMuPair",
68
70
72 // true for Gamma only.
73
74 void BuildPhysicsTable(const G4ParticleDefinition&) override;
75 // here dummy, the total cross section parametrization is used rather
76 // than tables, just calling PrintInfoDefinition
77
79 // Print few lines of informations about the process: validity range,
80 // origine ..etc..
81 // Invoked by BuildThePhysicsTable().
82
84 // Set the factor to artificially increase the crossSection (default 1)
85
86 inline G4double GetCrossSecFactor() const { return CrossSecFactor;}
87 // Get the factor to artificially increase the cross section
88
89 G4double GetMeanFreePath(const G4Track& aTrack,
90 G4double previousStepSize,
91 G4ForceCondition* condition) override;
92 // It returns the MeanFreePath of the process for the current track :
93 // (energy, material)
94 // The previousStepSize and G4ForceCondition* are not used.
95 // This function overloads a virtual function of the base class.
96 // It is invoked by the ProcessManager of the Particle.
97
99 const G4Element* anElement);
100 // It returns the total CrossSectionPerAtom of the process,
101 // for the current DynamicGamma (energy), in anElement.
102
104 const G4Step& aStep) override;
105 // It computes the final state of the process (at end of step),
106 // returned as a ParticleChange object.
107 // This function overloads a virtual function of the base class.
108 // It is invoked by the ProcessManager of the Particle.
109
111
113 const G4Material* aMaterial);
114
115 // hide assignment operator as private
117 operator=(const G4GammaConversionToMuons &right) = delete;
119
120private:
121
122 const G4Element* SelectRandomAtom(const G4DynamicParticle* aDynamicGamma,
123 const G4Material* aMaterial);
124
125 G4double Mmuon;
126 G4double Rc;
127 G4double LimitEnergy; // energy limit for accurate x-section
128 G4double LowestEnergyLimit; // low energy limit of the model
129 G4double HighestEnergyLimit; // high energy limit of the model
130 G4double Energy5DLimit = 0.0; // high energy limit for 5D final state sampling
131
132 G4double MeanFreePath = DBL_MAX;// actual MeanFreePath (current medium)
133 G4double CrossSecFactor = 1.0; // factor to artificially increase
134 // the cross section
135
136 G4LossTableManager* fManager;
137 G4BetheHeitler5DModel* f5Dmodel = nullptr;
138 const G4ParticleDefinition* theGamma;
139 const G4ParticleDefinition* theMuonPlus;
140 const G4ParticleDefinition* theMuonMinus;
141 std::vector<G4double> temp;
142};
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
146#endif
147
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double ComputeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
G4double ComputeCrossSectionPerAtom(G4double GammaEnergy, G4int Z)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
G4GammaConversionToMuons & operator=(const G4GammaConversionToMuons &right)=delete
G4bool IsApplicable(const G4ParticleDefinition &) override
G4double GetCrossSectionPerAtom(const G4DynamicParticle *aDynamicGamma, const G4Element *anElement)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4GammaConversionToMuons(const G4String &processName="GammaToMuPair", G4ProcessType type=fElectromagnetic)
G4GammaConversionToMuons(const G4GammaConversionToMuons &)=delete
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
#define DBL_MAX
Definition templates.hh:62