Geant4 11.4.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// Cerenkov Radiation Class Definition
28////////////////////////////////////////////////////////////////////////
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//
41////////////////////////////////////////////////////////////////////////
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 "G4VDiscreteProcess.hh"
52
53#include <map>
54
55class G4Material;
57class G4PhysicsTable;
58class G4Step;
59class G4Track;
61
63{
64 public:
65 explicit G4Cerenkov(const G4String& processName = "Cerenkov",
67 ~G4Cerenkov() override;
68
69 G4Cerenkov(const G4Cerenkov& right) = delete;
70 G4Cerenkov& operator=(const G4Cerenkov& right) = delete;
71
72 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
73 // Returns true -> 'is applicable', for all charged particles
74 // except short-lived particles.
75
76 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
77 // Build table at a right time
78
79 void PreparePhysicsTable(const G4ParticleDefinition& part) override;
80
81 G4double GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*) override;
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 void SetTrackSecondariesFirst(const G4bool state);
95 // If set, the primary particle tracking is interrupted and any
96 // produced Cerenkov photons are tracked next. When all have
97 // been tracked, the tracking of the primary resumes.
98
100 // Returns the boolean flag for tracking secondaries first.
101
102 void SetMaxBetaChangePerStep(const G4double d);
103 // Set the maximum allowed change in beta = v/c in % (perCent) per step.
104
106 // Returns the maximum allowed change in beta = v/c in % (perCent)
107
108 void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
109 // Set the maximum number of Cerenkov photons allowed to be generated during
110 // a tracking step. This is an average ONLY; the actual number will vary
111 // around this average. If invoked, the maximum photon stack will roughly be
112 // of the size set. If not called, the step is not limited by the number of
113 // photons generated.
114
116 // Returns the maximum number of Cerenkov photons allowed to be
117 // generated during a tracking step.
118
119 void SetStackPhotons(const G4bool);
120 // Call by the user to set the flag for stacking the Cerenkov photons
121
122 G4bool GetStackPhotons() const;
123 // Return the boolean for whether or not the Cerenkov photons are stacked
124
125 G4int GetNumPhotons() const;
126 // Returns the current number of scint. photons (after PostStepDoIt)
127
129 // Returns the address of the physics table.
130
131 void DumpPhysicsTable() const;
132 // Prints the physics table.
133
134 G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta,
135 const G4Material* aMaterial,
136 G4MaterialPropertyVector* Rindex) const;
137
138 void DumpInfo() const override;
139 void ProcessDescription(std::ostream& out) const override;
140
142 // sets verbosity
143
144 protected:
146 std::map<std::size_t, std::size_t> fIndexMPT;
147
148 private:
149 void Initialise();
150
151 G4double fMaxBetaChange;
152
153 G4int fMaxPhotons;
154 G4int fNumPhotons;
155
156 G4bool fStackingFlag;
157 G4bool fTrackSecondariesFirst;
158
159 G4int secID = -1; // creator modelID
160};
161
163{
164 return fTrackSecondariesFirst;
165}
166
168{
169 return fMaxBetaChange;
170}
171
172inline G4int G4Cerenkov::GetMaxNumPhotonsPerStep() const { return fMaxPhotons; }
173
174inline G4bool G4Cerenkov::GetStackPhotons() const { return fStackingFlag; }
175
176inline G4int G4Cerenkov::GetNumPhotons() const { return fNumPhotons; }
177
179{
180 return thePhysicsTable;
181}
182
183#endif /* G4Cerenkov_h */
G4ForceCondition
G4PhysicsFreeVector G4MaterialPropertyVector
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void ProcessDescription(std::ostream &out) const override
void SetMaxBetaChangePerStep(const G4double d)
G4bool GetStackPhotons() const
G4PhysicsTable * thePhysicsTable
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4bool GetTrackSecondariesFirst() const
void SetTrackSecondariesFirst(const G4bool state)
void DumpPhysicsTable() const
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
void SetVerboseLevel(G4int)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
void PreparePhysicsTable(const G4ParticleDefinition &part) override
G4double GetMaxBetaChangePerStep() const
G4int GetNumPhotons() const
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition G4Cerenkov.cc:79
G4int GetMaxNumPhotonsPerStep() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4PhysicsTable * GetPhysicsTable() const
G4Cerenkov & operator=(const G4Cerenkov &right)=delete
void SetStackPhotons(const G4bool)
G4Cerenkov(const G4Cerenkov &right)=delete
void DumpInfo() const override
std::map< std::size_t, std::size_t > fIndexMPT
~G4Cerenkov() override
Definition G4Cerenkov.cc:95
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)