Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4alphaIonisation.cc
Go to the documentation of this file.
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 // $Id: G4alphaIonisation.cc 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4alphaIonisation
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 28.10.2009 created from G4ionIonisation
38 //
39 // Modifications:
40 //
41 //
42 //
43 // -------------------------------------------------------------------
44 //
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "G4alphaIonisation.hh"
49 
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "G4Electron.hh"
53 #include "G4Alpha.hh"
54 #include "G4BraggIonModel.hh"
55 #include "G4BetheBlochModel.hh"
56 #include "G4UnitsTable.hh"
57 #include "G4LossTableManager.hh"
58 #include "G4IonFluctuations.hh"
60 #include "G4EmParameters.hh"
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
64 using namespace std;
65 
67  : G4VEnergyLossProcess(name),
68  theParticle(nullptr),
69  isInitialised(false)
70 {
71  G4Exception("G4alphaIonisation::G4alphaIonisation","em0007",JustWarning,
72  " The process is not ready for use - incorrect results are expected");
73  SetLinearLossLimit(0.02);
75  mass = 0.0;
76  ratio = 0.0;
77  eth = 8*MeV;
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83 {}
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88 {
89  return (!p.IsShortLived() &&
90  std::fabs(p.GetPDGCharge()/CLHEP::eplus - 2) < 0.01);
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96  const G4Material*,
97  G4double cut)
98 {
99  G4double x = 0.5*cut/electron_mass_c2;
100  G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
101  return mass*(gam - 1.0);
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
107  const G4ParticleDefinition* part,
108  const G4ParticleDefinition* bpart)
109 {
110  if(!isInitialised) {
111 
112  theParticle = part;
113  G4String pname = part->GetParticleName();
114 
115  // define base particle
116  const G4ParticleDefinition* theBaseParticle = nullptr;
117  if(bpart == 0) {
118  if(pname != "alpha") { theBaseParticle = G4Alpha::Alpha(); }
119  } else { theBaseParticle = bpart; }
120 
121  mass = part->GetPDGMass();
122  ratio = mass/electron_mass_c2/mass;
123 
124  SetBaseParticle(theBaseParticle);
126 
127  if (!EmModel(1)) { SetEmModel(new G4BraggIonModel(), 1); }
128 
130  G4double emin = param->MinKinEnergy();
131  EmModel(1)->SetLowEnergyLimit(emin);
132 
133  // model limit defined for alpha
134  eth = (EmModel(1)->HighEnergyLimit())*ratio;
135  EmModel(1)->SetHighEnergyLimit(eth);
136  AddEmModel(1, EmModel(1), new G4IonFluctuations());
137 
139 
140  if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
141  EmModel(2)->SetLowEnergyLimit(eth);
143  AddEmModel(2, EmModel(2), FluctModel());
144 
145  isInitialised = true;
146  }
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150 
152 {}
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
const XML_Char * name
Definition: expat.h:151
G4double MaxKinEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
virtual G4bool IsApplicable(const G4ParticleDefinition &p) final
const char * p
Definition: xmltok.h:285
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
void SetLinearLossLimit(G4double val)
tuple x
Definition: test.py:50
void SetFluctModel(G4VEmFluctuationModel *)
virtual void PrintInfo() override
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4VEmFluctuationModel * FluctModel()
bool G4bool
Definition: G4Types.hh:79
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static constexpr double eplus
float electron_mass_c2
Definition: hepunit.py:274
G4double MinKinEnergy() const
void SetSecondaryParticle(const G4ParticleDefinition *p)
string pname
Definition: eplot.py:33
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetPDGMass() const
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) final
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
virtual ~G4alphaIonisation()
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
G4double GetPDGCharge() const
G4alphaIonisation(const G4String &name="alphaIoni")
void SetBaseParticle(const G4ParticleDefinition *p)
G4VEmModel * EmModel(G4int index=1) const