| 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: G4eBremsstrahlungModel.hh,v 1.25 2008/11/13 19:28:58 vnivanch Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-02 $
|
|---|
| 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class header file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4eBremsstrahlungModel
|
|---|
| 35 | //
|
|---|
| 36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code
|
|---|
| 37 | //
|
|---|
| 38 | // Creation date: 07.01.2002
|
|---|
| 39 | //
|
|---|
| 40 | // Modifications:
|
|---|
| 41 | //
|
|---|
| 42 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
|
|---|
| 43 | // 24-01-03 Make models region aware (V.Ivanchenko)
|
|---|
| 44 | // 13-02-03 Add name (V.Ivanchenko)
|
|---|
| 45 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
|
|---|
| 46 | // 07-02-06 public function ComputeCrossSectionPerAtom() (mma)
|
|---|
| 47 | //
|
|---|
| 48 | //
|
|---|
| 49 | // Class Description:
|
|---|
| 50 | //
|
|---|
| 51 | // Implementation of energy loss for gamma emission by electrons and
|
|---|
| 52 | // positrons
|
|---|
| 53 |
|
|---|
| 54 | // -------------------------------------------------------------------
|
|---|
| 55 | //
|
|---|
| 56 |
|
|---|
| 57 | #ifndef G4eBremsstrahlungModel_h
|
|---|
| 58 | #define G4eBremsstrahlungModel_h 1
|
|---|
| 59 |
|
|---|
| 60 | #include "G4VEmModel.hh"
|
|---|
| 61 |
|
|---|
| 62 | class G4Element;
|
|---|
| 63 | class G4ParticleChangeForLoss;
|
|---|
| 64 |
|
|---|
| 65 | class G4eBremsstrahlungModel : public G4VEmModel
|
|---|
| 66 | {
|
|---|
| 67 |
|
|---|
| 68 | public:
|
|---|
| 69 |
|
|---|
| 70 | G4eBremsstrahlungModel(const G4ParticleDefinition* p = 0,
|
|---|
| 71 | const G4String& nam = "eBrem");
|
|---|
| 72 |
|
|---|
| 73 | virtual ~G4eBremsstrahlungModel();
|
|---|
| 74 |
|
|---|
| 75 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
|
|---|
| 76 |
|
|---|
| 77 | G4double MinEnergyCut(const G4ParticleDefinition*,
|
|---|
| 78 | const G4MaterialCutsCouple*);
|
|---|
| 79 |
|
|---|
| 80 | virtual G4double ComputeDEDXPerVolume(const G4Material*,
|
|---|
| 81 | const G4ParticleDefinition*,
|
|---|
| 82 | G4double kineticEnergy,
|
|---|
| 83 | G4double cutEnergy);
|
|---|
| 84 |
|
|---|
| 85 | virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
|
|---|
| 86 | G4double tkin,
|
|---|
| 87 | G4double Z, G4double,
|
|---|
| 88 | G4double cut,
|
|---|
| 89 | G4double maxE = DBL_MAX);
|
|---|
| 90 |
|
|---|
| 91 | virtual G4double CrossSectionPerVolume(const G4Material*,
|
|---|
| 92 | const G4ParticleDefinition*,
|
|---|
| 93 | G4double kineticEnergy,
|
|---|
| 94 | G4double cutEnergy,
|
|---|
| 95 | G4double maxEnergy);
|
|---|
| 96 |
|
|---|
| 97 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
|
|---|
| 98 | const G4MaterialCutsCouple*,
|
|---|
| 99 | const G4DynamicParticle*,
|
|---|
| 100 | G4double tmin,
|
|---|
| 101 | G4double maxEnergy);
|
|---|
| 102 |
|
|---|
| 103 | protected:
|
|---|
| 104 |
|
|---|
| 105 | inline G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
|
|---|
| 106 | G4double kineticEnergy);
|
|---|
| 107 |
|
|---|
| 108 | const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple);
|
|---|
| 109 |
|
|---|
| 110 | private:
|
|---|
| 111 |
|
|---|
| 112 | void SetParticle(const G4ParticleDefinition* p);
|
|---|
| 113 |
|
|---|
| 114 | G4double ComputeBremLoss(G4double Z, G4double tkin, G4double cut);
|
|---|
| 115 |
|
|---|
| 116 | G4double PositronCorrFactorLoss(G4double Z, G4double tkin, G4double cut);
|
|---|
| 117 |
|
|---|
| 118 | G4double PositronCorrFactorSigma(G4double Z, G4double tkin, G4double cut);
|
|---|
| 119 |
|
|---|
| 120 | G4DataVector* ComputePartialSumSigma(const G4Material* material,
|
|---|
| 121 | G4double tkin, G4double cut);
|
|---|
| 122 |
|
|---|
| 123 | G4double SupressionFunction(const G4Material* material, G4double tkin,
|
|---|
| 124 | G4double gammaEnergy);
|
|---|
| 125 |
|
|---|
| 126 | inline G4double ScreenFunction1(G4double ScreenVariable);
|
|---|
| 127 |
|
|---|
| 128 | inline G4double ScreenFunction2(G4double ScreenVariable);
|
|---|
| 129 |
|
|---|
| 130 | // hide assignment operator
|
|---|
| 131 | G4eBremsstrahlungModel & operator=(const G4eBremsstrahlungModel &right);
|
|---|
| 132 | G4eBremsstrahlungModel(const G4eBremsstrahlungModel&);
|
|---|
| 133 |
|
|---|
| 134 | protected:
|
|---|
| 135 |
|
|---|
| 136 | const G4ParticleDefinition* particle;
|
|---|
| 137 | G4ParticleDefinition* theGamma;
|
|---|
| 138 | G4ParticleChangeForLoss* fParticleChange;
|
|---|
| 139 |
|
|---|
| 140 | G4double minThreshold;
|
|---|
| 141 | G4bool isElectron;
|
|---|
| 142 |
|
|---|
| 143 | private:
|
|---|
| 144 |
|
|---|
| 145 | G4double highKinEnergy;
|
|---|
| 146 | G4double lowKinEnergy;
|
|---|
| 147 | G4double probsup;
|
|---|
| 148 | G4double MigdalConstant;
|
|---|
| 149 | G4double LPMconstant;
|
|---|
| 150 | G4bool isInitialised;
|
|---|
| 151 |
|
|---|
| 152 | std::vector<G4DataVector*> partialSumSigma;
|
|---|
| 153 |
|
|---|
| 154 | };
|
|---|
| 155 |
|
|---|
| 156 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 157 |
|
|---|
| 158 | inline G4double G4eBremsstrahlungModel::ScreenFunction1(G4double ScreenVariable)
|
|---|
| 159 |
|
|---|
| 160 | // compute the value of the screening function 3*PHI1 - PHI2
|
|---|
| 161 |
|
|---|
| 162 | {
|
|---|
| 163 | G4double screenVal;
|
|---|
| 164 |
|
|---|
| 165 | if (ScreenVariable > 1.)
|
|---|
| 166 | screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
|
|---|
| 167 | else
|
|---|
| 168 | screenVal = 42.392 - ScreenVariable* (7.796 - 1.961*ScreenVariable);
|
|---|
| 169 |
|
|---|
| 170 | return screenVal;
|
|---|
| 171 | }
|
|---|
| 172 |
|
|---|
| 173 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 174 |
|
|---|
| 175 | inline
|
|---|
| 176 | G4double G4eBremsstrahlungModel::ScreenFunction2(G4double ScreenVariable)
|
|---|
| 177 |
|
|---|
| 178 | // compute the value of the screening function 1.5*PHI1 - 0.5*PHI2
|
|---|
| 179 |
|
|---|
| 180 | {
|
|---|
| 181 | G4double screenVal;
|
|---|
| 182 |
|
|---|
| 183 | if (ScreenVariable > 1.)
|
|---|
| 184 | screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
|
|---|
| 185 | else
|
|---|
| 186 | screenVal = 41.734 - ScreenVariable* (6.484 - 1.250*ScreenVariable);
|
|---|
| 187 |
|
|---|
| 188 | return screenVal;
|
|---|
| 189 | }
|
|---|
| 190 |
|
|---|
| 191 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 192 |
|
|---|
| 193 | inline
|
|---|
| 194 | G4double G4eBremsstrahlungModel::MaxSecondaryEnergy(
|
|---|
| 195 | const G4ParticleDefinition*,
|
|---|
| 196 | G4double kineticEnergy)
|
|---|
| 197 | {
|
|---|
| 198 | return kineticEnergy;
|
|---|
| 199 | }
|
|---|
| 200 |
|
|---|
| 201 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 202 |
|
|---|
| 203 | #endif
|
|---|