Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4BetheBlochModel.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:     G4BetheBlochModel
 33 //
 34 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
 35 //
 36 // Creation date: 03.01.2002
 37 //
 38 // Modifications:
 39 //
 40 // Modified by Michel Maire and Vladimir Ivanchenko
 41 //
 42 // Class Description:
 43 //
 44 // Implementation of Bethe-Bloch model of energy loss and
 45 // delta-electron production by heavy charged particles
 46 
 47 // -------------------------------------------------------------------
 48 //
 49 
 50 #ifndef G4BetheBlochModel_h
 51 #define G4BetheBlochModel_h 1
 52 
 53 #include "G4VEmModel.hh"
 54 
 55 class G4EmCorrections;
 56 class G4ParticleChangeForLoss;
 57 class G4ICRU90StoppingData;
 58 class G4NistManager;
 59 
 60 class G4BetheBlochModel : public G4VEmModel
 61 {
 62 
 63 public:
 64 
 65   explicit G4BetheBlochModel(const G4ParticleDefinition* p = nullptr,
 66            const G4String& nam = "BetheBloch");
 67 
 68   ~G4BetheBlochModel() override;
 69 
 70   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 71 
 72   G4double MinEnergyCut(const G4ParticleDefinition*,
 73       const G4MaterialCutsCouple* couple) override;
 74 
 75   virtual G4double ComputeCrossSectionPerElectron(
 76          const G4ParticleDefinition*,
 77          G4double kineticEnergy,
 78          G4double cutEnergy,
 79          G4double maxEnergy);
 80          
 81   G4double ComputeCrossSectionPerAtom(
 82          const G4ParticleDefinition*,
 83          G4double kineticEnergy,
 84          G4double Z, G4double A,
 85          G4double cutEnergy,
 86          G4double maxEnergy) override;
 87                  
 88   G4double CrossSectionPerVolume(const G4Material*,
 89          const G4ParticleDefinition*,
 90          G4double kineticEnergy,
 91          G4double cutEnergy,
 92          G4double maxEnergy) override;
 93          
 94   G4double ComputeDEDXPerVolume(const G4Material*,
 95         const G4ParticleDefinition*,
 96         G4double kineticEnergy,
 97         G4double cutEnergy) override;
 98 
 99   G4double GetChargeSquareRatio(const G4ParticleDefinition* p,
100         const G4Material* mat,
101         G4double kineticEnergy) override;
102 
103   G4double GetParticleCharge(const G4ParticleDefinition* p,
104            const G4Material* mat,
105            G4double kineticEnergy) override;
106 
107   void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
108           const G4DynamicParticle* dp,
109           const G4double& length,
110           G4double& eloss) override;
111 
112   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
113        const G4MaterialCutsCouple*,
114        const G4DynamicParticle*,
115        G4double tmin,
116        G4double maxEnergy) override;
117 
118   // hide assignment operator
119   G4BetheBlochModel & operator=(const G4BetheBlochModel &right) = delete;
120   G4BetheBlochModel(const G4BetheBlochModel&) = delete;
121 
122 protected:
123 
124   G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
125             G4double kinEnergy) override;
126 
127   inline G4double GetChargeSquareRatio() const;
128 
129   inline void SetChargeSquareRatio(G4double val);
130 
131 private:
132 
133   void SetupParameters(const G4ParticleDefinition* p);
134 
135   const G4ParticleDefinition* particle = nullptr;
136   const G4ParticleDefinition* theElectron;
137   G4EmCorrections*            corr;
138   G4ParticleChangeForLoss*    fParticleChange = nullptr;
139   G4NistManager*              nist;
140   G4ICRU90StoppingData*       fICRU90 = nullptr;
141   const G4Material*           currentMaterial = nullptr;
142   const G4Material*           baseMaterial = nullptr;
143 
144   G4double mass = 0.0;
145   G4double tlimit = DBL_MAX;
146   G4double spin = 0.0;
147   G4double magMoment2 = 0.0;
148   G4double chargeSquare = 1.0;
149   G4double ratio = 1.0;
150   G4double formfact = 0.0;
151   G4double twoln10;
152   G4double fAlphaTlimit;
153   G4double fProtonTlimit;
154 
155   G4int    iICRU90 = -1;
156   G4bool   isIon = false;
157   G4bool   isAlpha = false;
158 };
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161 
162 inline G4double G4BetheBlochModel::GetChargeSquareRatio() const
163 {
164   return chargeSquare;
165 }
166 
167 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168 
169 inline void G4BetheBlochModel::SetChargeSquareRatio(G4double val)
170 {
171   chargeSquare = val;
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175 
176 #endif
177