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

#include <G4Cerenkov.hh>

Inheritance diagram for G4Cerenkov:

Public Member Functions

 G4Cerenkov (const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
 ~G4Cerenkov () override
 G4Cerenkov (const G4Cerenkov &right)=delete
G4Cerenkovoperator= (const G4Cerenkov &right)=delete
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType) override
void PreparePhysicsTable (const G4ParticleDefinition &part) override
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double, G4ForceCondition *) override
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
void SetTrackSecondariesFirst (const G4bool state)
G4bool GetTrackSecondariesFirst () const
void SetMaxBetaChangePerStep (const G4double d)
G4double GetMaxBetaChangePerStep () const
void SetMaxNumPhotonsPerStep (const G4int NumPhotons)
G4int GetMaxNumPhotonsPerStep () const
void SetStackPhotons (const G4bool)
G4bool GetStackPhotons () const
G4int GetNumPhotons () const
G4PhysicsTableGetPhysicsTable () const
void DumpPhysicsTable () const
G4double GetAverageNumberOfPhotons (const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
void DumpInfo () const override
void ProcessDescription (std::ostream &out) const override
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 &)

Protected Attributes

G4PhysicsTablethePhysicsTable
std::map< std::size_t, std::size_t > fIndexMPT
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

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 ()

Detailed Description

Definition at line 62 of file G4Cerenkov.hh.

Constructor & Destructor Documentation

◆ G4Cerenkov() [1/2]

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

Definition at line 79 of file G4Cerenkov.cc.

80 : G4VDiscreteProcess(processName, type)
81 , fNumPhotons(0)
82{
83 secID = G4PhysicsModelCatalog::GetModelID("model_Cerenkov");
85
86 thePhysicsTable = nullptr;
87 Initialise();
88
89 if (verboseLevel > 0) {
90 G4cout << GetProcessName() << " is created." << G4endl;
91 }
92}
@ fCerenkov
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * thePhysicsTable
static G4int GetModelID(const G4int modelIndex)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

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

◆ ~G4Cerenkov()

G4Cerenkov::~G4Cerenkov ( )
override

Definition at line 95 of file G4Cerenkov.cc.

96{
97 if(thePhysicsTable != nullptr)
98 {
99 thePhysicsTable->clearAndDestroy();
100 delete thePhysicsTable;
101 }
102}

◆ G4Cerenkov() [2/2]

G4Cerenkov::G4Cerenkov ( const G4Cerenkov & right)
delete

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 141 of file G4Cerenkov.cc.

142{
144 return;
145
146 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
147 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
148
149 // Find the number of materials that have non-empty material property tables
150 std::size_t numOfMaterialsWithMPT = 0;
151 for(std::size_t i = 0; i < numOfMaterials; ++i)
152 {
153 if(((*theMaterialTable)[i])->GetMaterialPropertiesTable())
154 {
155 ++numOfMaterialsWithMPT;
156 }
157 }
158
159 thePhysicsTable = new G4PhysicsTable(numOfMaterialsWithMPT);
160
161 // loop over materials
162 std::size_t indexMPT = 0;
163 for(std::size_t i = 0; i < numOfMaterials; ++i)
164 {
165 G4PhysicsFreeVector* cerenkovIntegral = nullptr;
166
167 // Retrieve vector of refraction indices for the material
168 // from the material's optical properties table
169 G4Material* aMaterial = (*theMaterialTable)[i];
170 G4MaterialPropertiesTable* MPT = aMaterial->GetMaterialPropertiesTable();
171
172 if(MPT)
173 {
174 cerenkovIntegral = new G4PhysicsFreeVector();
175 G4MaterialPropertyVector* refractiveIndex = MPT->GetProperty(kRINDEX);
176
177 if(refractiveIndex)
178 {
179 // Retrieve the first refraction index in vector
180 // of (photon energy, refraction index) pairs
181 G4double currentRI = (*refractiveIndex)[0];
182 if(currentRI > 1.0)
183 {
184 // Create first (photon energy, Cerenkov Integral) pair
185 G4double currentPM = refractiveIndex->Energy(0);
186 G4double currentCAI = 0.0;
187
188 cerenkovIntegral->InsertValues(currentPM, currentCAI);
189
190 // Set previous values to current ones prior to loop
191 G4double prevPM = currentPM;
192 G4double prevCAI = currentCAI;
193 G4double prevRI = currentRI;
194
195 // loop over all (photon energy, refraction index)
196 // pairs stored for this material
197 for(std::size_t ii = 1; ii < refractiveIndex->GetVectorLength(); ++ii)
198 {
199 currentRI = (*refractiveIndex)[ii];
200 currentPM = refractiveIndex->Energy(ii);
201 currentCAI = prevCAI + (currentPM - prevPM) * 0.5 *
202 (1.0 / (prevRI * prevRI) +
203 1.0 / (currentRI * currentRI));
204
205 cerenkovIntegral->InsertValues(currentPM, currentCAI);
206
207 prevPM = currentPM;
208 prevCAI = currentCAI;
209 prevRI = currentRI;
210 }
211 }
212 }
213 // The Cerenkov integral for a given material will be inserted in
214 // thePhysicsTable according to the position of the material in
215 // the material table.
216 thePhysicsTable->insertAt(indexMPT, cerenkovIntegral);
217 fIndexMPT.insert(std::make_pair(i, indexMPT));
218 ++indexMPT;
219 }
220 }
221}
G4PhysicsFreeVector G4MaterialPropertyVector
std::vector< G4Material * > G4MaterialTable
double G4double
Definition G4Types.hh:83
std::map< std::size_t, std::size_t > fIndexMPT
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
void InsertValues(const G4double energy, const G4double value)
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const

