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