Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Cerenkov.hh
이 파일의 문서화 페이지로 가기
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 108423 2018-02-13 11:18:13Z 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 
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  explicit G4Cerenkov(const G4String& processName = "Cerenkov",
92  ~G4Cerenkov();
93 
94  explicit G4Cerenkov(const G4Cerenkov &right);
95 
96 private:
97 
99  // Operators
101 
102  G4Cerenkov& operator=(const G4Cerenkov &right) = delete;
103 
104 public:
105 
107  // Methods
109 
110  G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
111  // Returns true -> 'is applicable', for all charged particles
112  // except short-lived particles.
113 
114  void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
115  // Build table at a right time
116 
117  G4double GetMeanFreePath(const G4Track& aTrack,
119  // Returns the discrete step limit and sets the 'StronglyForced'
120  // condition for the DoIt to be invoked at every step.
121 
123  G4double ,
124  G4ForceCondition* ) override;
125  // Returns the discrete step limit and sets the 'StronglyForced'
126  // condition for the DoIt to be invoked at every step.
127 
128  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
129  const G4Step& aStep) override;
130  // This is the method implementing the Cerenkov process.
131 
132  // no operation in AtRestDoIt and AlongStepDoIt
134  G4double ,
135  G4double ,
136  G4double& ,
138  ) override { return -1.0; };
139 
142  ) override { return -1.0; };
143 
144  // no operation in AtRestDoIt and AlongStepDoIt
145  virtual G4VParticleChange* AtRestDoIt(const G4Track& , const G4Step& )
146  override {return nullptr;};
147 
148  virtual G4VParticleChange* AlongStepDoIt(const G4Track& , const G4Step&)
149  override {return nullptr;};
150 
151  void SetTrackSecondariesFirst(const G4bool state);
152  // If set, the primary particle tracking is interrupted and any
153  // produced Cerenkov photons are tracked next. When all have
154  // been tracked, the tracking of the primary resumes.
155 
157  // Returns the boolean flag for tracking secondaries first.
158 
159  void SetMaxBetaChangePerStep(const G4double d);
160  // Set the maximum allowed change in beta = v/c in % (perCent)
161  // per step.
162 
164  // Returns the maximum allowed change in beta = v/c in % (perCent)
165 
166  void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
167  // Set the maximum number of Cerenkov photons allowed to be
168  // generated during a tracking step. This is an average ONLY;
169  // the actual number will vary around this average. If invoked,
170  // the maximum photon stack will roughly be of the size set.
171  // If not called, the step is not limited by the number of
172  // photons generated.
173 
175  // Returns the maximum number of Cerenkov photons allowed to be
176  // generated during a tracking step.
177 
178  void SetStackPhotons(const G4bool );
179  // Call by the user to set the flag for stacking the scint. photons
180 
181  G4bool GetStackPhotons() const;
182  // Return the boolean for whether or not the scint. photons are stacked
183 
184  G4int GetNumPhotons() const;
185  // Returns the current number of scint. photons (after PostStepDoIt)
186 
188  // Returns the address of the physics table.
189 
190  void DumpPhysicsTable() const;
191  // Prints the physics table.
192 
193 private:
194 
195  void BuildThePhysicsTable();
196 
198  // Helper Functions
200 
202  const G4double beta,
203  const G4Material *aMaterial,
204  G4MaterialPropertyVector* Rindex) const;
205 
207  // Class Data Members
209 
210 protected:
211 
213  // A Physics Table can be either a cross-sections table or
214  // an energy table (or can be used for other specific
215  // purposes).
216 
217 private:
218 
222 
224 
226 };
227 
229  // Inline methods
231 
232 inline
234 {
235  return fTrackSecondariesFirst;
236 }
237 
238 inline
240 {
241  return fMaxBetaChange;
242 }
243 
244 inline
246 {
247  return fMaxPhotons;
248 }
249 
250 inline
251 void G4Cerenkov::SetStackPhotons(const G4bool stackingFlag)
252 {
253  fStackingFlag = stackingFlag;
254 }
255 
256 inline
258 {
259  return fStackingFlag;
260 }
261 
262 inline
264 {
265  return fNumPhotons;
266 }
267 
268 inline
270 {
271  return thePhysicsTable;
272 }
273 
274 #endif /* G4Cerenkov_h */
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Cerenkov.cc:169
G4bool fStackingFlag
Definition: G4Cerenkov.hh:223
G4PhysicsTable * GetPhysicsTable() const
Definition: G4Cerenkov.hh:269
G4int fMaxPhotons
Definition: G4Cerenkov.hh:221
G4double GetMaxBetaChangePerStep() const
Definition: G4Cerenkov.hh:239
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
Definition: G4Cerenkov.hh:148
void SetStackPhotons(const G4bool)
Definition: G4Cerenkov.hh:251
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:150
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
Definition: G4Cerenkov.hh:140
G4double fMaxBetaChange
Definition: G4Cerenkov.hh:220
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:155
Double_t beta
void BuildThePhysicsTable()
Definition: G4Cerenkov.cc:391
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:145
G4int fNumPhotons
Definition: G4Cerenkov.hh:225
G4bool GetStackPhotons() const
Definition: G4Cerenkov.hh:257
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
G4ProcessType
G4GPILSelection
Float_t d
G4Cerenkov & operator=(const G4Cerenkov &right)=delete
Definition: G4Step.hh:76
G4bool fTrackSecondariesFirst
Definition: G4Cerenkov.hh:219
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition: G4Cerenkov.cc:492
G4int GetNumPhotons() const
Definition: G4Cerenkov.hh:263
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition: G4Cerenkov.cc:96
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
Definition: G4Cerenkov.hh:133
G4bool GetTrackSecondariesFirst() const
Definition: G4Cerenkov.hh:233
int G4int
Definition: G4Types.hh:78
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:160
void DumpPhysicsTable() const
Definition: G4Cerenkov.cc:683
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:485
G4ForceCondition
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
Definition: G4Cerenkov.cc:603
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
Definition: G4Cerenkov.hh:145
G4int GetMaxNumPhotonsPerStep() const
Definition: G4Cerenkov.hh:245
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:133