|
||||
File indexing completed on 2025-01-18 09:57:55
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 0028 // ------------ G4AnnihiToMuPair physics process ------ 0029 // by H.Burkhardt, S. Kelner and R. Kokoulin, November 2002 0030 // ----------------------------------------------------------------------------- 0031 0032 // class description 0033 // 0034 // (high energy) e+ (atomic) e- ---> mu+ mu- 0035 // inherit from G4VDiscreteProcess 0036 // 0037 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......// 0038 // 0039 // 04.02.03 : cosmetic simplifications (mma) 0040 // 27.01.03 : first implementation (hbu) 0041 // 0042 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0043 0044 #ifndef G4AnnihiToMuPair_h 0045 #define G4AnnihiToMuPair_h 1 0046 0047 #include "G4VDiscreteProcess.hh" 0048 #include "globals.hh" 0049 0050 class G4LossTableManager; 0051 class G4ParticleDefinition; 0052 class G4Track; 0053 class G4Step; 0054 0055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0056 0057 class G4AnnihiToMuPair : public G4VDiscreteProcess 0058 { 0059 public: // with description 0060 0061 explicit G4AnnihiToMuPair(const G4String& processName ="AnnihiToMuPair", 0062 G4ProcessType type = fElectromagnetic); 0063 0064 ~G4AnnihiToMuPair() override; 0065 0066 G4bool IsApplicable(const G4ParticleDefinition&) override; 0067 // true for positron only. 0068 0069 void BuildPhysicsTable(const G4ParticleDefinition&) override; 0070 // here dummy, just calling PrintInfoDefinition 0071 // the total cross section is calculated analytically 0072 0073 void PrintInfoDefinition(); 0074 // Print few lines of informations about the process: validity range, 0075 // origine ..etc.. 0076 // Invoked by BuildPhysicsTable(). 0077 0078 void SetCrossSecFactor(G4double fac); 0079 // Set the factor to artificially increase the crossSection (default 1) 0080 0081 G4double GetCrossSecFactor() {return fCrossSecFactor;}; 0082 // Get the factor to artificially increase the cross section 0083 0084 G4double CrossSectionPerVolume(G4double positronEnergy, const G4Material*); 0085 0086 G4double ComputeCrossSectionPerElectron(const G4double energy); 0087 0088 G4double ComputeCrossSectionPerAtom(const G4double energy, const G4double Z); 0089 0090 G4double GetMeanFreePath(const G4Track& aTrack, 0091 G4double previousStepSize, 0092 G4ForceCondition* ) override; 0093 // It returns the MeanFreePath of the process for the current track : 0094 // (energy, material) 0095 // The previousStepSize and G4ForceCondition* are not used. 0096 // This function overloads a virtual function of the base class. 0097 // It is invoked by the ProcessManager of the Particle. 0098 0099 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 0100 const G4Step& aStep) override; 0101 // It computes the final state of the process (at end of step), 0102 // returned as a ParticleChange object. 0103 // This function overloads a virtual function of the base class. 0104 // It is invoked by the ProcessManager of the Particle. 0105 0106 // hide assignment operator as private 0107 G4AnnihiToMuPair& operator=(const G4AnnihiToMuPair &right) = delete; 0108 G4AnnihiToMuPair(const G4AnnihiToMuPair& ) = delete; 0109 0110 private: 0111 0112 G4LossTableManager* fManager; 0113 const G4ParticleDefinition* part1; 0114 const G4ParticleDefinition* part2; 0115 G4double fMass; 0116 0117 G4double fLowEnergyLimit; // Energy threshold of e+ 0118 G4double fHighEnergyLimit; // Limit of validity of the model 0119 G4double fCurrentSigma; // the last value of cross section per volume 0120 G4double fCrossSecFactor; // factor to increase the cross section 0121 G4String fInfo = "e+e->mu+mu-"; 0122 }; 0123 0124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0125 0126 #endif 0127
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |