Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TransitionRadiation.hh
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// Class for description of transition radiation generated
27// by charged particle crossed interface between material 1
28// and material 2 (1 -> 2). Transition radiation could be of kind:
29// - optical back
30// - optical forward
31// - X-ray forward (for relativistic case Tkin/mass >= 10^2)
32//
33// GEANT 4 class header file --- Copyright CERN 1995
34//
35// History:
36// 18.12.97, V. Grichine (Vladimir.Grichine@cern.ch)
37// 02.02.00, V.Grichine, new data fEnergy and fVarAngle for double
38// numerical integration in inherited classes
39// 03.06.03, V.Ivanchenko fix compilation warnings
40// 28.07.05, P.Gumplinger add G4ProcessType to constructor
41
42#ifndef G4TransitionRadiation_h
43#define G4TransitionRadiation_h
44
45#include "globals.hh"
47#include "G4Step.hh"
48#include "G4Track.hh"
49#include "G4VDiscreteProcess.hh"
50#include "G4VParticleChange.hh"
51
53{
54 public:
55 explicit G4TransitionRadiation(const G4String& processName = "TR",
57
59
62
63 // Methods
64
65 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
66
68 G4ForceCondition* condition) override;
69
71 const G4Step&) override;
72
73 virtual void ProcessDescription(std::ostream&) const override;
74 virtual void DumpInfo() const override { ProcessDescription(G4cout); };
75
77 G4double varAngle) const = 0;
78
80 G4double varAngle) const;
81
83 G4double varAngle2) const;
84
86 G4double varAngle2) const;
87
89
90 protected:
91 // Local constants
92 // Accuracy of Sympson integration
93 static constexpr G4int fSympsonNumber = 100;
94 static constexpr G4int fGammaNumber = 15;
95 static constexpr G4int fPointNumber = 100;
96
100
101 G4double fMinEnergy; // min TR energy
102 G4double fMaxEnergy; // max TR energy
103 G4double fMaxTheta; // max theta of TR quanta
104
105 G4double fSigma1; // plasma energy Sq of matter1
106 G4double fSigma2; // plasma energy Sq of matter2
107
108 G4int fMatIndex1; // index of the 1st material
109 G4int fMatIndex2; // index of the 2nd material
110};
111
112#endif // G4TransitionRadiation_h
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:62
static constexpr G4int fSympsonNumber
G4double EnergyIntegralDistribution(G4double energy1, G4double energy2) const
G4TransitionRadiation(const G4TransitionRadiation &right)=delete
G4double IntegralOverAngle(G4double energy, G4double varAngle1, G4double varAngle2) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
virtual void DumpInfo() const override
virtual ~G4TransitionRadiation()
virtual G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const =0
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
virtual void ProcessDescription(std::ostream &) const override
G4double IntegralOverEnergy(G4double energy1, G4double energy2, G4double varAngle) const
static constexpr G4int fPointNumber
G4TransitionRadiation & operator=(const G4TransitionRadiation &right)=delete
static constexpr G4int fGammaNumber
G4double AngleIntegralDistribution(G4double varAngle1, G4double varAngle2) const
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override