Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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$
28//
29//
30////////////////////////////////////////////////////////////////////////
31// Cerenkov Radiation Class Definition
32////////////////////////////////////////////////////////////////////////
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//
46////////////////////////////////////////////////////////////////////////
47
48#ifndef G4Cerenkov_h
49#define G4Cerenkov_h 1
50
51/////////////
52// Includes
53/////////////
54
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
77/////////////////////
78// Class Definition
79/////////////////////
80
81class G4Cerenkov : public G4VProcess
82{
83
84public:
85
86 ////////////////////////////////
87 // Constructors and Destructor
88 ////////////////////////////////
89
90 G4Cerenkov(const G4String& processName = "Cerenkov",
93
94 G4Cerenkov(const G4Cerenkov &right);
95
96private:
97
98 //////////////
99 // Operators
100 //////////////
101
102 G4Cerenkov& operator=(const G4Cerenkov &right);
103
104public:
105
106 ////////////
107 // Methods
108 ////////////
109
110 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
111 // Returns true -> 'is applicable', for all charged particles
112 // except short-lived particles.
113
114 G4double GetMeanFreePath(const G4Track& aTrack,
115 G4double ,
117 // Returns the discrete step limit and sets the 'StronglyForced'
118 // condition for the DoIt to be invoked at every step.
119
121 G4double ,
123 // Returns the discrete step limit and sets the 'StronglyForced'
124 // condition for the DoIt to be invoked at every step.
125
126 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
127 const G4Step& aStep);
128 // This is the method implementing the Cerenkov process.
129
130 // no operation in AtRestDoIt and AlongStepDoIt
132 const G4Track&,
133 G4double ,
134 G4double ,
135 G4double& ,
137 ) { return -1.0; };
138
140 const G4Track& ,
142 ) { return -1.0; };
143
144 // no operation in AtRestDoIt and AlongStepDoIt
146 const G4Track& ,
147 const G4Step&
148 ) {return 0;};
149
151 const G4Track& ,
152 const G4Step&
153 ) {return 0;};
154
155 void SetTrackSecondariesFirst(const G4bool state);
156 // If set, the primary particle tracking is interrupted and any
157 // produced Cerenkov photons are tracked next. When all have
158 // been tracked, the tracking of the primary resumes.
159
160 void SetMaxBetaChangePerStep(const G4double d);
161 // Set the maximum allowed change in beta = v/c in % (perCent)
162 // per step.
163
164 void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
165 // Set the maximum number of Cerenkov photons allowed to be
166 // generated during a tracking step. This is an average ONLY;
167 // the actual number will vary around this average. If invoked,
168 // the maximum photon stack will roughly be of the size set.
169 // If not called, the step is not limited by the number of
170 // photons generated.
171
173 // Returns the address of the physics table.
174
175 void DumpPhysicsTable() const;
176 // Prints the physics table.
177
178private:
179
180 void BuildThePhysicsTable();
181
182 /////////////////////
183 // Helper Functions
184 /////////////////////
185
186 G4double GetAverageNumberOfPhotons(const G4double charge,
187 const G4double beta,
188 const G4Material *aMaterial,
189 G4MaterialPropertyVector* Rindex) const;
190
191 ///////////////////////
192 // Class Data Members
193 ///////////////////////
194
195protected:
196
198 // A Physics Table can be either a cross-sections table or
199 // an energy table (or can be used for other specific
200 // purposes).
201
202private:
203
204 G4bool fTrackSecondariesFirst;
205 G4double fMaxBetaChange;
206 G4int fMaxPhotons;
207};
208
209////////////////////
210// Inline methods
211////////////////////
212
213inline
215{
216 if (aParticleType.GetParticleName() == "chargedgeantino") return false;
217 if (aParticleType.IsShortLived()) return false;
218
219 return (aParticleType.GetPDGCharge() != 0);
220}
221
222inline
224{
225 fTrackSecondariesFirst = state;
226}
227
228inline
230{
231 fMaxBetaChange = value*CLHEP::perCent;
232}
233
234inline
236{
237 fMaxPhotons = NumPhotons;
238}
239
240inline
242{
243 G4int PhysicsTableSize = thePhysicsTable->entries();
245
246 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
247 {
249 v->DumpValues();
250 }
251}
252
253inline
255{
256 return thePhysicsTable;
257}
258
259#endif /* G4Cerenkov_h */
G4ForceCondition
G4GPILSelection
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.hh:229
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:197
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Cerenkov.cc:134
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.hh:223
void DumpPhysicsTable() const
Definition: G4Cerenkov.hh:241
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
Definition: G4Cerenkov.hh:131
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
Definition: G4Cerenkov.hh:145
G4PhysicsTable * GetPhysicsTable() const
Definition: G4Cerenkov.hh:254
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.hh:214
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:461
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
Definition: G4Cerenkov.hh:150
G4Cerenkov(const G4Cerenkov &right)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
Definition: G4Cerenkov.hh:139
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:454
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.hh:235
G4double GetPDGCharge() const
const G4String & GetParticleName() const
size_t entries() const
Definition: G4Step.hh:78