Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/highenergy/include/G4Vee2hadrons.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4Vee2hadrons
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 //
 36 // Creation date: 12.08.2004
 37 //
 38 // Modifications: 14 July 2014 N. Chikuma revised interfaces  
 39 //
 40 
 41 //
 42 // Class Description: base class to compute partial cross sections
 43 //                    of e+e- annihilation into hadrons and 
 44 //                    sample of final state in the centre mass frame
 45 
 46 // -------------------------------------------------------------------
 47 //
 48 #ifndef G4Vee2hadrons_h
 49 #define G4Vee2hadrons_h 1
 50 
 51 #include <vector>
 52 #include <CLHEP/Units/SystemOfUnits.h>
 53 
 54 #include "globals.hh"
 55 #include "G4ThreeVector.hh"
 56 #include "G4eeCrossSections.hh"
 57 #include "G4PhysicsLinearVector.hh"
 58 
 59 class G4DynamicParticle;
 60 class G4PhysicsVector;
 61 
 62 class G4Vee2hadrons 
 63 {
 64 
 65 public:
 66 
 67   explicit G4Vee2hadrons(G4eeCrossSections* cr,
 68     G4double vlowEnergy,
 69     G4double vhighEnergy,
 70     G4double vdelta) : cross(cr)
 71   {
 72     lowEnergy  = vlowEnergy;
 73     highEnergy = vhighEnergy;
 74     delta      = vdelta;
 75   };
 76 
 77   virtual ~G4Vee2hadrons() {};
 78 
 79   virtual G4double PeakEnergy() const = 0;
 80 
 81   virtual G4double ComputeCrossSection(G4double) const = 0;
 82 
 83   G4PhysicsVector* PhysicsVector() const
 84   {
 85     G4int nbins = std::max(3, G4int((highEnergy - lowEnergy)/delta) );
 86     G4PhysicsVector* pp = new G4PhysicsLinearVector(lowEnergy,highEnergy,nbins);
 87     return pp;
 88   };
 89 
 90   virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 91          G4double, const G4ThreeVector&) = 0;
 92 
 93   G4double LowEnergy() const {return lowEnergy;};
 94 
 95   G4double HighEnergy() const {return highEnergy;};
 96   
 97   // hide assignment operator
 98   G4Vee2hadrons & operator=(const  G4Vee2hadrons &right) = delete;
 99   G4Vee2hadrons(const  G4Vee2hadrons&) = delete;
100 
101 private:
102 
103   // parameters of the table
104   G4double lowEnergy;
105   G4double highEnergy;
106   G4double delta;       
107 
108 protected:
109 
110    G4eeCrossSections* cross;  // class to compute cross section
111 
112 };
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115 
116 #endif
117