|
||||
File indexing completed on 2025-01-18 09:59:20
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 // 0029 // GEANT4 Class header file 0030 // 0031 // 0032 // File name: G4Vee2hadrons 0033 // 0034 // Author: Vladimir Ivanchenko 0035 // 0036 // Creation date: 12.08.2004 0037 // 0038 // Modifications: 14 July 2014 N. Chikuma revised interfaces 0039 // 0040 0041 // 0042 // Class Description: base class to compute partial cross sections 0043 // of e+e- annihilation into hadrons and 0044 // sample of final state in the centre mass frame 0045 0046 // ------------------------------------------------------------------- 0047 // 0048 #ifndef G4Vee2hadrons_h 0049 #define G4Vee2hadrons_h 1 0050 0051 #include <vector> 0052 #include <CLHEP/Units/SystemOfUnits.h> 0053 0054 #include "globals.hh" 0055 #include "G4ThreeVector.hh" 0056 #include "G4eeCrossSections.hh" 0057 #include "G4PhysicsLinearVector.hh" 0058 0059 class G4DynamicParticle; 0060 class G4PhysicsVector; 0061 0062 class G4Vee2hadrons 0063 { 0064 0065 public: 0066 0067 explicit G4Vee2hadrons(G4eeCrossSections* cr, 0068 G4double vlowEnergy, 0069 G4double vhighEnergy, 0070 G4double vdelta) : cross(cr) 0071 { 0072 lowEnergy = vlowEnergy; 0073 highEnergy = vhighEnergy; 0074 delta = vdelta; 0075 }; 0076 0077 virtual ~G4Vee2hadrons() {}; 0078 0079 virtual G4double PeakEnergy() const = 0; 0080 0081 virtual G4double ComputeCrossSection(G4double) const = 0; 0082 0083 G4PhysicsVector* PhysicsVector() const 0084 { 0085 G4int nbins = std::max(3, G4int((highEnergy - lowEnergy)/delta) ); 0086 G4PhysicsVector* pp = new G4PhysicsLinearVector(lowEnergy,highEnergy,nbins); 0087 return pp; 0088 }; 0089 0090 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, 0091 G4double, const G4ThreeVector&) = 0; 0092 0093 G4double LowEnergy() const {return lowEnergy;}; 0094 0095 G4double HighEnergy() const {return highEnergy;}; 0096 0097 // hide assignment operator 0098 G4Vee2hadrons & operator=(const G4Vee2hadrons &right) = delete; 0099 G4Vee2hadrons(const G4Vee2hadrons&) = delete; 0100 0101 private: 0102 0103 // parameters of the table 0104 G4double lowEnergy; 0105 G4double highEnergy; 0106 G4double delta; 0107 0108 protected: 0109 0110 G4eeCrossSections* cross; // class to compute cross section 0111 0112 }; 0113 0114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 0115 0116 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |