Geant4  10.00.p01
G4AnnihiToMuPair.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 // $Id: G4AnnihiToMuPair.hh 66241 2012-12-13 18:34:42Z gunter $
28 //
29 // ------------ G4AnnihiToMuPair physics process ------
30 // by H.Burkhardt, S. Kelner and R. Kokoulin, November 2002
31 // -----------------------------------------------------------------------------
32 
33 // class description
34 //
35 // (high energy) e+ (atomic) e- ---> mu+ mu-
36 // inherit from G4VDiscreteProcess
37 //
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......//
39 //
40 // 04.02.03 : cosmetic simplifications (mma)
41 // 27.01.03 : first implementation (hbu)
42 //
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44 
45 #ifndef G4AnnihiToMuPair_h
46 #define G4AnnihiToMuPair_h 1
47 
48 #include "G4VDiscreteProcess.hh"
49 #include "globals.hh"
50 
52 class G4Track;
53 class G4Step;
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 
58 {
59  public: // with description
60 
61  G4AnnihiToMuPair(const G4String& processName ="AnnihiToMuPair",
63 
65 
67  // true for positron only.
68 
70  // here dummy, just calling PrintInfoDefinition
71  // the total cross section is calculated analytically
72 
73  void PrintInfoDefinition();
74  // Print few lines of informations about the process: validity range,
75  // origine ..etc..
76  // Invoked by BuildPhysicsTable().
77 
79  // Set the factor to artificially increase the crossSection (default 1)
80 
82  // Get the factor to artificially increase the cross section
83 
84  G4double CrossSectionPerVolume(G4double PositronEnergy,
85  const G4Material*);
86  // Compute total cross section
87 
89  G4double AtomicZ);
90  // Compute total cross section
91 
92  G4double GetMeanFreePath(const G4Track& aTrack,
93  G4double previousStepSize,
95  // It returns the MeanFreePath of the process for the current track :
96  // (energy, material)
97  // The previousStepSize and G4ForceCondition* are not used.
98  // This function overloads a virtual function of the base class.
99  // It is invoked by the ProcessManager of the Particle.
100 
101  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
102  const G4Step& aStep);
103  // It computes the final state of the process (at end of step),
104  // returned as a ParticleChange object.
105  // This function overloads a virtual function of the base class.
106  // It is invoked by the ProcessManager of the Particle.
107 
108  private:
109 
110  // hide assignment operator as private
113 
114  G4double LowestEnergyLimit; // Energy threshold of e+
115  G4double HighestEnergyLimit; // Limit of validity of the model
116 
117  G4double CurrentSigma; // the last value of cross section per volume
118 
119  G4double CrossSecFactor; // factor to artificially increase
120  // the cross section, static to make sure
121  // to have single value
122 };
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
126 #endif
127 
static const G4double fac
G4double ComputeCrossSectionPerAtom(G4double PositronEnergy, G4double AtomicZ)
G4bool IsApplicable(const G4ParticleDefinition &)
G4AnnihiToMuPair & operator=(const G4AnnihiToMuPair &right)
G4AnnihiToMuPair(const G4String &processName="AnnihiToMuPair", G4ProcessType type=fElectromagnetic)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
G4double GetCrossSecFactor()
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *)
G4double CrossSectionPerVolume(G4double PositronEnergy, const G4Material *)
void SetCrossSecFactor(G4double fac)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4ProcessType