Geant4  10.02
G4Cerenkov.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: G4Cerenkov.hh 85355 2014-10-28 09:58:59Z gcosmo $
28 //
29 //
31 // Cerenkov Radiation Class Definition
33 //
34 // File: G4Cerenkov.hh
35 // Description: Discrete Process - Generation of Cerenkov Photons
36 // Version: 2.0
37 // Created: 1996-02-21
38 // Author: Juliet Armstrong
39 // Updated: 2007-09-30 change inheritance to G4VDiscreteProcess
40 // 2005-07-28 add G4ProcessType to constructor
41 // 1999-10-29 add method and class descriptors
42 // 1997-04-09 by Peter Gumplinger
43 // > G4MaterialPropertiesTable; new physics/tracking scheme
44 // mail: gum@triumf.ca
45 //
47 
48 #ifndef G4Cerenkov_h
49 #define G4Cerenkov_h 1
50 
52 // Includes
54 
55 #include <CLHEP/Units/SystemOfUnits.h>
56 
57 #include "globals.hh"
58 #include "templates.hh"
59 #include "Randomize.hh"
60 #include "G4ThreeVector.hh"
61 #include "G4ParticleMomentum.hh"
62 #include "G4Step.hh"
63 #include "G4VProcess.hh"
64 #include "G4OpticalPhoton.hh"
65 #include "G4DynamicParticle.hh"
66 #include "G4Material.hh"
67 #include "G4PhysicsTable.hh"
71 
72 // Class Description:
73 // Discrete Process -- Generation of Cerenkov Photons.
74 // Class inherits publicly from G4VDiscreteProcess.
75 // Class Description - End:
76 
78 // Class Definition
80 
81 class G4Cerenkov : public G4VProcess
82 {
83 
84 public:
85 
87  // Constructors and Destructor
89 
90  G4Cerenkov(const G4String& processName = "Cerenkov",
92  ~G4Cerenkov();
93 
94  G4Cerenkov(const G4Cerenkov &right);
95 
96 private:
97 
99  // Operators
101 
102  G4Cerenkov& operator=(const G4Cerenkov &right);
103 
104 public:
105 
107  // Methods
109 
110  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
111  // Returns true -> 'is applicable', for all charged particles
112  // except short-lived particles.
113 
114  void BuildPhysicsTable(const G4ParticleDefinition& aParticleType);
115  // Build table at a right time
116 
117  G4double GetMeanFreePath(const G4Track& aTrack,
118  G4double ,
119  G4ForceCondition* );
120  // Returns the discrete step limit and sets the 'StronglyForced'
121  // condition for the DoIt to be invoked at every step.
122 
124  G4double ,
125  G4ForceCondition* );
126  // Returns the discrete step limit and sets the 'StronglyForced'
127  // condition for the DoIt to be invoked at every step.
128 
129  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
130  const G4Step& aStep);
131  // This is the method implementing the Cerenkov process.
132 
133  // no operation in AtRestDoIt and AlongStepDoIt
135  const G4Track&,
136  G4double ,
137  G4double ,
138  G4double& ,
140  ) { return -1.0; };
141 
143  const G4Track& ,
145  ) { return -1.0; };
146 
147  // no operation in AtRestDoIt and AlongStepDoIt
149  const G4Track& ,
150  const G4Step&
151  ) {return 0;};
152 
154  const G4Track& ,
155  const G4Step&
156  ) {return 0;};
157 
158  void SetTrackSecondariesFirst(const G4bool state);
159  // If set, the primary particle tracking is interrupted and any
160  // produced Cerenkov photons are tracked next. When all have
161  // been tracked, the tracking of the primary resumes.
162 
164  // Returns the boolean flag for tracking secondaries first.
165 
166  void SetMaxBetaChangePerStep(const G4double d);
167  // Set the maximum allowed change in beta = v/c in % (perCent)
168  // per step.
169 
171  // Returns the maximum allowed change in beta = v/c in % (perCent)
172 
173  void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
174  // Set the maximum number of Cerenkov photons allowed to be
175  // generated during a tracking step. This is an average ONLY;
176  // the actual number will vary around this average. If invoked,
177  // the maximum photon stack will roughly be of the size set.
178  // If not called, the step is not limited by the number of
179  // photons generated.
180 
182  // Returns the maximum number of Cerenkov photons allowed to be
183  // generated during a tracking step.
184 
186  // Returns the address of the physics table.
187 
188  void DumpPhysicsTable() const;
189  // Prints the physics table.
190 
191 private:
192 
193  void BuildThePhysicsTable();
194 
196  // Helper Functions
198 
200  const G4double beta,
201  const G4Material *aMaterial,
202  G4MaterialPropertyVector* Rindex) const;
203 
205  // Class Data Members
207 
208 protected:
209 
211  // A Physics Table can be either a cross-sections table or
212  // an energy table (or can be used for other specific
213  // purposes).
214 
215 private:
216 
220 };
221 
223 // Inline methods
225 
226 inline
228 {
229  return fTrackSecondariesFirst;
230 }
231 
232 inline
234 {
235  return fMaxBetaChange;
236 }
237 
238 inline
240 {
241  return fMaxPhotons;
242 }
243 
244 inline
246 {
247  G4int PhysicsTableSize = thePhysicsTable->entries();
249 
250  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
251  {
253  v->DumpValues();
254  }
255 }
256 
257 inline
259 {
260  return thePhysicsTable;
261 }
262 
263 #endif /* G4Cerenkov_h */
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:151
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
Definition: G4Cerenkov.hh:148
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
Definition: G4Cerenkov.hh:142
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:146
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
Definition: G4Cerenkov.cc:615
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:499
G4int GetMaxNumPhotonsPerStep() const
Definition: G4Cerenkov.hh:239
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:210
G4Cerenkov & operator=(const G4Cerenkov &right)
int G4int
Definition: G4Types.hh:78
G4bool GetTrackSecondariesFirst() const
Definition: G4Cerenkov.hh:227
G4int fMaxPhotons
Definition: G4Cerenkov.hh:219
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
Definition: G4Cerenkov.hh:153
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:492
bool G4bool
Definition: G4Types.hh:79
G4double fMaxBetaChange
Definition: G4Cerenkov.hh:218
G4PhysicsTable * GetPhysicsTable() const
Definition: G4Cerenkov.hh:258
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:156
Definition: G4Step.hh:76
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Cerenkov.cc:170
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition: G4Cerenkov.cc:100
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:161
G4ForceCondition
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
Definition: G4Cerenkov.hh:134
size_t entries() const
G4double GetMaxBetaChangePerStep() const
Definition: G4Cerenkov.hh:233
G4GPILSelection
void DumpPhysicsTable() const
Definition: G4Cerenkov.hh:245
void BuildThePhysicsTable()
Definition: G4Cerenkov.cc:393
G4ProcessType
G4bool fTrackSecondariesFirst
Definition: G4Cerenkov.hh:217
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:135