Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28//
29////////////////////////////////////////////////////////////////////////
30// Cerenkov Radiation Class Definition
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4Cerenkov.hh
34// Description: Discrete Process - Generation of Cerenkov Photons
35// Version: 2.0
36// Created: 1996-02-21
37// Author: Juliet Armstrong
38// Updated: 2007-09-30 change inheritance to G4VDiscreteProcess
39// 2005-07-28 add G4ProcessType to constructor
40// 1999-10-29 add method and class descriptors
41// 1997-04-09 by Peter Gumplinger
42// > G4MaterialPropertiesTable; new physics/tracking scheme
43//
44////////////////////////////////////////////////////////////////////////
45
46#ifndef G4Cerenkov_h
47#define G4Cerenkov_h 1
48
50
51#include "globals.hh"
52#include "templates.hh"
53#include "Randomize.hh"
54#include "G4ThreeVector.hh"
55#include "G4ParticleMomentum.hh"
56#include "G4Step.hh"
57#include "G4VProcess.hh"
58#include "G4OpticalPhoton.hh"
59#include "G4DynamicParticle.hh"
60#include "G4Material.hh"
61#include "G4PhysicsTable.hh"
65
66class G4Cerenkov : public G4VProcess
67{
68 public:
69 explicit G4Cerenkov(const G4String& processName = "Cerenkov",
72
73 explicit G4Cerenkov(const G4Cerenkov& right);
74
75 private:
76 G4Cerenkov& operator=(const G4Cerenkov& right) = delete;
77
78 public:
79 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
80 // Returns true -> 'is applicable', for all charged particles
81 // except short-lived particles.
82
83 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
84 // Build table at a right time
85
86 void PreparePhysicsTable(const G4ParticleDefinition& part) override;
87 void Initialise();
88
90 // Returns the discrete step limit and sets the 'StronglyForced'
91 // condition for the DoIt to be invoked at every step.
92
94 G4ForceCondition*) override;
95 // Returns the discrete step limit and sets the 'StronglyForced'
96 // condition for the DoIt to be invoked at every step.
97
99 const G4Step& aStep) override;
100 // This is the method implementing the Cerenkov process.
101
102 // no operation in AtRestDoIt and AlongStepDoIt
104 const G4Track&, G4double, G4double, G4double&, G4GPILSelection*) override
105 {
106 return -1.0;
107 };
108
110 const G4Track&, G4ForceCondition*) override
111 {
112 return -1.0;
113 };
114
115 // no operation in AtRestDoIt and AlongStepDoIt
116 virtual G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&) override
117 {
118 return nullptr;
119 };
120
122 const G4Step&) override
123 {
124 return nullptr;
125 };
126
127 void SetTrackSecondariesFirst(const G4bool state);
128 // If set, the primary particle tracking is interrupted and any
129 // produced Cerenkov photons are tracked next. When all have
130 // been tracked, the tracking of the primary resumes.
131
133 // Returns the boolean flag for tracking secondaries first.
134
135 void SetMaxBetaChangePerStep(const G4double d);
136 // Set the maximum allowed change in beta = v/c in % (perCent) per step.
137
139 // Returns the maximum allowed change in beta = v/c in % (perCent)
140
141 void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
142 // Set the maximum number of Cerenkov photons allowed to be generated during
143 // a tracking step. This is an average ONLY; the actual number will vary
144 // around this average. If invoked, the maximum photon stack will roughly be
145 // of the size set. If not called, the step is not limited by the number of
146 // photons generated.
147
149 // Returns the maximum number of Cerenkov photons allowed to be
150 // generated during a tracking step.
151
152 void SetStackPhotons(const G4bool);
153 // Call by the user to set the flag for stacking the scint. photons
154
155 G4bool GetStackPhotons() const;
156 // Return the boolean for whether or not the scint. photons are stacked
157
158 G4int GetNumPhotons() const;
159 // Returns the current number of scint. photons (after PostStepDoIt)
160
162 // Returns the address of the physics table.
163
164 void DumpPhysicsTable() const;
165 // Prints the physics table.
166
167 G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta,
168 const G4Material* aMaterial,
169 G4MaterialPropertyVector* Rindex) const;
170
171 protected:
173
174 private:
175 G4bool fTrackSecondariesFirst;
176 G4double fMaxBetaChange;
177 G4int fMaxPhotons;
178
179 G4bool fStackingFlag;
180
181 G4int fNumPhotons;
182};
183
185{
186 return fTrackSecondariesFirst;
187}
188
190{
191 return fMaxBetaChange;
192}
193
194inline G4int G4Cerenkov::GetMaxNumPhotonsPerStep() const { return fMaxPhotons; }
195
196inline void G4Cerenkov::SetStackPhotons(const G4bool stackingFlag)
197{
198 fStackingFlag = stackingFlag;
199}
200
201inline G4bool G4Cerenkov::GetStackPhotons() const { return fStackingFlag; }
202
203inline G4int G4Cerenkov::GetNumPhotons() const { return fNumPhotons; }
204
206{
207 return thePhysicsTable;
208}
209
210#endif /* G4Cerenkov_h */
G4ForceCondition
G4GPILSelection
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
Definition: G4Cerenkov.hh:103
void Initialise()
Definition: G4Cerenkov.cc:113
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:572
G4bool GetStackPhotons() const
Definition: G4Cerenkov.hh:201
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:172
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition: G4Cerenkov.cc:391
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Cerenkov.cc:199
G4bool GetTrackSecondariesFirst() const
Definition: G4Cerenkov.hh:184
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:566
void DumpPhysicsTable() const
Definition: G4Cerenkov.cc:584
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
Definition: G4Cerenkov.hh:121
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
Definition: G4Cerenkov.cc:492
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:124
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
Definition: G4Cerenkov.hh:116
void PreparePhysicsTable(const G4ParticleDefinition &part) override
Definition: G4Cerenkov.cc:378
G4double GetMaxBetaChangePerStep() const
Definition: G4Cerenkov.hh:189
G4int GetNumPhotons() const
Definition: G4Cerenkov.hh:203
G4int GetMaxNumPhotonsPerStep() const
Definition: G4Cerenkov.hh:194
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
Definition: G4Cerenkov.hh:109
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:102
G4PhysicsTable * GetPhysicsTable() const
Definition: G4Cerenkov.hh:205
void SetStackPhotons(const G4bool)
Definition: G4Cerenkov.hh:196
G4Cerenkov(const G4Cerenkov &right)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:384
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:578
Definition: G4Step.hh:62