Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GeneralCerenkov Class Reference

#include <G4GeneralCerenkov.hh>

Inheritance diagram for G4GeneralCerenkov:

Public Member Functions

 G4GeneralCerenkov (const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
 ~G4GeneralCerenkov () override
 G4GeneralCerenkov (const G4GeneralCerenkov &right)=delete
G4GeneralCerenkovoperator= (const G4GeneralCerenkov &right)=delete
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
void PreparePhysicsTable (const G4ParticleDefinition &part) override
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double, G4ForceCondition *) override
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
void AddModelForVolume (G4VXRayModel *, const G4String &nameLogVolume)
void DumpInfo () const override
void ProcessDescription (std::ostream &out) const override
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *) override
void SetTrackSecondariesFirst (const G4bool state)
void SetMaxBetaChangePerStep (const G4double d)
void SetMaxNumPhotonsPerStep (const G4int NumPhotons)
void SetStackPhotons (const G4bool)
void SetVerboseLevel (G4int)
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 G4VProcess (const G4VProcess &right)
virtual ~G4VProcess ()
G4VProcessoperator= (const G4VProcess &)=delete
G4bool operator== (const G4VProcess &right) const
G4bool operator!= (const G4VProcess &right) const
G4double GetCurrentInteractionLength () const
void SetPILfactor (G4double value)
G4double GetPILfactor () const
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
const G4StringGetProcessName () const
G4ProcessType GetProcessType () const
void SetProcessType (G4ProcessType)
G4int GetProcessSubType () const
void SetProcessSubType (G4int)
virtual const G4VProcessGetCreatorProcess () const
virtual void StartTracking (G4Track *)
virtual void EndTracking ()
virtual void SetProcessManager (const G4ProcessManager *)
virtual const G4ProcessManagerGetProcessManager ()
virtual void ResetNumberOfInteractionLengthLeft ()
G4double GetNumberOfInteractionLengthLeft () const
G4double GetTotalNumberOfInteractionLengthTraversed () const
G4bool isAtRestDoItIsEnabled () const
G4bool isAlongStepDoItIsEnabled () const
G4bool isPostStepDoItIsEnabled () const
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const
virtual void SetMasterProcess (G4VProcess *masterP)
const G4VProcessGetMasterProcess () const
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()
Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
G4VParticleChangepParticleChange = nullptr
G4ParticleChange aParticleChange
G4double theNumberOfInteractionLengthLeft = -1.0
G4double currentInteractionLength = -1.0
G4double theInitialNumberOfInteractionLength = -1.0
G4String theProcessName
G4String thePhysicsTableFileName
G4ProcessType theProcessType = fNotDefined
G4int theProcessSubType = -1
G4double thePILfactor = 1.0
G4int verboseLevel = 0
G4bool enableAtRestDoIt = true
G4bool enableAlongStepDoIt = true
G4bool enablePostStepDoIt = true

Detailed Description

Definition at line 58 of file G4GeneralCerenkov.hh.

Constructor & Destructor Documentation

◆ G4GeneralCerenkov() [1/2]

G4GeneralCerenkov::G4GeneralCerenkov ( const G4String & processName = "Cerenkov",
G4ProcessType type = fElectromagnetic )
explicit

Definition at line 47 of file G4GeneralCerenkov.cc.

48 : G4VDiscreteProcess(nam, type)
49{
50 secID = G4PhysicsModelCatalog::GetModelID("model_Cerenkov");
52 if (nullptr == fLV) {
53 // initialise static data
54 fSharedModels = new std::vector<G4VXRayModel*>;
55 fLV = new std::vector<std::vector<const G4LogicalVolume*>* >;
56 fLVNames = new std::vector<G4String>;
57
58 isInitializer = true;
59 }
60}
@ fCerenkov
static G4int GetModelID(const G4int modelIndex)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
void SetProcessSubType(G4int)

Referenced by G4GeneralCerenkov(), and operator=().

◆ ~G4GeneralCerenkov()

G4GeneralCerenkov::~G4GeneralCerenkov ( )
override

Definition at line 63 of file G4GeneralCerenkov.cc.

