Geant4  10.00.p02
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 // $Id: G4eBremsstrahlung.cc 75582 2013-11-04 12:13:01Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eBremsstrahlung
34 //
35 // Author: Michel Maire
36 //
37 // Creation date: 26.06.1996
38 //
39 // Modifications:
40 //
41 // 26-09-96 extension of the total crosssection above 100 GeV, M.Maire
42 // 1-10-96 new type G4OrderedTable; ComputePartialSumSigma(), M.Maire
43 // 16-10-96 DoIt() call to the non static GetEnergyCuts(), L.Urban
44 // 13-12-96 Sign corrected in grejmax and greject
45 // error definition of screenvar, L.Urban
46 // 20-03-97 new energy loss+ionisation+brems scheme, L.Urban
47 // 07-04-98 remove 'tracking cut' of the diffracted particle, MMa
48 // 13-08-98 new methods SetBining() PrintInfo()
49 // 03-03-99 Bug fixed in LPM effect, L.Urban
50 // 10-02-00 modifications , new e.m. structure, L.Urban
51 // 07-08-00 new cross section/en.loss parametrisation, LPM flag , L.Urban
52 // 21-09-00 corrections in the LPM implementation, L.Urban
53 // 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
54 // 09-08-01 new methods Store/Retrieve PhysicsTable (mma)
55 // 17-09-01 migration of Materials to pure STL (mma)
56 // 21-09-01 completion of RetrievePhysicsTable() (mma)
57 // 29-10-01 all static functions no more inlined (mma)
58 // 08-11-01 particleMass becomes a local variable
59 // 30-04-02 V.Ivanchenko update to new design
60 // 23-12-02 Change interface in order to move to cut per region (VI)
61 // 26-12-02 Secondary production moved to derived classes (VI)
62 // 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko)
63 // 08-08-03 STD substitute standard (V.Ivanchenko)
64 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
65 // 04-11-04 add gamma threshold (V.Ivanchenko)
66 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
67 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
68 // 22-05-06 Use gammaThreshold from manager (V.Ivantchenko)
69 // 15-01-07 use SetEmModel() from G4VEnergyLossProcess (mma)
70 // use RelEmModel above 1GeV (AS & VI)
71 // 13-11-08 reenable LPM switch (A.Schaelicke)
72 // -------------------------------------------------------------------
73 //
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
77 #include "G4eBremsstrahlung.hh"
78 #include "G4SystemOfUnits.hh"
79 #include "G4Gamma.hh"
80 #include "G4SeltzerBergerModel.hh"
82 #include "G4UnitsTable.hh"
83 #include "G4LossTableManager.hh"
84 
85 #include "G4ProductionCutsTable.hh"
86 #include "G4MaterialCutsCouple.hh"
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
90 using namespace std;
91 
93  G4VEnergyLossProcess(name),
94  isInitialised(false)
95 {
98  SetIonisation(false);
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102 
104 {}
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109 {
110  return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
115 void
117  const G4ParticleDefinition*)
118 {
119  if(!isInitialised) {
120 
121  // if (!EmModel(1)) { SetEmModel(new G4eBremsstrahlungModel(), 1); }
122  if (!EmModel(1)) { SetEmModel(new G4SeltzerBergerModel(), 1); }
123  if (!EmModel(2)) { SetEmModel(new G4eBremsstrahlungRelModel(), 2); }
124 
125  G4double energyLimit = 1*GeV;
126 
128  EmModel(1)->SetHighEnergyLimit(energyLimit);
129  EmModel(2)->SetLowEnergyLimit(energyLimit);
131 
133  AddEmModel(1, EmModel(1), fm);
134  AddEmModel(2, EmModel(2), fm);
135  isInitialised = true;
136  }
138  G4double eth = man->BremsstrahlungTh();
141 
142  // Only high energy model LMP flag is ON/OFF
143  EmModel(1)->SetLPMFlag(false);
144  EmModel(2)->SetLPMFlag(man->LPMFlag());
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
150 {
151  if(EmModel(1)) {
153  G4double eth = man->BremsstrahlungTh();
154  G4cout << " LPM flag: " << man->LPMFlag() << " for E > "
155  << EmModel(1)->HighEnergyLimit()/GeV << " GeV";
156  if(eth < DBL_MAX) {
157  G4cout << ", HighEnergyThreshold(GeV)= " << eth/GeV;
158  }
159  G4cout << G4endl;
160  }
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
static G4LossTableManager * Instance()
void SetIonisation(G4bool val)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4String name
Definition: TRTMaterials.hh:40
virtual void PrintInfo()
void SetSecondaryThreshold(G4double)
Definition: G4VEmModel.hh:725
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:196
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double BremsstrahlungTh() const
void SetSecondaryParticle(const G4ParticleDefinition *p)
#define fm
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:732
G4eBremsstrahlung(const G4String &name="eBrem")
void SetEmModel(G4VEmModel *, G4int index=1)
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)
#define DBL_MAX
Definition: templates.hh:83
G4VEmModel * EmModel(G4int index=1) const