Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4BraggModel.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:     G4BraggModel
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 //
 36 // Creation date: 03.01.2002
 37 //
 38 // Modifications:
 39 // 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
 40 // 24-01-03 Make models region aware (V.Ivanchenko)
 41 // 13-02-03 Add name (V.Ivanchenko)
 42 // 12-11-03 Fix for GenericIons (V.Ivanchenko)
 43 // 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
 44 // 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
 45 // 25-04-06 Added stopping data from PSTAR (V.Ivanchenko)
 46 // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio, 
 47 //          CorrectionsAlongStep needed for ions(V.Ivanchenko)
 48 
 49 //
 50 // Class Description:
 51 //
 52 // Implementation of energy loss and delta-electron production
 53 // by heavy slow charged particles using ICRU'49 and NIST evaluated data 
 54 // for protons
 55 
 56 // -------------------------------------------------------------------
 57 //
 58 
 59 #ifndef G4BraggModel_h
 60 #define G4BraggModel_h 1
 61 
 62 #include "G4VEmModel.hh"
 63 #include "G4PSTARStopping.hh"
 64 
 65 class G4ParticleChangeForLoss;
 66 class G4EmCorrections;
 67 class G4ICRU90StoppingData;
 68 class G4PSTARStopping;
 69 
 70 class G4BraggModel : public G4VEmModel
 71 {
 72 
 73 public:
 74 
 75   explicit G4BraggModel(const G4ParticleDefinition* p = nullptr,
 76       const G4String& nam = "Bragg");
 77 
 78   ~G4BraggModel() override;
 79 
 80   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 81 
 82   G4double MinEnergyCut(const G4ParticleDefinition*,
 83       const G4MaterialCutsCouple* couple) override;
 84 
 85   G4double ComputeCrossSectionPerElectron(
 86          const G4ParticleDefinition*,
 87          G4double kineticEnergy,
 88          G4double cutEnergy,
 89          G4double maxEnergy);
 90          
 91   G4double ComputeCrossSectionPerAtom(
 92          const G4ParticleDefinition*,
 93          G4double kineticEnergy,
 94          G4double Z, G4double A,
 95          G4double cutEnergy,
 96          G4double maxEnergy) override;
 97                  
 98   G4double CrossSectionPerVolume(const G4Material*,
 99          const G4ParticleDefinition*,
100          G4double kineticEnergy,
101          G4double cutEnergy,
102          G4double maxEnergy) override;
103          
104   G4double ComputeDEDXPerVolume(const G4Material*,
105         const G4ParticleDefinition*,
106         G4double kineticEnergy,
107         G4double cutEnergy) override;
108 
109   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
110        const G4MaterialCutsCouple*,
111        const G4DynamicParticle*,
112        G4double tmin,
113        G4double maxEnergy) override;
114 
115   // Compute ion charge 
116   G4double GetChargeSquareRatio(const G4ParticleDefinition*,
117         const G4Material*,
118         G4double kineticEnergy) override;
119 
120   G4double GetParticleCharge(const G4ParticleDefinition* p,
121            const G4Material* mat,
122            G4double kineticEnergy) override;
123 
124   // hide assignment operator
125   G4BraggModel & operator=(const  G4BraggModel &right) = delete;
126   G4BraggModel(const  G4BraggModel&) = delete;
127 
128 protected:
129 
130   void SetParticle(const G4ParticleDefinition* p);
131 
132   G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
133             G4double kinEnergy) final;
134 
135   inline void SetChargeSquareRatio(G4double val);
136 
137 private:
138 
139   void HasMaterial(const G4Material* material);
140 
141   G4double StoppingPower(const G4Material* material,
142        G4double kineticEnergy);
143 
144   G4double ElectronicStoppingPower(G4double z,
145                                    G4double kineticEnergy) const;
146 
147   G4double DEDX(const G4Material* material, G4double kineticEnergy);
148 
149   G4bool MolecIsInZiegler1988(const G4Material* material);
150 
151   G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
152 
153 protected:
154 
155   const G4ParticleDefinition* particle = nullptr;
156   G4ParticleDefinition* theElectron = nullptr;
157   G4ParticleChangeForLoss* fParticleChange = nullptr;
158 
159   const G4Material* currentMaterial = nullptr;
160   const G4Material* baseMaterial = nullptr;
161 
162   G4EmCorrections* corr = nullptr;
163 
164   static G4ICRU90StoppingData* fICRU90;
165   static G4PSTARStopping* fPSTAR;
166 
167   G4double mass = 0.0;
168   G4double spin = 0.0;
169   G4double chargeSquare = 1.0;
170   G4double massRate = 1.0;
171   G4double ratio = 1.0;
172   G4double protonMassAMU = 1.007276;
173   G4double lowestKinEnergy;
174   G4double theZieglerFactor;
175   G4double expStopPower125;  // Experimental Stopping power at 125keV
176 
177   G4int iMolecula = -1;   // index in the molecula's table
178   G4int iPSTAR = -1;      // index in NIST PSTAR
179   G4int iICRU90 = -1;
180 
181 private:
182 
183   G4bool isIon = false;
184   G4bool isFirst = false;
185 };
186 
187 inline void G4BraggModel::SetChargeSquareRatio(G4double val)
188 {
189   chargeSquare = val;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
194 #endif
195