Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/cuts/include/G4VRangeToEnergyConverter.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 // G4VRangeToEnergyConverter
 27 //
 28 // Class description:
 29 //
 30 // Base class for Range to Energy Converters.
 31 // Cut in energy corresponding to given cut value in range
 32 // is calculated for a material by using Convert() method.
 33 
 34 // Author: H.Kurashige, 05 October 2002 - First implementation
 35 // --------------------------------------------------------------------
 36 #ifndef G4VRangeToEnergyConverter_hh
 37 #define G4VRangeToEnergyConverter_hh 1
 38 
 39 #include <vector>
 40 
 41 #include "globals.hh"
 42 #include "G4ParticleDefinition.hh"
 43 #include "G4Material.hh"
 44 
 45 class G4VRangeToEnergyConverter
 46 {
 47 public:
 48 
 49   explicit G4VRangeToEnergyConverter();
 50 
 51   virtual ~G4VRangeToEnergyConverter();
 52 
 53   // operators are not used
 54   G4VRangeToEnergyConverter(const G4VRangeToEnergyConverter& r) = delete;
 55   G4VRangeToEnergyConverter& operator=
 56   (const G4VRangeToEnergyConverter &r) = delete;
 57   G4bool operator==(const G4VRangeToEnergyConverter& r) const = delete;
 58   G4bool operator!=(const G4VRangeToEnergyConverter& r) const = delete;
 59 
 60   // Calculate energy cut from given range cut for the material
 61   virtual G4double Convert(const G4double rangeCut, const G4Material* material);
 62 
 63   // Set energy range for all particle type
 64   // if highedge > 10 GeV, highedge value is not changed
 65   static void SetEnergyRange(const G4double lowedge, const G4double highedge);
 66 
 67   // Get energy range for all particle type
 68   static G4double GetLowEdgeEnergy();
 69   static G4double GetHighEdgeEnergy();
 70 
 71   // Get/set max cut energy for all particle type
 72   // No check on the value
 73   static G4double GetMaxEnergyCut();
 74   static void SetMaxEnergyCut(const G4double value);
 75   
 76   // Return pointer to the particle type which this converter takes care of
 77   inline const G4ParticleDefinition* GetParticleType() const;
 78 
 79   inline void SetVerboseLevel(G4int value);
 80   inline G4int GetVerboseLevel() const;
 81       // control flag for output message
 82       //  0: Silent
 83       //  1: Warning message
 84       //  2: More
 85 
 86 protected:
 87 
 88   virtual G4double ComputeValue(const G4int Z, const G4double kinEnergy) = 0;
 89 
 90 private:
 91 
 92   static void FillEnergyVector(const G4double emin, const G4double emax);
 93 
 94   G4double ConvertForGamma(const G4double rangeCut, const G4Material* material);
 95 
 96   G4double ConvertForElectron(const G4double rangeCut, 
 97                               const G4Material* material);
 98 
 99   inline G4double LiniearInterpolation(const G4double e1, const G4double e2, 
100                                        const G4double r1, const G4double r2, 
101                                        const G4double r);
102 
103 protected:
104 
105   const G4ParticleDefinition* theParticle = nullptr;
106   G4int fPDG = 0;
107 
108 private:
109 
110   static G4double sEmin;
111   static G4double sEmax; 
112   static std::vector<G4double>* sEnergy;
113   static G4int sNbinPerDecade;
114   static G4int sNbin;
115 
116   G4int verboseLevel = 1;
117   G4bool isFirstInstance = false;
118 };
119 
120 // ------------------
121 // Inline methods
122 // ------------------
123 
124 inline 
125 void G4VRangeToEnergyConverter::SetVerboseLevel(G4int value)
126 {
127   verboseLevel = value;
128 }
129 
130 inline 
131 G4int G4VRangeToEnergyConverter::GetVerboseLevel() const
132 {
133   return verboseLevel;
134 }
135 
136 inline 
137 const G4ParticleDefinition* G4VRangeToEnergyConverter::GetParticleType() const
138 {
139   return theParticle;
140 }
141 
142 inline G4double G4VRangeToEnergyConverter::LiniearInterpolation(
143                 const G4double e1, const G4double e2, 
144                 const G4double r1, const G4double r2, const G4double r)
145 {
146   return (r1 == r2) ? e1 : e1 + (e2 - e1)*(r - r1)/(r2 - r1);
147 }
148 
149 #endif
150