Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
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
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63using namespace std;
64
67 theParticle(0),
68 isInitialised(false)
69{
70 G4Exception("G4alphaIonisation::G4alphaIonisation","em0007",JustWarning,
71 " The process is not ready for use - incorrect results are expected");
73 SetStepFunction(0.2, 0.01*mm);
74 // SetIntegral(true);
76 // SetVerboseLevel(1);
77 mass = 0.0;
78 ratio = 0.0;
79 eth = 8*MeV;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85{}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90{
91 return (!p.IsShortLived() &&
92 std::fabs(p.GetPDGCharge()/CLHEP::eplus - 2) < 0.01);
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
98 const G4Material*,
99 G4double cut)
100{
101 G4double x = 0.5*cut/electron_mass_c2;
102 G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
103 return mass*(gam - 1.0);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109 const G4ParticleDefinition* part,
110 const G4ParticleDefinition* bpart)
111{
112 if(!isInitialised) {
113
114 theParticle = part;
115 G4String pname = part->GetParticleName();
116
117 // define base particle
118 const G4ParticleDefinition* theBaseParticle = 0;
119 if(bpart == 0) {
120 if(pname != "alpha") { theBaseParticle = G4Alpha::Alpha(); }
121 } else { theBaseParticle = bpart; }
122
123 mass = part->GetPDGMass();
124 ratio = electron_mass_c2/mass;
125
126 SetBaseParticle(theBaseParticle);
128
129 if (!EmModel(1)) { SetEmModel(new G4BraggIonModel(), 1); }
131
132 // model limit defined for alpha
133 eth = (EmModel(1)->HighEnergyLimit())*mass/proton_mass_c2;
135
138
139 if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
140 EmModel(2)->SetLowEnergyLimit(eth);
142 AddEmModel(2, EmModel(2), FluctModel());
143
144 isInitialised = true;
145 }
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
151{}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fIonisation
@ JustWarning
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetStepFunction(G4double v1, G4double v2)
void SetBaseParticle(const G4ParticleDefinition *p)
void SetLinearLossLimit(G4double val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4VEmFluctuationModel * FluctModel()
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
virtual void PrintInfo()
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
G4alphaIonisation(const G4String &name="alphaIoni")
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)
virtual ~G4alphaIonisation()
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41