◆ DumpInfo()

void G4Cerenkov::DumpInfo ( ) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 651 of file G4Cerenkov.cc.

652{
654}
void ProcessDescription(std::ostream &out) const override

◆ DumpPhysicsTable()

void G4Cerenkov::DumpPhysicsTable ( ) const

Definition at line 634 of file G4Cerenkov.cc.

635{
636 G4cout << "Dump Physics Table!" << G4endl;
637 for(std::size_t i = 0; i < thePhysicsTable->entries(); ++i)
638 {
639 (*thePhysicsTable)[i]->DumpValues();
640 }
641}

◆ GetAverageNumberOfPhotons()

G4double G4Cerenkov::GetAverageNumberOfPhotons ( const G4double charge,
const G4double beta,
const G4Material * aMaterial,
G4MaterialPropertyVector * Rindex ) const

Definition at line 517 of file G4Cerenkov.cc.

522{
523 constexpr G4double Rfact = 369.81 / (eV * cm);
524 if(beta <= 0.0)
525 return 0.0;
526 G4double BetaInverse = 1. / beta;
527
528 // Vectors used in computation of Cerenkov Angle Integral:
529 // - Refraction Indices for the current material
530 // - new G4PhysicsFreeVector allocated to hold CAI's
531 std::size_t materialIndex = aMaterial->GetIndex();
532
533 // Retrieve the Cerenkov Angle Integrals for this material
534 auto it = fIndexMPT.find(materialIndex);
535
536 std::size_t indexMPT = 0;
537 if(it != fIndexMPT.end())
538 {
539 indexMPT = it->second;
540 }
541 else
542 {
544 ed << "G4MaterialPropertiesTable for " << aMaterial->GetName()
545 << " is not found!" << G4endl;
546 G4Exception("G4Cerenkov::GetAverageNumberOfPhotons", "Cerenkov01",
547 FatalException, ed);
548 }
549
550 G4PhysicsVector* CerenkovAngleIntegrals = ((*thePhysicsTable)(indexMPT));
551
552 std::size_t length = CerenkovAngleIntegrals->GetVectorLength();
553 if(0 == length)
554 return 0.0;
555
556 // Min and Max photon energies
557 G4double Pmin = Rindex->Energy(0);
558 G4double Pmax = Rindex->GetMaxEnergy();
559
560 // Min and Max Refraction Indices
561 G4double nMin = Rindex->GetMinValue();
562 G4double nMax = Rindex->GetMaxValue();
563
564 // Max Cerenkov Angle Integral
565 G4double CAImax = (*CerenkovAngleIntegrals)[length - 1];
566
567 G4double dp, ge;
568 // If n(Pmax) < 1/Beta -- no photons generated
569 if(nMax < BetaInverse)
570 {
571 dp = 0.0;
572 ge = 0.0;
573 }
574 // otherwise if n(Pmin) >= 1/Beta -- photons generated
575 else if(nMin > BetaInverse)
576 {
577 dp = Pmax - Pmin;
578 ge = CAImax;
579 }
580 // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then we need to find a P such
581 // that the value of n(P) == 1/Beta. Interpolation is performed by the
582 // GetEnergy() and Value() methods of the G4MaterialPropertiesTable and
583 // the Value() method of G4PhysicsVector.
584 else
585 {
586 Pmin = Rindex->GetEnergy(BetaInverse);
587 dp = Pmax - Pmin;
588
589 G4double CAImin = CerenkovAngleIntegrals->Value(Pmin);
590 ge = CAImax - CAImin;
591
592 if(verboseLevel > 1)
593 {
594 G4cout << "CAImin = " << CAImin << G4endl << "ge = " << ge << G4endl;
595 }
596 }
597
598 // Calculate number of photons
599 G4double NumPhotons = Rfact * charge / eplus * charge / eplus *
600 (dp - ge * BetaInverse * BetaInverse);
601
602 return NumPhotons;
603}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
std::size_t GetIndex() const
const G4String & GetName() const
G4double GetEnergy(const G4double value) const
G4double GetMinValue() const
G4double GetMaxEnergy() const
G4double GetMaxValue() const
G4double Value(const G4double energy, std::size_t &lastidx) const