64{
65 if (isInitializer) {
66 delete fSharedModels;
67 fSharedModels = nullptr;
68 delete fLVNames;
69 for (auto const & p : *fLV) {
70 delete p;
71 }
72 delete fLV;
73 fLV = nullptr;
74 }
75}

◆ G4GeneralCerenkov() [2/2]

G4GeneralCerenkov::G4GeneralCerenkov ( const G4GeneralCerenkov & right)
delete

Member Function Documentation

◆ AddModelForVolume()

void G4GeneralCerenkov::AddModelForVolume ( G4VXRayModel * model,
const G4String & nameLogVolume )

Definition at line 84 of file G4GeneralCerenkov.cc.

86{
87 if (isPrepared || !isInitializer || nullptr == model) {
89 G4String nam;
90 if (model != nullptr) { nam = model->GetName(); }
91 ed << " Attempt to add Cerenkov model <" << nam << "> for LogicalVolume "
92 << nameLogVolume << " is failed!\n isPrepared:" << isPrepared
93 << " isInitilizer:" << isInitializer;
94 G4Exception("G4GeneralCerenkov::AddModelForVolume", "em0304",
95 FatalException, ed, "");
96 return;
97 }
98 fSharedModels->push_back(model);
99 fLVNames->push_back(nameLogVolume);
100}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
const G4String & GetName() const

◆ BuildPhysicsTable()

void G4GeneralCerenkov::BuildPhysicsTable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 181 of file G4GeneralCerenkov.cc.

182{
183 // Initialisation of models is done only once
184 if (isBuilt) { return; }
185 isBuilt = true;
186 if (nModels == 0) { return; }
187
188 // worker thread
189 if (!isInitializer) {
190 // worker initialisation - clone master models
191 fModels.reserve(nModels);
192 for (G4int i=0; i<nModels; ++i) {
193 auto newmod = new G4VXRayModel(*((*fSharedModels)[i]));
194 fModels.push_back(newmod);
195 G4double b = newmod->Initialise((*fLV)[i]);
196 fBetaMin = std::min(fBetaMin, b);
197 }
198 } else if (verboseLevel > 0) {
199 // needed for printout
200 std::size_t nn = 0;
201 for (G4int i=0; i<nModels; ++i) {
202 G4double b = (*fSharedModels)[i]->Initialise((*fLV)[i]);
203 fBetaMin = std::min(fBetaMin, b);
204 nn += ((*fLV)[i])->size();
205 }
206 G4long pres = G4cout.precision();
207 G4cout.precision(6);
208 G4cout << " " << GetProcessName() << std::setw(20) << " fBetaMin=" << fBetaMin
209 << " fMaxBetaChange=" << fMaxBetaChange << G4endl;
210 G4cout << std::setw(20) << "fMaxNphot=" << fMaxPhotons
211 << " Nlv=" << nn << " fStackingFlag:" << fStackingFlag
212 << " fTrackSecondariesFirst:" << fTrackSecondariesFirst
213 << G4endl;
214 for (G4int i=0; i<nModels; ++i) {
215 G4int n = (G4int)((*fLV)[i]->size());
216 G4cout << std::setw(10) << (*fSharedModels)[i]->GetName()
217 << std::setw(30) << "Nvolumes=" << n << " Volumes:" << G4endl;
218 G4cout << std::setw(12);
219 for (G4int j=0; j<n; ++j) {
220 G4cout << (*((*fLV)[i]))[j]->GetName() << " ";
221 if (0 != j && (j/5)*5 == j) {
222 G4cout << G4endl << std::setw(12);
223 }
224 }
225 G4cout << G4endl;
226 }
227 G4cout.precision(pres);
228 }
229}
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
const G4String & GetProcessName() const

◆ DumpInfo()

void G4GeneralCerenkov::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4GeneralCerenkov.hh.

void ProcessDescription(std::ostream &out) const override

◆ GetMeanFreePath()

G4double G4GeneralCerenkov::GetMeanFreePath ( const G4Track & ,
G4double ,
G4ForceCondition *  )
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 345 of file G4GeneralCerenkov.cc.

346{
347 return DBL_MAX;
348}
#define DBL_MAX
Definition templates.hh:62

◆ IsApplicable()

G4bool G4GeneralCerenkov::IsApplicable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 78 of file G4GeneralCerenkov.cc.

