Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4hBremsstrahlung.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: G4hBremsstrahlung
34//
35// Author: Vladimir Ivanchenko on base of model for muons
36//
37// Creation date: 01.03.2008
38//
39// Modifications:
40//
41//
42// -------------------------------------------------------------------
43//
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
47#include "G4hBremsstrahlung.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4Gamma.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55using namespace std;
56
59 theParticle(0),
60 theBaseParticle(0),
61 lowestKinEnergy(1.*GeV),
62 isInitialised(false)
63{
66 SetIonisation(false);
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72{}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77{
78 return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 110.0*MeV);
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
84 const G4Material*,
86{
87 return lowestKinEnergy;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
93 const G4ParticleDefinition* part,
95{
96 if(!isInitialised) {
97
98 isInitialised = true;
99
100 theParticle = part;
101
102 if (!EmModel()) { SetEmModel(new G4hBremsstrahlungModel()); }
103
104 G4VEmFluctuationModel* fm = 0;
107 AddEmModel(1, EmModel(), fm);
108 }
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112
114{}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117
@ fBremsstrahlung
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetPDGCharge() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetIonisation(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
G4hBremsstrahlung(const G4String &processName="hBrems")
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)
virtual ~G4hBremsstrahlung()
virtual void PrintInfo()
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut)
virtual G4bool IsApplicable(const G4ParticleDefinition &p)