Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4hhIonisation.cc
이 파일의 문서화 페이지로 가기
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: G4hhIonisation.cc 106715 2017-10-20 09:39:06Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4hhIonisation
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 30.09.2005
38 //
39 // Modifications:
40 // 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
41 // 27-10-06 Add maxKinEnergy (V.Ivantchenko)
42 //
43 //
44 // -------------------------------------------------------------------
45 //
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
49 #include "G4hhIonisation.hh"
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "G4BraggNoDeltaModel.hh"
54 #include "G4ICRU73NoDeltaModel.hh"
56 #include "G4BohrFluctuations.hh"
57 #include "G4IonFluctuations.hh"
58 #include "G4UnitsTable.hh"
59 #include "G4Electron.hh"
60 #include "G4EmParameters.hh"
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
65  : G4VEnergyLossProcess(name),
66  theParticle(nullptr),
67  //theBaseParticle(nullptr),
68  isInitialised(false)
69 {
70  SetStepFunction(0.1, 0.1*mm);
71  SetVerboseLevel(1);
74  mass = 0.0;
75  ratio = 0.0;
76  flucModel = nullptr;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82 {}
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 100.0*MeV &&
89  !p.IsShortLived());
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95  const G4Material*,
96  G4double cut)
97 {
98  G4double x = 0.5*cut/electron_mass_c2;
100  G4double gam = x*y + std::sqrt((1. + x)*(1. + x*y*y));
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) { return; }
111 
112  theParticle = part;
113  if(bpart) {
114  G4cout << "G4hhIonisation::InitialiseEnergyLossProcess WARNING: no "
115  << "base particle should be defined for the process "
116  << GetProcessName() << G4endl;
117  }
118  SetBaseParticle(0);
123 
125  G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
126  G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
127 
128  SetMinKinEnergy(emin);
129  SetMaxKinEnergy(emax);
130  G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
131  SetDEDXBinning(bin);
132 
133  G4VEmModel* em = nullptr;
134  if(part->GetPDGCharge() > 0.0) { em = new G4BraggNoDeltaModel(); }
135  else { em = new G4ICRU73NoDeltaModel(); }
136  em->SetLowEnergyLimit(emin);
137  em->SetHighEnergyLimit(eth);
138  AddEmModel(1, em, flucModel);
139 
140  em = new G4BetheBlochNoDeltaModel();
141  em->SetLowEnergyLimit(eth);
142  em->SetHighEnergyLimit(emax);
143  SetEmModel(em);
144  AddEmModel(1, em, flucModel);
145 
146  if(verboseLevel>1) {
147  G4cout << "G4hhIonisation is initialised" << G4endl;
148  }
149  isInitialised = true;
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153 
155 {
156  G4cout << " Delta-ray will not be produced; "
157  << G4endl;
158 }
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161 
162 void G4hhIonisation::ProcessDescription(std::ostream& out) const
163 {
164  out << "No description available.";
165  out << "<br>\n";
167 }
168 
169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Float_t x
Definition: compare.C:6
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const XML_Char * name
Definition: expat.h:151
virtual ~G4hhIonisation()
void SetMinKinEnergy(G4double e)
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
G4int NumberOfBinsPerDecade() const
void SetEmModel(G4VEmModel *, G4int index=0)
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:707
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
G4double MaxKinEnergy() const
G4double GetPDGCharge() const
static const G4double emax
void SetDEDXBinning(G4int nbins)
G4double GetPDGMass() const
static constexpr double proton_mass_c2
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:440
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
static constexpr double electron_mass_c2
virtual void PrintInfo() override
virtual void ProcessDescription(std::ostream &) const override
const G4String & GetProcessName() const
Definition: G4VProcess.hh:411
void SetMaxKinEnergy(G4double e)
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override
static G4Electron * Electron()
Definition: G4Electron.cc:94
float bin[41]
Definition: plottest35.C:14
G4double MinKinEnergy() const
const G4ParticleDefinition * theParticle
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetBaseParticle(const G4ParticleDefinition *p)
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:714
G4int verboseLevel
Definition: G4VProcess.hh:371
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) override
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4GLOB_DLL std::ostream G4cout
G4hhIonisation(const G4String &name="hhIoni")
static G4EmParameters * Instance()
G4VEmFluctuationModel * flucModel
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual void ProcessDescription(std::ostream &outFile) const override