Referenced by PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ GetMaxBetaChangePerStep()

G4double G4Cerenkov::GetMaxBetaChangePerStep ( ) const
inline

Definition at line 167 of file G4Cerenkov.hh.

168{
169 return fMaxBetaChange;
170}

◆ GetMaxNumPhotonsPerStep()

G4int G4Cerenkov::GetMaxNumPhotonsPerStep ( ) const
inline

Definition at line 172 of file G4Cerenkov.hh.

172{ return fMaxPhotons; }

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 404 of file G4Cerenkov.cc.

406{
407 return DBL_MAX;
408}
#define DBL_MAX
Definition templates.hh:62

◆ GetNumPhotons()

G4int G4Cerenkov::GetNumPhotons ( ) const
inline

Definition at line 176 of file G4Cerenkov.hh.

176{ return fNumPhotons; }

◆ GetPhysicsTable()

G4PhysicsTable * G4Cerenkov::GetPhysicsTable ( ) const
inline

Definition at line 178 of file G4Cerenkov.hh.

179{
180 return thePhysicsTable;
181}

◆ GetStackPhotons()

G4bool G4Cerenkov::GetStackPhotons ( ) const
inline

Definition at line 174 of file G4Cerenkov.hh.

174{ return fStackingFlag; }

◆ GetTrackSecondariesFirst()

G4bool G4Cerenkov::GetTrackSecondariesFirst ( ) const
inline

Definition at line 162 of file G4Cerenkov.hh.

163{
164 return fTrackSecondariesFirst;
165}

◆ IsApplicable()

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

Reimplemented from G4VProcess.

Definition at line 132 of file G4Cerenkov.cc.

133{
134 return ((aParticleType.GetPDGCharge() != 0.0 &&
135 aParticleType.GetPDGMass() != 0.0 &&
136 !aParticleType.IsShortLived()) ||
137 aParticleType.GetParticleName() == "unknown");
138}
const G4String & GetParticleName() const

◆ operator=()

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 224 of file G4Cerenkov.cc.

