Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Cerenkov.cc
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 Implementation
28////////////////////////////////////////////////////////////////////////
29//
30// File: G4Cerenkov.cc
31// Description: Discrete Process -- Generation of Cerenkov Photons
32// Version: 2.1
33// Created: 1996-02-21
34// Author: Juliet Armstrong
35// Updated: 2007-09-30 by Peter Gumplinger
36// > change inheritance to G4VDiscreteProcess
37// GetContinuousStepLimit -> GetMeanFreePath (StronglyForced)
38// AlongStepDoIt -> PostStepDoIt
39// 2005-08-17 by Peter Gumplinger
40// > change variable name MeanNumPhotons -> MeanNumberOfPhotons
41// 2005-07-28 by Peter Gumplinger
42// > add G4ProcessType to constructor
43// 2001-09-17, migration of Materials to pure STL (mma)
44// 2000-11-12 by Peter Gumplinger
45// > add check on CerenkovAngleIntegrals->IsFilledVectorExist()
46// in method GetAverageNumberOfPhotons
47// > and a test for MeanNumberOfPhotons <= 0.0 in DoIt
48// 2000-09-18 by Peter Gumplinger
49// > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition();
50// aSecondaryTrack->SetTouchable(0);
51// 1999-10-29 by Peter Gumplinger
52// > change: == into <= in GetContinuousStepLimit
53// 1997-08-08 by Peter Gumplinger
54// > add protection against /0
55// > G4MaterialPropertiesTable; new physics/tracking scheme
56//
57////////////////////////////////////////////////////////////////////////
58
59#include "G4Cerenkov.hh"
60
61#include "G4ios.hh"
62#include "G4LossTableManager.hh"
63#include "G4Material.hh"
67#include "G4OpticalPhoton.hh"
69#include "G4ParticleMomentum.hh"
72#include "G4Poisson.hh"
73#include "G4SystemOfUnits.hh"
74#include "G4ThreeVector.hh"
75#include "Randomize.hh"
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96{
97 if(thePhysicsTable != nullptr)
98 {
99 thePhysicsTable->clearAndDestroy();
100 delete thePhysicsTable;
101 }
102}
103
104void G4Cerenkov::ProcessDescription(std::ostream& out) const
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
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}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121void G4Cerenkov::Initialise()
122{
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133{
134 return ((aParticleType.GetPDGCharge() != 0.0 &&
135 aParticleType.GetPDGMass() != 0.0 &&
136 !aParticleType.IsShortLived()) ||
137 aParticleType.GetParticleName() == "unknown");
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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];
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}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225 const G4Step& aStep)
226// This routine is called for each tracking Step of a charged particle
227// in a radiator. A Poisson-distributed number of photons is generated
228// according to the Cerenkov formula, distributed evenly along the track
229// segment and uniformly azimuth w.r.t. the particle direction. The
230// parameters are then transformed into the Master Reference System, and
231// they are added to the particle change.
232
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
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 =
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}
396
397//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
399{
400 Initialise();
401}
402
403//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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}
515
516//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
518 const G4double charge, const G4double beta, const G4Material* aMaterial,
519 G4MaterialPropertyVector* Rindex) const
520// This routine computes the number of Cerenkov photons produced per
521// Geant4-unit (millimeter) in the current medium.
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}
604
605//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
607{
608 fTrackSecondariesFirst = state;
610 fTrackSecondariesFirst);
611}
612
613//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
615{
616 fMaxBetaChange = value * CLHEP::perCent;
618}
619
620//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
622{
623 fMaxPhotons = NumPhotons;
625}
626
627void G4Cerenkov::SetStackPhotons(const G4bool stackingFlag)
628{
629 fStackingFlag = stackingFlag;
631}
632
633//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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}
642
643//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
649
650//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
652{
654}
655
656//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fCerenkov
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ StronglyForced
@ NotForced
G4PhysicsFreeVector G4MaterialPropertyVector
std::vector< G4Material * > G4MaterialTable
G4ThreeVector G4ParticleMomentum
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
G4ProcessType
CLHEP::Hep3Vector G4ThreeVector
@ fSuspend
@ fAlive
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
static double getTolerance()
Hep3Vector & rotateUz(const Hep3Vector &)
void ProcessDescription(std::ostream &out) const override
void SetMaxBetaChangePerStep(const G4double d)
G4PhysicsTable * thePhysicsTable
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
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
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition G4Cerenkov.cc:79
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
void SetStackPhotons(const G4bool)
void DumpInfo() const override
std::map< std::size_t, std::size_t > fIndexMPT
~G4Cerenkov() override
Definition G4Cerenkov.cc:95
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
G4double GetCharge() const
G4ParticleDefinition * GetDefinition() 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)
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static std::size_t GetNumberOfMaterials()
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
void SetCerenkovMaxBetaChange(G4double)
void SetCerenkovMaxPhotonsPerStep(G4int)
G4int GetCerenkovVerboseLevel() const
G4int GetCerenkovMaxPhotonsPerStep() const
static G4OpticalParameters * Instance()
G4double GetCerenkovMaxBetaChange() const
void SetCerenkovStackPhotons(G4bool)
void SetCerenkovTrackSecondariesFirst(G4bool)
G4bool GetCerenkovTrackSecondariesFirst() const
G4bool GetCerenkovStackPhotons() const
static G4OpticalPhoton * OpticalPhoton()
const G4String & GetParticleName() const
void InsertValues(const G4double energy, const G4double value)
static G4int GetModelID(const G4int modelIndex)
G4double GetMaxEnergy() const
G4double GetMaxValue() const
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
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)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetParentID(const G4int aValue)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4ParticleChange aParticleChange
G4int verboseLevel
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
virtual void DumpInfo() const
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62