Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GeneralCerenkov.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// G4GeneralCerenkov
29//
30// Class description:
31// The Cerenkov process using model approach when an object G4VXRayModel
32// is assign to G4LogicalVolume. A model is fully responsible for
33// Cerenkov gamma production. This process class performing only tecnical
34// operation and interaction with Geant4 kernel.
35//
36// Created 25.05.2025 V.Ivanchenko on base of G4Cerenkov class
37//
38// --------------------------------------------------------------------
39
40#ifndef G4GeneralCerenkov_h
41#define G4GeneralCerenkov_h 1
42
43#include "globals.hh"
44#include "G4ForceCondition.hh"
45#include "G4LogicalVolume.hh"
46#include "G4VXRayModel.hh"
47#include "G4VDiscreteProcess.hh"
48
49#include <vector>
50
51class G4Material;
53class G4PhysicsTable;
54class G4Step;
55class G4Track;
57
59{
60 public:
61 explicit G4GeneralCerenkov(const G4String& processName = "Cerenkov",
63 ~G4GeneralCerenkov() override;
64
65 G4GeneralCerenkov(const G4GeneralCerenkov& right) = delete;
67
68 G4bool IsApplicable(const G4ParticleDefinition& aParticleType) override;
69
70 void PreparePhysicsTable(const G4ParticleDefinition& part) override;
71
72 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) override;
73
75 G4ForceCondition*) override;
76 // Returns the discrete step limit and sets the 'StronglyForced'
77 // condition for the DoIt to be invoked at every step.
78
80 const G4Step& aStep) override;
81 // This is the method implementing the Cerenkov process.
82
83 void AddModelForVolume(G4VXRayModel*, const G4String& nameLogVolume);
84 // explicit addition of a custom Cerenkov model to a logical volume
85
86 void DumpInfo() const override { ProcessDescription(G4cout); };
87 void ProcessDescription(std::ostream& out) const override;
88
90 G4ForceCondition*) override;
91
92 // Obsolete methods to be removed for the next major release
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 // Set the maximum allowed change in beta = v/c in % (perCent) per step.
101
102 void SetMaxNumPhotonsPerStep(const G4int NumPhotons);
103 // Set the maximum number of Cerenkov photons allowed to be generated during
104 // a tracking step. This is an average ONLY; the actual number will vary
105 // around this average. If invoked, the maximum photon stack will roughly be
106 // of the size set. If not called, the step is not limited by the number of
107 // photons generated.
108
109 void SetStackPhotons(const G4bool);
110 // Call by the user to set the flag for stacking the Cerenkov photons
111
113
114private:
115
116 const G4LogicalVolume* fCurrentLV{nullptr};
117 G4VXRayModel* fCurrentModel{nullptr};
118
119 G4double fMaxBetaChange{0.1};
120 G4double fBetaMin{1.0};
121 G4double fPreStepBeta{0.0};
122
123 G4int fMaxPhotons{100};
124
125 G4bool fStackingFlag{true};
126 G4bool fTrackSecondariesFirst{true};
127 G4bool isInitializer{false};
128 G4bool isPrepared{false};
129 G4bool isBuilt{false};
130
131 G4int secID{-1}; // creator modelID
132 G4int nModels{0};
133
134 // map includes logical volume pointer and index of the model
135 static std::vector<std::vector<const G4LogicalVolume*>* >* fLV;
136
137 // vector is used only at initialisation
138 // these models are destructed by G4LossTableManager
139 static std::vector<G4VXRayModel*>* fSharedModels;
140
141 // vector of names of logical volumes for master used for initilisation
142 // not filled for a worker thread
143 std::vector<G4String>* fLVNames{nullptr};
144
145 // models used in run time - they are thread local, are
146 // instantiated in worker thread, and are cloned from fSharedModels
147 // these models are destructed by G4LossTableManager
148 std::vector<G4VXRayModel*> fModels;
149
150 // buffer for X-Rays
151 std::vector<G4Track*> fSecondaries;
152};
153
154#endif
G4ForceCondition
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
G4GeneralCerenkov & operator=(const G4GeneralCerenkov &right)=delete
void PreparePhysicsTable(const G4ParticleDefinition &part) override
G4GeneralCerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
~G4GeneralCerenkov() override
void SetTrackSecondariesFirst(const G4bool state)
void AddModelForVolume(G4VXRayModel *, const G4String &nameLogVolume)
void SetMaxBetaChangePerStep(const G4double d)
void DumpInfo() const override
void SetStackPhotons(const G4bool)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *) override
G4GeneralCerenkov(const G4GeneralCerenkov &right)=delete
void ProcessDescription(std::ostream &out) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)