Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4hCoulombScatteringModel.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:     G4hCoulombScatteringModel
 33 //
 34 // Author:        Vladimir Ivanchenko 
 35 //
 36 // Creation date: 08.06.2012 from G4eCoulombScatteringModel
 37 //
 38 // Modifications:
 39 //
 40 // Class Description:
 41 //
 42 // Implementation of Coulomb Scattering of a charge particle 
 43 // on Atomic Nucleus for interval of scattering anles in Lab system 
 44 // thetaMin - ThetaMax.
 45 //   The model based on analysis of J.M.Fernandez-Varea et al. 
 46 // NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with 
 47 // screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
 48 // 
 49 
 50 // -------------------------------------------------------------------
 51 //
 52 
 53 #ifndef G4hCoulombScatteringModel_h
 54 #define G4hCoulombScatteringModel_h 1
 55 
 56 #include "G4VEmModel.hh"
 57 #include "globals.hh"
 58 #include "G4MaterialCutsCouple.hh"
 59 #include "G4WentzelVIRelXSection.hh"
 60 
 61 class G4ParticleChangeForGamma;
 62 class G4ParticleDefinition;
 63 class G4IonTable;
 64 class G4NistManager;
 65 
 66 class G4hCoulombScatteringModel : public G4VEmModel
 67 {
 68 
 69 public:
 70 
 71   explicit G4hCoulombScatteringModel(G4bool combined = true);
 72  
 73   ~G4hCoulombScatteringModel() override;
 74 
 75   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 76 
 77   void InitialiseLocal(const G4ParticleDefinition*, 
 78                        G4VEmModel* masterModel) override;
 79 
 80   G4double ComputeCrossSectionPerAtom(
 81                                 const G4ParticleDefinition*,
 82         G4double kinEnergy, 
 83         G4double Z, 
 84         G4double A, 
 85         G4double cut,
 86         G4double emax) override;
 87 
 88   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 89        const G4MaterialCutsCouple*,
 90        const G4DynamicParticle*,
 91        G4double tmin,
 92        G4double maxEnergy) override;
 93 
 94   G4double MinPrimaryEnergy(const G4Material*,
 95           const G4ParticleDefinition*,
 96           G4double) final;
 97 
 98   // defines low energy limit of the model
 99   inline void SetLowEnergyThreshold(G4double val);
100 
101   // user definition of low-energy threshold of recoil
102   inline void SetRecoilThreshold(G4double eth);
103 
104   // defines low energy limit on energy transfer to atomic electron
105   inline void SetFixedCut(G4double);
106 
107   // low energy limit on energy transfer to atomic electron
108   inline G4double GetFixedCut() const;
109 
110   // hide assignment operator
111   G4hCoulombScatteringModel & operator=
112   (const G4hCoulombScatteringModel &right) = delete;
113   G4hCoulombScatteringModel(const  G4hCoulombScatteringModel&) = delete;
114 
115 protected:
116 
117   inline void DefineMaterial(const G4MaterialCutsCouple*);
118 
119   inline void SetupParticle(const G4ParticleDefinition*);
120 
121 private:
122  
123   G4IonTable*               theIonTable;
124   G4ParticleChangeForGamma* fParticleChange;
125   G4WentzelVIRelXSection*   wokvi;
126   G4NistManager*            fNistManager;
127 
128   const G4ParticleDefinition*  particle;
129   const G4ParticleDefinition*  theProton;
130   const std::vector<G4double>* pCuts;
131 
132   const G4MaterialCutsCouple* currentCouple;
133   const G4Material*           currentMaterial;
134   G4int                       currentMaterialIndex;
135 
136   G4double                  cosThetaMin;
137   G4double                  cosThetaMax;
138   G4double                  recoilThreshold;
139   G4double                  elecRatio;
140   G4double                  mass;
141 
142   G4double                  fixedCut;
143 
144   G4bool                    isCombined;  
145 };
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
149 inline
150 void G4hCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup) 
151 { 
152   if(cup != currentCouple) {
153     currentCouple = cup;
154     currentMaterial = cup->GetMaterial();
155     currentMaterialIndex = currentCouple->GetIndex(); 
156   }
157 }
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160 
161 inline 
162 void G4hCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
163 {
164   // Initialise mass and charge
165   if(p != particle) {
166     particle = p;
167     mass = particle->GetPDGMass();
168     wokvi->SetupParticle(p);
169   }
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173 
174 inline void G4hCoulombScatteringModel::SetRecoilThreshold(G4double eth)
175 {
176   recoilThreshold = eth;
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180 
181 inline void G4hCoulombScatteringModel::SetFixedCut(G4double val)
182 {
183   fixedCut = val;
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187 
188 inline G4double G4hCoulombScatteringModel::GetFixedCut() const
189 {
190   return fixedCut;
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
195 #endif
196