Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eBremsstrahlung.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4eBremsstrahlung
33//
34// Author: Michel Maire
35//
36// Creation date: 26.06.1996
37//
38// Modified by Michel Maire, Vladimir Ivanchenko, Andreas Schaelicke
39//
40// -------------------------------------------------------------------
41//
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
45#include "G4eBremsstrahlung.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4Gamma.hh"
50#include "G4UnitsTable.hh"
51#include "G4LossTableManager.hh"
52
55#include "G4EmParameters.hh"
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59using namespace std;
60
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
83void
86{
87 if(!isInitialised) {
89
90 G4double emax = param->MaxKinEnergy();
91 G4VEmFluctuationModel* fm = nullptr;
92
93 if (nullptr == EmModel(0)) { SetEmModel(new G4SeltzerBergerModel()); }
94 G4double energyLimit = std::min(EmModel(0)->HighEnergyLimit(), CLHEP::GeV);
95 EmModel(0)->SetHighEnergyLimit(energyLimit);
97 AddEmModel(1, EmModel(0), fm);
98
99 if(emax > energyLimit) {
100 if (nullptr == EmModel(1)) {
102 }
103 EmModel(1)->SetLowEnergyLimit(energyLimit);
104 EmModel(1)->SetHighEnergyLimit(emax);
106 AddEmModel(1, EmModel(1), fm);
107 }
108 isInitialised = true;
109 }
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
114void G4eBremsstrahlung::StreamProcessInfo(std::ostream& out) const
115{
116 if(nullptr != EmModel(0)) {
118 G4double eth = param->BremsstrahlungTh();
119 out << " LPM flag: " << param->LPM() << " for E > "
120 << EmModel(0)->HighEnergyLimit()/GeV << " GeV";
121 if(eth < DBL_MAX) {
122 out << ", VertexHighEnergyTh(GeV)= " << eth/GeV;
123 }
124 out << G4endl;
125 }
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
129
130void G4eBremsstrahlung::ProcessDescription(std::ostream& out) const
131{
132 out << " Bremsstrahlung";
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fBremsstrahlung
@ fEmTwoPeaks
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4double BremsstrahlungTh() const
G4double MaxKinEnergy() const
G4bool LPM() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetHighEnergyLimit(G4double)
void SetSecondaryThreshold(G4double)
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void ProcessDescription(std::ostream &outFile) const override
G4VEmModel * EmModel(std::size_t index=0) const
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionType(G4CrossSectionType val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
G4bool IsApplicable(const G4ParticleDefinition &p) final
void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
void ProcessDescription(std::ostream &) const override
G4eBremsstrahlung(const G4String &name="eBrem")
~G4eBremsstrahlung() override
void StreamProcessInfo(std::ostream &outFile) const override
#define DBL_MAX
Definition templates.hh:62