79{
80 return true;
81}

◆ operator=()

G4GeneralCerenkov & G4GeneralCerenkov::operator= ( const G4GeneralCerenkov & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4GeneralCerenkov::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 272 of file G4GeneralCerenkov.cc.

274{
275 aParticleChange.Initialize(aTrack);
276 if (fCurrentModel == nullptr) { return &aParticleChange; }
277 fCurrentModel->SampleXRays(fSecondaries, aStep);
278 if (!fSecondaries.empty()) {
279
280 // X-rays
281 auto touch = aStep.GetPreStepPoint()->GetTouchableHandle();
282 G4int parent = aTrack.GetTrackID();
283 for (auto & t : fSecondaries) {
284 t->SetTouchableHandle(touch);
285 t->SetParentID(parent);
286 t->SetCreatorModelID(secID);
287 aParticleChange.AddSecondary(t);
288 }
289 fSecondaries.clear();
290
291 // primary track suspended
292 if (fTrackSecondariesFirst && aTrack.GetTrackStatus() == fAlive) {
293 aParticleChange.ProposeTrackStatus(fSuspend);
294 }
295 }
296 return &aParticleChange;
297}
@ fSuspend
@ fAlive
const G4TouchableHandle & GetTouchableHandle() const
G4StepPoint * GetPreStepPoint() const
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
G4ParticleChange aParticleChange

◆ PostStepGetPhysicalInteractionLength()

G4double G4GeneralCerenkov::PostStepGetPhysicalInteractionLength ( const G4Track & aTrack,
G4double ,
G4ForceCondition * cond )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 233 of file G4GeneralCerenkov.cc.

236{
237 *cond = NotForced;
238 G4double limit = DBL_MAX;
239 auto const dp = aTrack.GetDynamicParticle();
240 if (dp->GetCharge() == 0.0) { return limit; }
241 fCurrentModel = nullptr;
242 fPreStepBeta = dp->GetBeta();
243 if (fPreStepBeta <= fBetaMin) { return limit; }
244
245 auto volume = aTrack.GetVolume();
246 if (nullptr == volume) { return limit; }
247
248 fCurrentLV = volume->GetLogicalVolume();
249 auto const MPT = fCurrentLV->GetMaterial()->GetMaterialPropertiesTable();
250 if (nullptr == MPT) { return limit; }
251
252 G4bool ok{false};
253 for (G4int i=0; i<nModels; ++i) {
254 auto const v = (*fLV)[i];
255 std::size_t nn = v->size();
256 for (std::size_t j = 0; j < nn; ++j) {
257 if ((*v)[j] == fCurrentLV) {
258 fCurrentModel = fModels[i];
259 if (fCurrentModel->StepLimit(j, aTrack, fPreStepBeta, limit)) {
260 *cond = StronglyForced;
261 }
262 ok = true;
263 break;
264 }
265 }
266 if (ok) { break; }
267 }
268 return limit;
269}
@ StronglyForced
@ NotForced
bool G4bool
Definition G4Types.hh:86
G4double GetBeta() const
G4Material * GetMaterial() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4VPhysicalVolume * GetVolume() const
const G4DynamicParticle * GetDynamicParticle() const
G4LogicalVolume * GetLogicalVolume() const

◆ PreparePhysicsTable()

void G4GeneralCerenkov::PreparePhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 103 of file G4GeneralCerenkov.cc.

104{
105 // definition of models is done only once
106 if (isPrepared) { return; }
107 isPrepared = true;
108
109 const G4OpticalParameters* params = G4OpticalParameters::Instance();
110 fMaxBetaChange = params->GetCerenkovMaxBetaChange();
111 fMaxPhotons = params->GetCerenkovMaxPhotonsPerStep();
112 fStackingFlag = params->GetCerenkovStackPhotons();
113 fTrackSecondariesFirst = params->GetCerenkovTrackSecondariesFirst();
115
116 auto nmod = fSharedModels->size();
117 if (0 == nmod) {
118 // the default model is added without association with a logical volume
119 G4VXRayModel* mod = new G4StandardCerenkovModel();
120 fSharedModels->push_back(mod);
121 nmod = 1;
122 }
123
124 fSecondaries.reserve(fMaxPhotons);
125 nModels = (G4int)nmod;
126
127 // fill static data structures
128 if (isInitializer) {
129 const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
130 const auto & modAndVol = params->ActiveVolumes();
131
132 // preparation of logical volume vector per model
133 fLV->reserve(nmod);
134 for (G4int i=0; i<nModels; ++i) {
135 auto v = new std::vector<const G4LogicalVolume*>;
136 fLV->push_back(v);
137 }
138
139 for (auto const & lv : *lvs) {
140 // only volumes with material property defined are considered
141 auto const MPT = lv->GetMaterial()->GetMaterialPropertiesTable();
142 if (nullptr == MPT) { continue; }
143
144 const G4String& lvname = lv->GetName();
145 G4bool ok{false};
146 // search for the default model in the list
147 if (!modAndVol.empty()) {
148 for (auto const & it : modAndVol) {
149 if (it.second == lvname) {
150 if (kCerenkovDefault == it.first) {
151 (*fLV)[0]->push_back(lv);
152 fLVNames->push_back(lvname);
153 ok = true;
154 break;
155 }
156 }
157 }
158 }
159 if (ok) { continue; }
160
161 // search in external models
162 for (G4int i=0; i<nModels; ++i) {
163 if (lvname == (*fLVNames)[i]) {
164 (*fLV)[i]->push_back(lv);
165 ok = true;
166 break;
167 }
168 }
169 if (ok) { continue; }
170
171 // temporary for backward compatibility search for a RINDEX of the volume
172 if (nullptr != MPT->GetProperty(kRINDEX)) {
173 (*fLV)[0]->push_back(lv);
174 fLVNames->push_back(lvname);
175 }
176 }
177 }
178}
@ kCerenkovDefault
static G4LogicalVolumeStore * GetInstance()
G4int GetCerenkovVerboseLevel() const
const std::vector< std::pair< G4XRayModelType, const G4String > > & ActiveVolumes() const
G4int GetCerenkovMaxPhotonsPerStep() const
static G4OpticalParameters * Instance()
G4double GetCerenkovMaxBetaChange() const
G4bool GetCerenkovTrackSecondariesFirst() const
G4bool GetCerenkovStackPhotons() const

◆ ProcessDescription()

void G4GeneralCerenkov::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 300 of file G4GeneralCerenkov.cc.

301{
302 out << "The Cerenkov effect simulates optical photons created by the\n";
303 out << "passage of charged particles through matter. Materials need\n";
304 out << "to have the property RINDEX (refractive index) defined." << G4endl;
305}

Referenced by DumpInfo().

◆ SetMaxBetaChangePerStep()

void G4GeneralCerenkov::SetMaxBetaChangePerStep ( const G4double d)

Definition at line 316 of file G4GeneralCerenkov.cc.

317{
318 fMaxBetaChange = value;
320}
void SetCerenkovMaxBetaChange(G4double)

◆ SetMaxNumPhotonsPerStep()

void G4GeneralCerenkov::SetMaxNumPhotonsPerStep ( const G4int NumPhotons)

Definition at line 323 of file G4GeneralCerenkov.cc.

324{
325 fMaxPhotons = NumPhotons;
327}
void SetCerenkovMaxPhotonsPerStep(G4int)

◆ SetStackPhotons()

void G4GeneralCerenkov::SetStackPhotons ( const G4bool stackingFlag)

Definition at line 330 of file G4GeneralCerenkov.cc.

331{
332 fStackingFlag = stackingFlag;
334}
void SetCerenkovStackPhotons(G4bool)

◆ SetTrackSecondariesFirst()

void G4GeneralCerenkov::SetTrackSecondariesFirst ( const G4bool state)

Definition at line 308 of file G4GeneralCerenkov.cc.

309{
310 fTrackSecondariesFirst = state;
312 fTrackSecondariesFirst);
313}
void SetCerenkovTrackSecondariesFirst(G4bool)

◆ SetVerboseLevel()

void G4GeneralCerenkov::SetVerboseLevel ( G4int verbose)

Definition at line 337 of file G4GeneralCerenkov.cc.


The documentation for this class was generated from the following files: