source: trunk/source/processes/electromagnetic/xrays/include/G4TransitionRadiation.hh@ 1066

Last change on this file since 1066 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 5.2 KB
Line 
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//
27// $Id: G4TransitionRadiation.hh,v 1.9 2006/06/29 19:55:47 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30// G4TransitionRadiation -- header file
31//
32// Class for description of transition radiation generated
33// by charged particle crossed interface between material 1
34// and material 2 (1 -> 2). Transition radiation could be of kind:
35// - optical back
36// - optical forward
37// - X-ray forward (for relativistic case Tkin/mass >= 10^2)
38//
39// GEANT 4 class header file --- Copyright CERN 1995
40// CERB Geneva Switzerland
41//
42// for information related to this code, please, contact
43// CERN, CN Division, ASD Group
44// History:
45// 18.12.97, V. Grichine (Vladimir.Grichine@cern.ch)
46// 02.02.00, V.Grichine, new data fEnergy and fVarAngle for double
47// numerical integration in inherited classes
48// 03.06.03, V.Ivanchenko fix compilation warnings
49// 28.07.05, P.Gumplinger add G4ProcessType to constructor
50
51#ifndef G4TransitionRadiation_h
52#define G4TransitionRadiation_h
53
54
55#include "G4VDiscreteProcess.hh"
56#include "G4Material.hh"
57// #include "G4OpBoundaryProcess.hh"
58
59class G4TransitionRadiation : public G4VDiscreteProcess
60{
61public:
62
63// Constructors
64
65
66 G4TransitionRadiation( const G4String& processName = "TR",
67 G4ProcessType type = fElectromagnetic) ;
68
69
70// G4TransitionRadiation(const G4TransitionRadiation& right) ;
71
72// Destructor
73
74 virtual ~G4TransitionRadiation() ;
75
76// Operators
77// G4TransitionRadiation& operator=(const G4TransitionRadiation& right) ;
78// G4int operator==(const G4TransitionRadiation& right)const ;
79// G4int operator!=(const G4TransitionRadiation& right)const ;
80
81// Methods
82
83 G4bool IsApplicable(const G4ParticleDefinition& aParticleType)
84 {
85 return ( aParticleType.GetPDGCharge() != 0.0 );
86 }
87
88 G4double GetMeanFreePath(const G4Track&,
89 G4double,
90 G4ForceCondition* condition)
91 {
92 *condition = Forced;
93 return DBL_MAX; // so TR doesn't limit mean free path
94 }
95
96 G4VParticleChange* PostStepDoIt(const G4Track&,
97 const G4Step&)
98 {
99 ClearNumberOfInteractionLengthLeft();
100 return &aParticleChange;
101 }
102
103
104
105
106virtual
107G4double SpectralAngleTRdensity( G4double energy,
108 G4double varAngle ) const = 0 ;
109
110G4double IntegralOverEnergy( G4double energy1,
111 G4double energy2,
112 G4double varAngle ) const ;
113
114G4double IntegralOverAngle( G4double energy,
115 G4double varAngle1,
116 G4double varAngle2 ) const ;
117
118G4double AngleIntegralDistribution( G4double varAngle1,
119 G4double varAngle2 ) const ;
120
121G4double EnergyIntegralDistribution( G4double energy1,
122 G4double energy2 ) const ;
123
124
125
126// Access functions
127
128
129protected :
130
131G4int fMatIndex1 ; // index of the 1st material
132G4int fMatIndex2 ; // index of the 2nd material
133
134// private :
135
136G4double fGamma ;
137G4double fEnergy ;
138G4double fVarAngle ;
139
140// Local constants
141static const G4int fSympsonNumber ; // Accuracy of Sympson integration 10
142static const G4int fGammaNumber ; // = 15
143static const G4int fPointNumber ; // = 100
144
145G4double fMinEnergy ; // min TR energy
146G4double fMaxEnergy ; // max TR energy
147G4double fMaxTheta ; // max theta of TR quanta
148
149G4double fSigma1 ; // plasma energy Sq of matter1
150G4double fSigma2 ; // plasma energy Sq of matter2
151
152
153} ;
154
155#endif // G4TransitionRadiation_h
Note: See TracBrowser for help on using the repository browser.