Geant4 11.1.1
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 : G4VProcess(processName, type)
81 , fNumPhotons(0)
82{
83 secID = G4PhysicsModelCatalog::GetModelID("model_Cerenkov");
85
86 thePhysicsTable = nullptr;
87
88 if(verboseLevel > 0)
89 {
90 G4cout << GetProcessName() << " is created." << G4endl;
91 }
92 Initialise();
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97{
98 if(thePhysicsTable != nullptr)
99 {
101 delete thePhysicsTable;
102 }
103}
104
105void G4Cerenkov::ProcessDescription(std::ostream& out) const
106{
107 out << "The Cerenkov effect simulates optical photons created by the\n";
108 out << "passage of charged particles through matter. Materials need\n";
109 out << "to have the property RINDEX (refractive index) defined.\n";
111
113 out << "Maximum beta change per step: " << params->GetCerenkovMaxBetaChange();
114 out << "Maximum photons per step: " << params->GetCerenkovMaxPhotonsPerStep();
115 out << "Track secondaries first: "
117 out << "Stack photons: " << params->GetCerenkovStackPhotons();
118 out << "Verbose level: " << params->GetCerenkovVerboseLevel();
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123{
124 return (aParticleType.GetPDGCharge() != 0.0 &&
125 aParticleType.GetPDGMass() != 0.0 &&
126 aParticleType.GetParticleName() != "chargedgeantino" &&
127 !aParticleType.IsShortLived())
128 ? true
129 : false;
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134{
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145{
147 return;
148
149 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
150 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
151
152 thePhysicsTable = new G4PhysicsTable(numOfMaterials);
153
154 // loop over materials
155 for(std::size_t i = 0; i < numOfMaterials; ++i)
156 {
157 G4PhysicsFreeVector* cerenkovIntegral = nullptr;
158
159 // Retrieve vector of refraction indices for the material
160 // from the material's optical properties table
161 G4Material* aMaterial = (*theMaterialTable)[i];
163
164 if(MPT)
165 {
166 cerenkovIntegral = new G4PhysicsFreeVector();
167 G4MaterialPropertyVector* refractiveIndex = MPT->GetProperty(kRINDEX);
168
169 if(refractiveIndex)
170 {
171 // Retrieve the first refraction index in vector
172 // of (photon energy, refraction index) pairs
173 G4double currentRI = (*refractiveIndex)[0];
174 if(currentRI > 1.0)
175 {
176 // Create first (photon energy, Cerenkov Integral) pair
177 G4double currentPM = refractiveIndex->Energy(0);
178 G4double currentCAI = 0.0;
179
180 cerenkovIntegral->InsertValues(currentPM, currentCAI);
181
182 // Set previous values to current ones prior to loop
183 G4double prevPM = currentPM;
184 G4double prevCAI = currentCAI;
185 G4double prevRI = currentRI;
186
187 // loop over all (photon energy, refraction index)
188 // pairs stored for this material
189 for(std::size_t ii = 1; ii < refractiveIndex->GetVectorLength(); ++ii)
190 {
191 currentRI = (*refractiveIndex)[ii];
192 currentPM = refractiveIndex->Energy(ii);
193 currentCAI = prevCAI + (currentPM - prevPM) * 0.5 *
194 (1.0 / (prevRI * prevRI) +
195 1.0 / (currentRI * currentRI));
196
197 cerenkovIntegral->InsertValues(currentPM, currentCAI);
198
199 prevPM = currentPM;
200 prevCAI = currentCAI;
201 prevRI = currentRI;
202 }
203 }
204 }
205 }
206
207 // The Cerenkov integral for a given material will be inserted in
208 // thePhysicsTable according to the position of the material in
209 // the material table.
210 thePhysicsTable->insertAt(i, cerenkovIntegral);
211 }
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216 const G4Step& aStep)
217// This routine is called for each tracking Step of a charged particle
218// in a radiator. A Poisson-distributed number of photons is generated
219// according to the Cerenkov formula, distributed evenly along the track
220// segment and uniformly azimuth w.r.t. the particle direction. The
221// parameters are then transformed into the Master Reference System, and
222// they are added to the particle change.
223
224{
226
227 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
228 const G4Material* aMaterial = aTrack.GetMaterial();
229
230 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
231 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
232
233 G4ThreeVector x0 = pPreStepPoint->GetPosition();
234 G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
235 G4double t0 = pPreStepPoint->GetGlobalTime();
236
238 if(!MPT)
239 return pParticleChange;
240
242 if(!Rindex)
243 return pParticleChange;
244
245 G4double charge = aParticle->GetDefinition()->GetPDGCharge();
246 G4double beta = (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta()) * 0.5;
247
248 G4double MeanNumberOfPhotons =
249 GetAverageNumberOfPhotons(charge, beta, aMaterial, Rindex);
250
251 if(MeanNumberOfPhotons <= 0.0)
252 {
253 // return unchanged particle and no secondaries
255 return pParticleChange;
256 }
257
258 G4double step_length = aStep.GetStepLength();
259 MeanNumberOfPhotons = MeanNumberOfPhotons * step_length;
260 fNumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons);
261
262 if(fNumPhotons <= 0 || !fStackingFlag)
263 {
264 // return unchanged particle and no secondaries
266 return pParticleChange;
267 }
268
269 ////////////////////////////////////////////////////////////////
271
272 if(fTrackSecondariesFirst)
273 {
274 if(aTrack.GetTrackStatus() == fAlive)
276 }
277
278 ////////////////////////////////////////////////////////////////
279 G4double Pmin = Rindex->Energy(0);
280 G4double Pmax = Rindex->GetMaxEnergy();
281 G4double dp = Pmax - Pmin;
282
283 G4double nMax = Rindex->GetMaxValue();
284 G4double BetaInverse = 1. / beta;
285
286 G4double maxCos = BetaInverse / nMax;
287 G4double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
288
289 G4double beta1 = pPreStepPoint->GetBeta();
290 G4double beta2 = pPostStepPoint->GetBeta();
291
292 G4double MeanNumberOfPhotons1 =
293 GetAverageNumberOfPhotons(charge, beta1, aMaterial, Rindex);
294 G4double MeanNumberOfPhotons2 =
295 GetAverageNumberOfPhotons(charge, beta2, aMaterial, Rindex);
296
297 for(G4int i = 0; i < fNumPhotons; ++i)
298 {
299 // Determine photon energy
300 G4double rand;
301 G4double sampledEnergy, sampledRI;
302 G4double cosTheta, sin2Theta;
303
304 // sample an energy
305 do
306 {
307 rand = G4UniformRand();
308 sampledEnergy = Pmin + rand * dp;
309 sampledRI = Rindex->Value(sampledEnergy);
310 cosTheta = BetaInverse / sampledRI;
311
312 sin2Theta = (1.0 - cosTheta) * (1.0 + cosTheta);
313 rand = G4UniformRand();
314
315 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
316 } while(rand * maxSin2 > sin2Theta);
317
318 // Create photon momentum direction vector. The momentum direction is still
319 // with respect to the coordinate system where the primary particle
320 // direction is aligned with the z axis
321 rand = G4UniformRand();
322 G4double phi = twopi * rand;
323 G4double sinPhi = std::sin(phi);
324 G4double cosPhi = std::cos(phi);
325 G4double sinTheta = std::sqrt(sin2Theta);
326 G4ParticleMomentum photonMomentum(sinTheta * cosPhi, sinTheta * sinPhi,
327 cosTheta);
328
329 // Rotate momentum direction back to global reference system
330 photonMomentum.rotateUz(p0);
331
332 // Determine polarization of new photon
333 G4ThreeVector photonPolarization(cosTheta * cosPhi, cosTheta * sinPhi,
334 -sinTheta);
335
336 // Rotate back to original coord system
337 photonPolarization.rotateUz(p0);
338
339 // Generate a new photon:
340 auto aCerenkovPhoton =
342
343 aCerenkovPhoton->SetPolarization(photonPolarization);
344 aCerenkovPhoton->SetKineticEnergy(sampledEnergy);
345
346 G4double NumberOfPhotons, N;
347
348 do
349 {
350 rand = G4UniformRand();
351 NumberOfPhotons = MeanNumberOfPhotons1 -
352 rand * (MeanNumberOfPhotons1 - MeanNumberOfPhotons2);
353 N =
354 G4UniformRand() * std::max(MeanNumberOfPhotons1, MeanNumberOfPhotons2);
355 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
356 } while(N > NumberOfPhotons);
357
358 G4double delta = rand * aStep.GetStepLength();
359 G4double deltaTime =
360 delta /
361 (pPreStepPoint->GetVelocity() +
362 rand * (pPostStepPoint->GetVelocity() - pPreStepPoint->GetVelocity()) *
363 0.5);
364
365 G4double aSecondaryTime = t0 + deltaTime;
366 G4ThreeVector aSecondaryPosition = x0 + rand * aStep.GetDeltaPosition();
367
368 // Generate new G4Track object:
369 G4Track* aSecondaryTrack =
370 new G4Track(aCerenkovPhoton, aSecondaryTime, aSecondaryPosition);
371
372 aSecondaryTrack->SetTouchableHandle(
374 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
375 aSecondaryTrack->SetCreatorModelID(secID);
376 aParticleChange.AddSecondary(aSecondaryTrack);
377 }
378
379 if(verboseLevel > 1)
380 {
381 G4cout << "\n Exiting from G4Cerenkov::DoIt -- NumberOfSecondaries = "
383 }
384
385 return pParticleChange;
386}
387
388//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
390{
391 Initialise();
392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397{
398 return 1.;
399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
404{
406 G4double StepLimit = DBL_MAX;
407 fNumPhotons = 0;
408
409 const G4Material* aMaterial = aTrack.GetMaterial();
410 std::size_t materialIndex = aMaterial->GetIndex();
411
412 // If Physics Vector is not defined no Cerenkov photons
413 if(!(*thePhysicsTable)[materialIndex])
414 {
415 return StepLimit;
416 }
417
418 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
419 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
420
421 G4double kineticEnergy = aParticle->GetKineticEnergy();
422 const G4ParticleDefinition* particleType = aParticle->GetDefinition();
423 G4double mass = particleType->GetPDGMass();
424
425 G4double beta = aParticle->GetTotalMomentum() / aParticle->GetTotalEnergy();
426 G4double gamma = aParticle->GetTotalEnergy() / mass;
427
428 G4MaterialPropertiesTable* aMaterialPropertiesTable =
429 aMaterial->GetMaterialPropertiesTable();
430
431 G4MaterialPropertyVector* Rindex = nullptr;
432
433 if(aMaterialPropertiesTable)
434 Rindex = aMaterialPropertiesTable->GetProperty(kRINDEX);
435
436 G4double nMax;
437 if(Rindex)
438 {
439 nMax = Rindex->GetMaxValue();
440 }
441 else
442 {
443 return StepLimit;
444 }
445
446 G4double BetaMin = 1. / nMax;
447 if(BetaMin >= 1.)
448 return StepLimit;
449
450 G4double GammaMin = 1. / std::sqrt(1. - BetaMin * BetaMin);
451 if(gamma < GammaMin)
452 return StepLimit;
453
454 G4double kinEmin = mass * (GammaMin - 1.);
455 G4double RangeMin =
456 G4LossTableManager::Instance()->GetRange(particleType, kinEmin, couple);
458 particleType, kineticEnergy, couple);
459 G4double Step = Range - RangeMin;
460
461 // If the step is smaller than G4ThreeVector::getTolerance(), it may happen
462 // that the particle does not move. See bug 1992.
463 static const G4double minAllowedStep = G4ThreeVector::getTolerance();
464 if(Step < minAllowedStep)
465 return StepLimit;
466
467 if(Step < StepLimit)
468 StepLimit = Step;
469
470 // If user has defined an average maximum number of photons to be generated in
471 // a Step, then calculate the Step length for that number of photons.
472 if(fMaxPhotons > 0)
473 {
474 const G4double charge = aParticle->GetDefinition()->GetPDGCharge();
475 G4double MeanNumberOfPhotons =
476 GetAverageNumberOfPhotons(charge, beta, aMaterial, Rindex);
477 Step = 0.;
478 if(MeanNumberOfPhotons > 0.0)
479 Step = fMaxPhotons / MeanNumberOfPhotons;
480 if(Step > 0. && Step < StepLimit)
481 StepLimit = Step;
482 }
483
484 // If user has defined an maximum allowed change in beta per step
485 if(fMaxBetaChange > 0.)
486 {
488 particleType, kineticEnergy, couple);
489 G4double deltaGamma =
490 gamma - 1. / std::sqrt(1. - beta * beta * (1. - fMaxBetaChange) *
491 (1. - fMaxBetaChange));
492
493 Step = mass * deltaGamma / dedx;
494 if(Step > 0. && Step < StepLimit)
495 StepLimit = Step;
496 }
497
499 return StepLimit;
500}
501
502//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
504 const G4double charge, const G4double beta, const G4Material* aMaterial,
505 G4MaterialPropertyVector* Rindex) const
506// This routine computes the number of Cerenkov photons produced per
507// Geant4-unit (millimeter) in the current medium.
508{
509 constexpr G4double Rfact = 369.81 / (eV * cm);
510 if(beta <= 0.0)
511 return 0.0;
512 G4double BetaInverse = 1. / beta;
513
514 // Vectors used in computation of Cerenkov Angle Integral:
515 // - Refraction Indices for the current material
516 // - new G4PhysicsFreeVector allocated to hold CAI's
517 std::size_t materialIndex = aMaterial->GetIndex();
518
519 // Retrieve the Cerenkov Angle Integrals for this material
520 G4PhysicsVector* CerenkovAngleIntegrals = ((*thePhysicsTable)(materialIndex));
521
522 std::size_t length = CerenkovAngleIntegrals->GetVectorLength();
523 if(0 == length)
524 return 0.0;
525
526 // Min and Max photon energies
527 G4double Pmin = Rindex->Energy(0);
528 G4double Pmax = Rindex->GetMaxEnergy();
529
530 // Min and Max Refraction Indices
531 G4double nMin = Rindex->GetMinValue();
532 G4double nMax = Rindex->GetMaxValue();
533
534 // Max Cerenkov Angle Integral
535 G4double CAImax = (*CerenkovAngleIntegrals)[length - 1];
536
537 G4double dp, ge;
538 // If n(Pmax) < 1/Beta -- no photons generated
539 if(nMax < BetaInverse)
540 {
541 dp = 0.0;
542 ge = 0.0;
543 }
544 // otherwise if n(Pmin) >= 1/Beta -- photons generated
545 else if(nMin > BetaInverse)
546 {
547 dp = Pmax - Pmin;
548 ge = CAImax;
549 }
550 // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then we need to find a P such
551 // that the value of n(P) == 1/Beta. Interpolation is performed by the
552 // GetEnergy() and Value() methods of the G4MaterialPropertiesTable and
553 // the Value() method of G4PhysicsVector.
554 else
555 {
556 Pmin = Rindex->GetEnergy(BetaInverse);
557 dp = Pmax - Pmin;
558
559 G4double CAImin = CerenkovAngleIntegrals->Value(Pmin);
560 ge = CAImax - CAImin;
561
562 if(verboseLevel > 1)
563 {
564 G4cout << "CAImin = " << CAImin << G4endl << "ge = " << ge << G4endl;
565 }
566 }
567
568 // Calculate number of photons
569 G4double NumPhotons = Rfact * charge / eplus * charge / eplus *
570 (dp - ge * BetaInverse * BetaInverse);
571
572 return NumPhotons;
573}
574
575//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
577{
578 fTrackSecondariesFirst = state;
580 fTrackSecondariesFirst);
581}
582
583//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
585{
586 fMaxBetaChange = value * CLHEP::perCent;
588}
589
590//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
592{
593 fMaxPhotons = NumPhotons;
595}
596
597void G4Cerenkov::SetStackPhotons(const G4bool stackingFlag)
598{
599 fStackingFlag = stackingFlag;
601}
602
603//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
605{
606 G4cout << "Dump Physics Table!" << G4endl;
607 for(std::size_t i = 0; i < thePhysicsTable->entries(); ++i)
608 {
609 (*thePhysicsTable)[i]->DumpValues();
610 }
611}
612
613//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
615{
616 verboseLevel = verbose;
618}
@ fCerenkov
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ StronglyForced
@ NotForced
std::vector< G4Material * > G4MaterialTable
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
G4ProcessType
@ 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:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
static double getTolerance()
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
void Initialise()
Definition: G4Cerenkov.cc:133
void ProcessDescription(std::ostream &out) const override
Definition: G4Cerenkov.cc:105
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:584
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:170
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition: G4Cerenkov.cc:402
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Cerenkov.cc:215
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:576
void DumpPhysicsTable() const
Definition: G4Cerenkov.cc:604
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
Definition: G4Cerenkov.cc:503
void SetVerboseLevel(G4int)
Definition: G4Cerenkov.cc:614
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:144
void PreparePhysicsTable(const G4ParticleDefinition &part) override
Definition: G4Cerenkov.cc:389
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition: G4Cerenkov.cc:79
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:122
void SetStackPhotons(const G4bool)
Definition: G4Cerenkov.cc:597
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:395
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:591
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
Definition: G4Material.hh:251
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
size_t GetIndex() const
Definition: G4Material.hh:255
void SetCerenkovMaxBetaChange(G4double)
void SetCerenkovMaxPhotonsPerStep(G4int)
G4int GetCerenkovVerboseLevel() const
G4int GetCerenkovMaxPhotonsPerStep() const
static G4OpticalParameters * Instance()
G4double GetCerenkovMaxBetaChange() const
void SetCerenkovStackPhotons(G4bool)
void SetCerenkovTrackSecondariesFirst(G4bool)
void SetCerenkovVerboseLevel(G4int)
G4bool GetCerenkovTrackSecondariesFirst() const
G4bool GetCerenkovStackPhotons() const
static G4OpticalPhoton * OpticalPhoton()
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void InsertValues(const G4double energy, const G4double value)
static G4int GetModelID(const G4int modelIndex)
void clearAndDestroy()
std::size_t entries() const
void insertAt(std::size_t, G4PhysicsVector *)
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
Definition: G4Step.hh:62
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)
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
virtual void DumpInfo() const
Definition: G4VProcess.cc:173
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define N
Definition: crc32.c:56
#define DBL_MAX
Definition: templates.hh:62