233{
234 aParticleChange.Initialize(aTrack);
235
236 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
237 const G4Material* aMaterial = aTrack.GetMaterial();
238
239 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
240 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
241
242 G4ThreeVector x0 = pPreStepPoint->GetPosition();
243 G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
244 G4double t0 = pPreStepPoint->GetGlobalTime();
245
246 G4MaterialPropertiesTable* MPT = aMaterial->GetMaterialPropertiesTable();
247 if(!MPT)
248 return pParticleChange;
249
251 if(!Rindex)
252 return pParticleChange;
253
254 G4double charge = aParticle->GetDefinition()->GetPDGCharge();
255
256 G4double beta1 = pPreStepPoint->GetBeta();
257 G4double beta2 = pPostStepPoint->GetBeta();
258 G4double beta = (beta1 + beta2) * 0.5;
259
260 G4double MeanNumberOfPhotons =
261 GetAverageNumberOfPhotons(charge, beta, aMaterial, Rindex);
262 G4double MeanNumberOfPhotons1 =
263 GetAverageNumberOfPhotons(charge, beta1, aMaterial, Rindex);
264 G4double MeanNumberOfPhotons2 =
265 GetAverageNumberOfPhotons(charge, beta2, aMaterial, Rindex);
266
267 if(MeanNumberOfPhotons <= 0.0)
268 {
269 // return unchanged particle and no secondaries
270 aParticleChange.SetNumberOfSecondaries(0);
271 return pParticleChange;
272 }
273
274 MeanNumberOfPhotons *= aStep.GetStepLength();
275 fNumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons);
276
277 // third condition added to prevent infinite loop in do-while below,
278 // see bugzilla 2555
279 if(fNumPhotons <= 0 || !fStackingFlag ||
280 std::max(MeanNumberOfPhotons1, MeanNumberOfPhotons2) < 1e-15)
281 {
282 // return unchanged particle and no secondaries
283 aParticleChange.SetNumberOfSecondaries(0);
284 return pParticleChange;
285 }
286
287 G4double deltaVelocity = pPostStepPoint->GetVelocity() -
288 pPreStepPoint->GetVelocity();
289 auto touchableHandle = aStep.GetPreStepPoint()->GetTouchableHandle();
290
291 ////////////////////////////////////////////////////////////////
292 aParticleChange.SetNumberOfSecondaries(fNumPhotons);
293
294 if(fTrackSecondariesFirst)
295 {
296 if(aTrack.GetTrackStatus() == fAlive)
297 aParticleChange.ProposeTrackStatus(fSuspend);
298 }
299
300 ////////////////////////////////////////////////////////////////
301 G4double Pmin = Rindex->Energy(0);
302 G4double Pmax = Rindex->GetMaxEnergy();
303 G4double dp = Pmax - Pmin;
304 G4double deltaNumberOfPhotons = MeanNumberOfPhotons1 - MeanNumberOfPhotons2;
305 G4double maxNumberOfPhotons =
306 std::max(MeanNumberOfPhotons1, MeanNumberOfPhotons2);
307
308 G4double nMax = Rindex->GetMaxValue();
309 G4double BetaInverse = 1. / beta;
310
311 G4double maxCos = BetaInverse / nMax;
312 G4double maxSin2 = 1.0 - maxCos * maxCos;
313
314 for(G4int i = 0; i < fNumPhotons; ++i)
315 {
316 // Determine photon energy
317 G4double rand;
318 G4double sampledEnergy;
319 G4double cosTheta, sin2Theta;
320
321 // sample an energy
322 do
323 {
324 rand = G4UniformRand();
325 sampledEnergy = Pmin + rand * dp;
326 cosTheta = BetaInverse / Rindex->Value(sampledEnergy);
327
328 sin2Theta = 1.0 - cosTheta * cosTheta;
329 rand = G4UniformRand();
330
331 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
332 } while(rand * maxSin2 > sin2Theta);
333
334 // Create photon momentum direction vector. The momentum direction is still
335 // with respect to the coordinate system where the primary particle
336 // direction is aligned with the z axis
337 rand = G4UniformRand();
338 G4double phi = twopi * rand;
339 G4double sinPhi = std::sin(phi);
340 G4double cosPhi = std::cos(phi);
341 G4double sinTheta = std::sqrt(sin2Theta);
342 G4ParticleMomentum photonMomentum(sinTheta * cosPhi, sinTheta * sinPhi,
343 cosTheta);
344
345 // Rotate momentum direction back to global reference system
346 photonMomentum.rotateUz(p0);
347
348 // Determine polarization of new photon
349 G4ThreeVector photonPolarization(cosTheta * cosPhi, cosTheta * sinPhi,
350 -sinTheta);
351
352 // Rotate back to original coord system
353 photonPolarization.rotateUz(p0);
354
355 // Generate a new photon:
356 auto aCerenkovPhoton =
357 new G4DynamicParticle(G4OpticalPhoton::OpticalPhoton(), photonMomentum);
358
359 aCerenkovPhoton->SetPolarization(photonPolarization);
360 aCerenkovPhoton->SetKineticEnergy(sampledEnergy);
361
362 G4double NumberOfPhotons;
363
364 do
365 {
366 rand = G4UniformRand();
367 NumberOfPhotons = MeanNumberOfPhotons1 - rand * deltaNumberOfPhotons;
368 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
369 } while(G4UniformRand() * maxNumberOfPhotons > NumberOfPhotons);
370
371 G4double delta = rand * aStep.GetStepLength();
372 G4double deltaTime =
373 delta / (pPreStepPoint->GetVelocity() + rand * deltaVelocity * 0.5);
374
375 G4double aSecondaryTime = t0 + deltaTime;
376 G4ThreeVector aSecondaryPosition = x0 + rand * aStep.GetDeltaPosition();
377
378 // Generate new G4Track object:
379 G4Track* aSecondaryTrack =
380 new G4Track(aCerenkovPhoton, aSecondaryTime, aSecondaryPosition);
381
382 aSecondaryTrack->SetTouchableHandle(touchableHandle);
383 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
384 aSecondaryTrack->SetCreatorModelID(secID);
385 aParticleChange.AddSecondary(aSecondaryTrack);
386 }
387
388 if(verboseLevel > 1)
389 {
390 G4cout << "\n Exiting from G4Cerenkov::DoIt -- NumberOfSecondaries = "
391 << aParticleChange.GetNumberOfSecondaries() << G4endl;
392 }
393
394 return pParticleChange;
395}
G4ThreeVector G4ParticleMomentum
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
CLHEP::Hep3Vector G4ThreeVector
@ fSuspend
@ fAlive
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
G4ParticleDefinition * GetDefinition() const
static G4OpticalPhoton * OpticalPhoton()
G4double GetVelocity() const
G4double GetBeta() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4ThreeVector GetDeltaPosition() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
void SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
G4ParticleChange aParticleChange
G4VParticleChange * pParticleChange

