Geant4  10.00.p01
G4VTransitionRadiation.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 //
27 // $Id: G4VTransitionRadiation.hh 68037 2013-03-13 14:15:08Z gcosmo $
28 //
29 // G4VTransitionRadiation -- header file
30 //
31 // Generic process of transition radiation
32 //
33 // History:
34 // 29.02.04, V.Ivanchenko created
35 // 28.07.05, P.Gumplinger add G4ProcessType to constructor
36 
37 #ifndef G4VTransitionRadiation_h
38 #define G4VTransitionRadiation_h
39 
40 
41 #include "G4VDiscreteProcess.hh"
42 #include "G4Track.hh"
43 #include "G4ForceCondition.hh"
44 #include "globals.hh"
45 #include <vector>
46 
47 class G4Material;
48 class G4Region;
49 class G4VTRModel;
50 class G4particleDefinition;
51 
53 {
54 public:
55 
56 // Constructors
57  G4VTransitionRadiation( const G4String& processName = "TR",
59 
60 
61 // Destructor
62  virtual ~G4VTransitionRadiation() ;
63 
64  virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
65 
66  virtual G4double GetMeanFreePath(const G4Track& track, G4double,
68 
69  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
70  const G4Step& step);
71 
72  virtual void PrintInfoDefinition();
73  // Print out of the class parameters
74 
75  void SetRegion(const G4Region* reg);
76 
77  void SetModel(G4VTRModel* m);
78 
79 // private :
80 
81  void Clear();
82 
83  // hide assignment operator
86 
87  std::vector<const G4Material*> materials;
88  std::vector<G4double> steps;
89  std::vector<G4ThreeVector> normals;
90 
93  const G4Region* region;
95 
97 
100 
101 };
102 
104  const G4Track& track, G4double,
106 {
107  if(nSteps > 0) {
108  *condition = StronglyForced;
109  } else {
110  *condition = NotForced;
111  if(track.GetKineticEnergy()/track.GetDefinition()->GetPDGMass() + 1.0 > gammaMin &&
112  track.GetVolume()->GetLogicalVolume()->GetRegion() == region) {
113  *condition = StronglyForced;
114  }
115  }
116  return DBL_MAX; // so TR doesn't limit mean free path
117 }
118 
119 
120 #endif // G4VTransitionRadiation_h
G4double condition(const G4ErrorSymMatrix &m)
G4ParticleDefinition * GetDefinition() const
CLHEP::Hep3Vector G4ThreeVector
G4Region * GetRegion() const
int G4int
Definition: G4Types.hh:78
static const G4double reg
G4VTransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4double GetKineticEnergy() const
std::vector< G4ThreeVector > normals
bool G4bool
Definition: G4Types.hh:79
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
Definition: G4Step.hh:76
virtual G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition)
G4VTransitionRadiation & operator=(const G4VTransitionRadiation &right)
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
G4VPhysicalVolume * GetVolume() const
static const double m
Definition: G4SIunits.hh:110
double G4double
Definition: G4Types.hh:76
std::vector< const G4Material * > materials
G4ForceCondition
std::vector< G4double > steps
#define DBL_MAX
Definition: templates.hh:83
void SetRegion(const G4Region *reg)
G4ProcessType