[819] | 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,v 1.48 2007/05/23 08:47:34 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-01-patch-02 $ |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4eBremsstrahlung |
---|
| 35 | // |
---|
| 36 | // Author: Michel Maire |
---|
| 37 | // |
---|
| 38 | // Creation date: 26.06.1996 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 26-09-96 extension of the total crosssection above 100 GeV, M.Maire |
---|
| 43 | // 1-10-96 new type G4OrderedTable; ComputePartialSumSigma(), M.Maire |
---|
| 44 | // 16-10-96 DoIt() call to the non static GetEnergyCuts(), L.Urban |
---|
| 45 | // 13-12-96 Sign corrected in grejmax and greject |
---|
| 46 | // error definition of screenvar, L.Urban |
---|
| 47 | // 20-03-97 new energy loss+ionisation+brems scheme, L.Urban |
---|
| 48 | // 07-04-98 remove 'tracking cut' of the diffracted particle, MMa |
---|
| 49 | // 13-08-98 new methods SetBining() PrintInfo() |
---|
| 50 | // 03-03-99 Bug fixed in LPM effect, L.Urban |
---|
| 51 | // 10-02-00 modifications , new e.m. structure, L.Urban |
---|
| 52 | // 07-08-00 new cross section/en.loss parametrisation, LPM flag , L.Urban |
---|
| 53 | // 21-09-00 corrections in the LPM implementation, L.Urban |
---|
| 54 | // 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation |
---|
| 55 | // 09-08-01 new methods Store/Retrieve PhysicsTable (mma) |
---|
| 56 | // 17-09-01 migration of Materials to pure STL (mma) |
---|
| 57 | // 21-09-01 completion of RetrievePhysicsTable() (mma) |
---|
| 58 | // 29-10-01 all static functions no more inlined (mma) |
---|
| 59 | // 08-11-01 particleMass becomes a local variable |
---|
| 60 | // 30-04-02 V.Ivanchenko update to new design |
---|
| 61 | // 23-12-02 Change interface in order to move to cut per region (VI) |
---|
| 62 | // 26-12-02 Secondary production moved to derived classes (VI) |
---|
| 63 | // 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko) |
---|
| 64 | // 08-08-03 STD substitute standard (V.Ivanchenko) |
---|
| 65 | // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko) |
---|
| 66 | // 04-11-04 add gamma threshold (V.Ivanchenko) |
---|
| 67 | // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko) |
---|
| 68 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
| 69 | // 22-05-06 Use gammaThreshold from manager (V.Ivantchenko) |
---|
| 70 | // 15-01-07 use SetEmModel() from G4VEnergyLossProcess (mma) |
---|
| 71 | // |
---|
| 72 | // ------------------------------------------------------------------- |
---|
| 73 | // |
---|
| 74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 75 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 76 | |
---|
| 77 | #include "G4eBremsstrahlung.hh" |
---|
| 78 | #include "G4Gamma.hh" |
---|
| 79 | #include "G4eBremsstrahlungModel.hh" |
---|
| 80 | #include "G4UniversalFluctuation.hh" |
---|
| 81 | #include "G4UnitsTable.hh" |
---|
| 82 | #include "G4LossTableManager.hh" |
---|
| 83 | |
---|
| 84 | #include "G4ProductionCutsTable.hh" |
---|
| 85 | #include "G4MaterialCutsCouple.hh" |
---|
| 86 | |
---|
| 87 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 88 | |
---|
| 89 | using namespace std; |
---|
| 90 | |
---|
| 91 | G4eBremsstrahlung::G4eBremsstrahlung(const G4String& name): |
---|
| 92 | G4VEnergyLossProcess(name), |
---|
| 93 | isInitialised(false) |
---|
| 94 | {} |
---|
| 95 | |
---|
| 96 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 97 | |
---|
| 98 | G4eBremsstrahlung::~G4eBremsstrahlung() |
---|
| 99 | {} |
---|
| 100 | |
---|
| 101 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 102 | |
---|
| 103 | void G4eBremsstrahlung::InitialiseEnergyLossProcess( |
---|
| 104 | const G4ParticleDefinition* p, |
---|
| 105 | const G4ParticleDefinition*) |
---|
| 106 | { |
---|
| 107 | if(!isInitialised) { |
---|
| 108 | particle = p; |
---|
| 109 | SetSecondaryParticle(G4Gamma::Gamma()); |
---|
| 110 | SetIonisation(false); |
---|
| 111 | if (!EmModel()) SetEmModel(new G4eBremsstrahlungModel()); |
---|
| 112 | EmModel()->SetLowEnergyLimit (100*eV); |
---|
| 113 | EmModel()->SetHighEnergyLimit(100*TeV); |
---|
| 114 | if (!FluctModel()) SetFluctModel(new G4UniversalFluctuation()); |
---|
| 115 | |
---|
| 116 | AddEmModel(1, EmModel(), FluctModel()); |
---|
| 117 | isInitialised = true; |
---|
| 118 | } |
---|
| 119 | G4LossTableManager* man = G4LossTableManager::Instance(); |
---|
| 120 | dynamic_cast<G4eBremsstrahlungModel*>(EmModel()) |
---|
| 121 | ->SetEnergyThreshold(man->BremsstrahlungTh()); |
---|
| 122 | dynamic_cast<G4eBremsstrahlungModel*>(EmModel()) |
---|
| 123 | ->SetLPMflag(man->LPMFlag()); |
---|
| 124 | } |
---|
| 125 | |
---|
| 126 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 127 | |
---|
| 128 | void G4eBremsstrahlung::PrintInfo() |
---|
| 129 | { |
---|
| 130 | if(EmModel()) { |
---|
| 131 | G4cout << " Total cross sections and sampling from " |
---|
| 132 | << EmModel()->GetName() << " model" |
---|
| 133 | << " (based on the EEDL data library) " |
---|
| 134 | << "\n Good description from 1 KeV to 100 GeV, " |
---|
| 135 | << "log scale extrapolation above 100 GeV." |
---|
| 136 | << " LPM flag " |
---|
| 137 | << dynamic_cast<G4eBremsstrahlungModel*>(EmModel())->LPMflag() |
---|
| 138 | << G4endl; |
---|
| 139 | G4double eth = dynamic_cast<G4eBremsstrahlungModel*>(EmModel())->EnergyThreshold(); |
---|
| 140 | if(eth < DBL_MIN) |
---|
| 141 | G4cout << " HighEnergyThreshold(GeV)= " << eth/GeV |
---|
| 142 | << G4endl; |
---|
| 143 | } |
---|
| 144 | } |
---|
| 145 | |
---|
| 146 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|