◆ PostStepGetPhysicalInteractionLength()

G4double G4Cerenkov::PostStepGetPhysicalInteractionLength ( const G4Track & aTrack,
G4double ,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 411 of file G4Cerenkov.cc.

413{
415 G4double StepLimit = DBL_MAX;
416 if (aTrack.GetDynamicParticle()->GetCharge() == 0.0) {
417 return StepLimit;
418 }
419 fNumPhotons = 0;
420
421 const G4Material* aMaterial = aTrack.GetMaterial();
422 std::size_t materialIndex = aMaterial->GetIndex();
423
424 // If Physics Vector is not defined no Cerenkov photons
425 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
426 auto const MPT =
427 ((*materialTable)[materialIndex])->GetMaterialPropertiesTable();
428 if (nullptr == MPT) {
429 return StepLimit;
430 }
431
432 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
433 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
434
435 G4double kineticEnergy = aParticle->GetKineticEnergy();
436 const G4ParticleDefinition* particleType = aParticle->GetDefinition();
437 G4double mass = particleType->GetPDGMass();
438
439 G4double beta = aParticle->GetTotalMomentum() / aParticle->GetTotalEnergy();
440 G4double gamma = aParticle->GetTotalEnergy() / mass;
441
442 G4MaterialPropertiesTable* aMaterialPropertiesTable =
443 aMaterial->GetMaterialPropertiesTable();
444
445 G4MaterialPropertyVector* Rindex = nullptr;
446
447 if(aMaterialPropertiesTable)
448 Rindex = aMaterialPropertiesTable->GetProperty(kRINDEX);
449
450 G4double nMax;
451 if(Rindex)
452 {
453 nMax = Rindex->GetMaxValue();
454 }
455 else
456 {
457 return StepLimit;
458 }
459
460 G4double BetaMin = 1. / nMax;
461 if(BetaMin >= 1.)
462 return StepLimit;
463
464 G4double GammaMin = 1. / std::sqrt(1. - BetaMin * BetaMin);
465 if(gamma < GammaMin)
466 return StepLimit;
467
468 G4double kinEmin = mass * (GammaMin - 1.);
469 G4double RangeMin =
470 G4LossTableManager::Instance()->GetRange(particleType, kinEmin, couple);
472 particleType, kineticEnergy, couple);
473 G4double Step = Range - RangeMin;
474
475 // If the step is smaller than G4ThreeVector::getTolerance(), it may happen
476 // that the particle does not move. See bug 1992.
477 static const G4double minAllowedStep = G4ThreeVector::getTolerance();
478 if(Step < minAllowedStep)
479 return StepLimit;
480
481 if(Step < StepLimit)
482 StepLimit = Step;
483
484 // If user has defined an average maximum number of photons to be generated in
485 // a Step, then calculate the Step length for that number of photons.
486 if(fMaxPhotons > 0)
487 {
488 const G4double charge = aParticle->GetDefinition()->GetPDGCharge();
489 G4double MeanNumberOfPhotons =
490 GetAverageNumberOfPhotons(charge, beta, aMaterial, Rindex);
491 Step = 0.;
492 if(MeanNumberOfPhotons > 0.0)
493 Step = fMaxPhotons / MeanNumberOfPhotons;
494 if(Step > 0. && Step < StepLimit)
495 StepLimit = Step;
496 }
497
498 // If user has defined an maximum allowed change in beta per step
499 if(fMaxBetaChange > 0.)
500 {
502 particleType, kineticEnergy, couple);
503 G4double deltaGamma =
504 gamma - 1. / std::sqrt(1. - beta * beta * (1. - fMaxBetaChange) *
505 (1. - fMaxBetaChange));
506
507 Step = mass * deltaGamma / dedx;
508 if(Step > 0. && Step < StepLimit)
509 StepLimit = Step;
510 }
511
513 return StepLimit;
514}
G4double condition(const G4ErrorSymMatrix &m)
@ StronglyForced
@ NotForced
static double getTolerance()
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
static G4LossTableManager * Instance()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

◆ PreparePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 398 of file G4Cerenkov.cc.

399{
400 Initialise();
401}

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 104 of file G4Cerenkov.cc.

105{
106 out << "The Cerenkov effect simulates optical photons created by the\n";
107 out << "passage of charged particles through matter. Materials need\n";
108 out << "to have the property RINDEX (refractive index) defined.\n";
110
111 G4OpticalParameters* params = G4OpticalParameters::Instance();
112 out << "Maximum beta change per step: " << params->GetCerenkovMaxBetaChange();
113 out << "Maximum photons per step: " << params->GetCerenkovMaxPhotonsPerStep();
114 out << "Track secondaries first: "
116 out << "Stack photons: " << params->GetCerenkovStackPhotons();
117 out << "Verbose level: " << params->GetCerenkovVerboseLevel();
118}
G4int GetCerenkovVerboseLevel() const
G4int GetCerenkovMaxPhotonsPerStep() const
static G4OpticalParameters * Instance()
G4double GetCerenkovMaxBetaChange() const
G4bool GetCerenkovTrackSecondariesFirst() const
G4bool GetCerenkovStackPhotons() const
virtual void DumpInfo() const

Referenced by DumpInfo().

◆ SetMaxBetaChangePerStep()

void G4Cerenkov::SetMaxBetaChangePerStep ( const G4double d)

Definition at line 614 of file G4Cerenkov.cc.

615{
616 fMaxBetaChange = value * CLHEP::perCent;
618}
void SetCerenkovMaxBetaChange(G4double)

◆ SetMaxNumPhotonsPerStep()

void G4Cerenkov::SetMaxNumPhotonsPerStep ( const G4int NumPhotons)

Definition at line 621 of file G4Cerenkov.cc.

622{
623 fMaxPhotons = NumPhotons;
625}
void SetCerenkovMaxPhotonsPerStep(G4int)

◆ SetStackPhotons()

void G4Cerenkov::SetStackPhotons ( const G4bool stackingFlag)

Definition at line 627 of file G4Cerenkov.cc.

628{
629 fStackingFlag = stackingFlag;
631}
void SetCerenkovStackPhotons(G4bool)

◆ SetTrackSecondariesFirst()

void G4Cerenkov::SetTrackSecondariesFirst ( const G4bool state)

Definition at line 606 of file G4Cerenkov.cc.

607{
608 fTrackSecondariesFirst = state;
610 fTrackSecondariesFirst);
611}
void SetCerenkovTrackSecondariesFirst(G4bool)

◆ SetVerboseLevel()

void G4Cerenkov::SetVerboseLevel ( G4int verbose)

Definition at line 644 of file G4Cerenkov.cc.

Member Data Documentation

◆ fIndexMPT

std::map<std::size_t, std::size_t> G4Cerenkov::fIndexMPT
protected

Definition at line 146 of file G4Cerenkov.hh.

Referenced by BuildPhysicsTable(), and GetAverageNumberOfPhotons().

◆ thePhysicsTable

G4PhysicsTable* G4Cerenkov::thePhysicsTable
protected

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