Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // $Id: G4Cerenkov.cc 98002 2016-06-30 13:03:36Z gcosmo $
28 //
30 // Cerenkov Radiation Class Implementation
32 //
33 // File: G4Cerenkov.cc
34 // Description: Discrete Process -- Generation of Cerenkov Photons
35 // Version: 2.1
36 // Created: 1996-02-21
37 // Author: Juliet Armstrong
38 // Updated: 2007-09-30 by Peter Gumplinger
39 // > change inheritance to G4VDiscreteProcess
40 // GetContinuousStepLimit -> GetMeanFreePath (StronglyForced)
41 // AlongStepDoIt -> PostStepDoIt
42 // 2005-08-17 by Peter Gumplinger
43 // > change variable name MeanNumPhotons -> MeanNumberOfPhotons
44 // 2005-07-28 by Peter Gumplinger
45 // > add G4ProcessType to constructor
46 // 2001-09-17, migration of Materials to pure STL (mma)
47 // 2000-11-12 by Peter Gumplinger
48 // > add check on CerenkovAngleIntegrals->IsFilledVectorExist()
49 // in method GetAverageNumberOfPhotons
50 // > and a test for MeanNumberOfPhotons <= 0.0 in DoIt
51 // 2000-09-18 by Peter Gumplinger
52 // > change: aSecondaryPosition=x0+rand*aStep.GetDeltaPosition();
53 // aSecondaryTrack->SetTouchable(0);
54 // 1999-10-29 by Peter Gumplinger
55 // > change: == into <= in GetContinuousStepLimit
56 // 1997-08-08 by Peter Gumplinger
57 // > add protection against /0
58 // > G4MaterialPropertiesTable; new physics/tracking scheme
59 //
60 // mail: gum@triumf.ca
61 //
63 
64 #include "G4ios.hh"
65 #include "G4PhysicalConstants.hh"
66 #include "G4SystemOfUnits.hh"
67 #include "G4Poisson.hh"
68 #include "G4EmProcessSubType.hh"
69 
70 #include "G4LossTableManager.hh"
71 #include "G4MaterialCutsCouple.hh"
72 #include "G4ParticleDefinition.hh"
73 
74 #include "G4Cerenkov.hh"
75 
77 // Class Implementation
79 
80 //G4bool G4Cerenkov::fTrackSecondariesFirst = false;
81 //G4double G4Cerenkov::fMaxBetaChange = 0.;
82 //G4int G4Cerenkov::fMaxPhotons = 0;
83 
85  // Operators
87 
88 // G4Cerenkov::operator=(const G4Cerenkov &right)
89 // {
90 // }
91 
93  // Constructors
95 
97  : G4VProcess(processName, type),
98  fTrackSecondariesFirst(false),
99  fMaxBetaChange(0.0),
100  fMaxPhotons(0),
101  fStackingFlag(true),
102  fNumPhotons(0)
103 {
105 
106  thePhysicsTable = nullptr;
107 
108  if (verboseLevel>0) {
109  G4cout << GetProcessName() << " is created " << G4endl;
110  }
111 }
112 
113 // G4Cerenkov::G4Cerenkov(const G4Cerenkov &right)
114 // {
115 // }
116 
118  // Destructors
120 
122 {
123  if (thePhysicsTable != nullptr) {
125  delete thePhysicsTable;
126  }
127 }
128 
130  // Methods
132 
134 {
135  G4bool result = false;
136 
137  if (aParticleType.GetPDGCharge() != 0.0 &&
138  aParticleType.GetPDGMass() != 0.0 &&
139  aParticleType.GetParticleName() != "chargedgeantino" &&
140  !aParticleType.IsShortLived() ) { result = true; }
141 
142  return result;
143 }
144 
146 {
147  fTrackSecondariesFirst = state;
148 }
149 
151 {
152  fMaxBetaChange = value*CLHEP::perCent;
153 }
154 
156 {
157  fMaxPhotons = NumPhotons;
158 }
159 
161 {
162  if (!thePhysicsTable) BuildThePhysicsTable();
163 }
164 
165 // PostStepDoIt
166 // -------------
167 //
169 G4Cerenkov::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
170 
171 // This routine is called for each tracking Step of a charged particle
172 // in a radiator. A Poisson-distributed number of photons is generated
173 // according to the Cerenkov formula, distributed evenly along the track
174 // segment and uniformly azimuth w.r.t. the particle direction. The
175 // parameters are then transformed into the Master Reference System, and
176 // they are added to the particle change.
177 
178 {
180  // Should we ensure that the material is dispersive?
182 
183  aParticleChange.Initialize(aTrack);
184 
185  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
186  const G4Material* aMaterial = aTrack.GetMaterial();
187 
188  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
189  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
190 
191  G4ThreeVector x0 = pPreStepPoint->GetPosition();
192  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
193  G4double t0 = pPreStepPoint->GetGlobalTime();
194 
195  G4MaterialPropertiesTable* aMaterialPropertiesTable =
196  aMaterial->GetMaterialPropertiesTable();
197  if (!aMaterialPropertiesTable) return pParticleChange;
198 
199  G4MaterialPropertyVector* Rindex =
200  aMaterialPropertiesTable->GetProperty("RINDEX");
201  if (!Rindex) return pParticleChange;
202 
203  // particle charge
204  G4double charge = aParticle->GetDefinition()->GetPDGCharge();
205 
206  // particle beta
207  G4double beta = (pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta())*0.5;
208 
209  fNumPhotons = 0;
210 
211  G4double MeanNumberOfPhotons =
212  GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex);
213 
214  if (MeanNumberOfPhotons <= 0.0) {
215 
216  // return unchanged particle and no secondaries
217 
219 
220  return pParticleChange;
221 
222  }
223 
224  G4double step_length = aStep.GetStepLength();
225 
226  MeanNumberOfPhotons = MeanNumberOfPhotons * step_length;
227 
228  fNumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons);
229 
230  if ( fNumPhotons <= 0 || !fStackingFlag ) {
231 
232  // return unchanged particle and no secondaries
233 
235 
236  return pParticleChange;
237 
238  }
239 
241 
243 
244  if (fTrackSecondariesFirst) {
245  if (aTrack.GetTrackStatus() == fAlive )
247  }
248 
250 
251  G4double Pmin = Rindex->GetMinLowEdgeEnergy();
252  G4double Pmax = Rindex->GetMaxLowEdgeEnergy();
253  G4double dp = Pmax - Pmin;
254 
255  G4double nMax = Rindex->GetMaxValue();
256 
257  G4double BetaInverse = 1./beta;
258 
259  G4double maxCos = BetaInverse / nMax;
260  G4double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
261 
262  G4double beta1 = pPreStepPoint ->GetBeta();
263  G4double beta2 = pPostStepPoint->GetBeta();
264 
265  G4double MeanNumberOfPhotons1 =
266  GetAverageNumberOfPhotons(charge,beta1,aMaterial,Rindex);
267  G4double MeanNumberOfPhotons2 =
268  GetAverageNumberOfPhotons(charge,beta2,aMaterial,Rindex);
269 
270  for (G4int i = 0; i < fNumPhotons; i++) {
271 
272  // Determine photon energy
273 
274  G4double rand;
275  G4double sampledEnergy, sampledRI;
276  G4double cosTheta, sin2Theta;
277 
278  // sample an energy
279 
280  do {
281  rand = G4UniformRand();
282  sampledEnergy = Pmin + rand * dp;
283  sampledRI = Rindex->Value(sampledEnergy);
284  cosTheta = BetaInverse / sampledRI;
285 
286  sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
287  rand = G4UniformRand();
288 
289  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
290  } while (rand*maxSin2 > sin2Theta);
291 
292  // Generate random position of photon on cone surface
293  // defined by Theta
294 
295  rand = G4UniformRand();
296 
297  G4double phi = twopi*rand;
298  G4double sinPhi = std::sin(phi);
299  G4double cosPhi = std::cos(phi);
300 
301  // calculate x,y, and z components of photon energy
302  // (in coord system with primary particle direction
303  // aligned with the z axis)
304 
305  G4double sinTheta = std::sqrt(sin2Theta);
306  G4double px = sinTheta*cosPhi;
307  G4double py = sinTheta*sinPhi;
308  G4double pz = cosTheta;
309 
310  // Create photon momentum direction vector
311  // The momentum direction is still with respect
312  // to the coordinate system where the primary
313  // particle direction is aligned with the z axis
314 
315  G4ParticleMomentum photonMomentum(px, py, pz);
316 
317  // Rotate momentum direction back to global reference
318  // system
319 
320  photonMomentum.rotateUz(p0);
321 
322  // Determine polarization of new photon
323 
324  G4double sx = cosTheta*cosPhi;
325  G4double sy = cosTheta*sinPhi;
326  G4double sz = -sinTheta;
327 
328  G4ThreeVector photonPolarization(sx, sy, sz);
329 
330  // Rotate back to original coord system
331 
332  photonPolarization.rotateUz(p0);
333 
334  // Generate a new photon:
335 
336  G4DynamicParticle* aCerenkovPhoton =
338 
339  aCerenkovPhoton->SetPolarization(photonPolarization.x(),
340  photonPolarization.y(),
341  photonPolarization.z());
342 
343  aCerenkovPhoton->SetKineticEnergy(sampledEnergy);
344 
345  // Generate new G4Track object:
346 
347  G4double NumberOfPhotons, N;
348 
349  do {
350  rand = G4UniformRand();
351  NumberOfPhotons = MeanNumberOfPhotons1 - rand *
352  (MeanNumberOfPhotons1-MeanNumberOfPhotons2);
353  N = G4UniformRand() *
354  std::max(MeanNumberOfPhotons1,MeanNumberOfPhotons2);
355  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
356  } while (N > NumberOfPhotons);
357 
358  G4double delta = rand * aStep.GetStepLength();
359 
360  G4double deltaTime = delta / (pPreStepPoint->GetVelocity()+
361  rand*(pPostStepPoint->GetVelocity()-
362  pPreStepPoint->GetVelocity())*0.5);
363 
364  G4double aSecondaryTime = t0 + deltaTime;
365 
366  G4ThreeVector aSecondaryPosition = x0 + rand * aStep.GetDeltaPosition();
367 
368  G4Track* aSecondaryTrack =
369  new G4Track(aCerenkovPhoton,aSecondaryTime,aSecondaryPosition);
370 
371  aSecondaryTrack->SetTouchableHandle(
373 
374  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
375 
376  aParticleChange.AddSecondary(aSecondaryTrack);
377  }
378 
379  if (verboseLevel>0) {
380  G4cout <<"\n Exiting from G4Cerenkov::DoIt -- NumberOfSecondaries = "
382  }
383 
384  return pParticleChange;
385 }
386 
387 // BuildThePhysicsTable for the Cerenkov process
388 // ---------------------------------------------
389 //
390 
391 void G4Cerenkov::BuildThePhysicsTable()
392 {
393  if (thePhysicsTable) return;
394 
395  const G4MaterialTable* theMaterialTable=
397  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
398 
399  // create new physics table
400 
401  thePhysicsTable = new G4PhysicsTable(numOfMaterials);
402 
403  // loop for materials
404 
405  for (G4int i=0 ; i < numOfMaterials; i++) {
406 
407  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector = 0;
408 
409  // Retrieve vector of refraction indices for the material
410  // from the material's optical properties table
411 
412  G4Material* aMaterial = (*theMaterialTable)[i];
413 
414  G4MaterialPropertiesTable* aMaterialPropertiesTable =
415  aMaterial->GetMaterialPropertiesTable();
416 
417  if (aMaterialPropertiesTable) {
418  aPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
419  G4MaterialPropertyVector* theRefractionIndexVector =
420  aMaterialPropertiesTable->GetProperty("RINDEX");
421 
422  if (theRefractionIndexVector) {
423 
424  // Retrieve the first refraction index in vector
425  // of (photon energy, refraction index) pairs
426 
427  G4double currentRI = (*theRefractionIndexVector)[0];
428 
429  if (currentRI > 1.0) {
430 
431  // Create first (photon energy, Cerenkov Integral)
432  // pair
433 
434  G4double currentPM = theRefractionIndexVector->Energy(0);
435  G4double currentCAI = 0.0;
436 
437  aPhysicsOrderedFreeVector->InsertValues(currentPM , currentCAI);
438 
439  // Set previous values to current ones prior to loop
440 
441  G4double prevPM = currentPM;
442  G4double prevCAI = currentCAI;
443  G4double prevRI = currentRI;
444 
445  // loop over all (photon energy, refraction index)
446  // pairs stored for this material
447 
448  for (size_t ii = 1;
449  ii < theRefractionIndexVector->GetVectorLength();
450  ++ii) {
451  currentRI = (*theRefractionIndexVector)[ii];
452  currentPM = theRefractionIndexVector->Energy(ii);
453 
454  currentCAI = 0.5*(1.0/(prevRI*prevRI) +
455  1.0/(currentRI*currentRI));
456 
457  currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
458 
459  aPhysicsOrderedFreeVector->
460  InsertValues(currentPM, currentCAI);
461 
462  prevPM = currentPM;
463  prevCAI = currentCAI;
464  prevRI = currentRI;
465  }
466 
467  }
468  }
469  }
470 
471  // The Cerenkov integral for a given material
472  // will be inserted in thePhysicsTable
473  // according to the position of the material in
474  // the material table.
475 
476  thePhysicsTable->insertAt(i,aPhysicsOrderedFreeVector);
477 
478  }
479 }
480 
481 // GetMeanFreePath
482 // ---------------
483 //
484 
486  G4double,
488 {
489  return 1.;
490 }
491 
493  const G4Track& aTrack,
494  G4double,
496 {
497  *condition = NotForced;
498  G4double StepLimit = DBL_MAX;
499 
500  const G4Material* aMaterial = aTrack.GetMaterial();
501  G4int materialIndex = aMaterial->GetIndex();
502 
503  // If Physics Vector is not defined no Cerenkov photons
504  // this check avoid string comparison below
505 
506  if(!(*thePhysicsTable)[materialIndex]) { return StepLimit; }
507 
508  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
509  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
510 
511  G4double kineticEnergy = aParticle->GetKineticEnergy();
512  const G4ParticleDefinition* particleType = aParticle->GetDefinition();
513  G4double mass = particleType->GetPDGMass();
514 
515  // particle beta
516  G4double beta = aParticle->GetTotalMomentum() /
517  aParticle->GetTotalEnergy();
518  // particle gamma
519  G4double gamma = aParticle->GetTotalEnergy()/mass;
520 
521  G4MaterialPropertiesTable* aMaterialPropertiesTable =
522  aMaterial->GetMaterialPropertiesTable();
523 
524  G4MaterialPropertyVector* Rindex = NULL;
525 
526  if (aMaterialPropertiesTable)
527  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
528 
529  G4double nMax;
530  if (Rindex) {
531  nMax = Rindex->GetMaxValue();
532  } else {
533  return StepLimit;
534  }
535 
536  G4double BetaMin = 1./nMax;
537  if ( BetaMin >= 1. ) return StepLimit;
538 
539  G4double GammaMin = 1./std::sqrt(1.-BetaMin*BetaMin);
540 
541  if (gamma < GammaMin ) return StepLimit;
542 
543  G4double kinEmin = mass*(GammaMin-1.);
544 
545  G4double RangeMin = G4LossTableManager::Instance()->GetRange(particleType,
546  kinEmin,
547  couple);
548  G4double Range = G4LossTableManager::Instance()->GetRange(particleType,
549  kineticEnergy,
550  couple);
551 
552  G4double Step = Range - RangeMin;
553 // if (Step < 1.*um ) return StepLimit;
554 
555  if (Step > 0. && Step < StepLimit) StepLimit = Step;
556 
557  // If user has defined an average maximum number of photons to
558  // be generated in a Step, then calculate the Step length for
559  // that number of photons.
560 
561  if (fMaxPhotons > 0) {
562 
563  // particle charge
564  const G4double charge = aParticle->GetDefinition()->GetPDGCharge();
565 
566  G4double MeanNumberOfPhotons =
567  GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex);
568 
569  Step = 0.;
570  if (MeanNumberOfPhotons > 0.0) Step = fMaxPhotons / MeanNumberOfPhotons;
571 
572  if (Step > 0. && Step < StepLimit) StepLimit = Step;
573  }
574 
575  // If user has defined an maximum allowed change in beta per step
576  if (fMaxBetaChange > 0.) {
577 
578  G4double dedx = G4LossTableManager::Instance()->GetDEDX(particleType,
579  kineticEnergy,
580  couple);
581 
582  G4double deltaGamma = gamma - 1./std::sqrt(1.-beta*beta*
583  (1.-fMaxBetaChange)*
584  (1.-fMaxBetaChange));
585 
586  Step = mass * deltaGamma / dedx;
587 
588  if (Step > 0. && Step < StepLimit) StepLimit = Step;
589 
590  }
591 
592  *condition = StronglyForced;
593  return StepLimit;
594 }
595 
596 // GetAverageNumberOfPhotons
597 // -------------------------
598 // This routine computes the number of Cerenkov photons produced per
599 // GEANT-unit (millimeter) in the current medium.
600 // ^^^^^^^^^^
601 
602 G4double
603  G4Cerenkov::GetAverageNumberOfPhotons(const G4double charge,
604  const G4double beta,
605  const G4Material* aMaterial,
606  G4MaterialPropertyVector* Rindex) const
607 {
608  const G4double Rfact = 369.81/(eV * cm);
609 
610  if(beta <= 0.0)return 0.0;
611 
612  G4double BetaInverse = 1./beta;
613 
614  // Vectors used in computation of Cerenkov Angle Integral:
615  // - Refraction Indices for the current material
616  // - new G4PhysicsOrderedFreeVector allocated to hold CAI's
617 
618  G4int materialIndex = aMaterial->GetIndex();
619 
620  // Retrieve the Cerenkov Angle Integrals for this material
621 
622  G4PhysicsOrderedFreeVector* CerenkovAngleIntegrals =
623  (G4PhysicsOrderedFreeVector*)((*thePhysicsTable)(materialIndex));
624 
625  if(!(CerenkovAngleIntegrals->IsFilledVectorExist()))return 0.0;
626 
627  // Min and Max photon energies
628  G4double Pmin = Rindex->GetMinLowEdgeEnergy();
629  G4double Pmax = Rindex->GetMaxLowEdgeEnergy();
630 
631  // Min and Max Refraction Indices
632  G4double nMin = Rindex->GetMinValue();
633  G4double nMax = Rindex->GetMaxValue();
634 
635  // Max Cerenkov Angle Integral
636  G4double CAImax = CerenkovAngleIntegrals->GetMaxValue();
637 
638  G4double dp, ge;
639 
640  // If n(Pmax) < 1/Beta -- no photons generated
641 
642  if (nMax < BetaInverse) {
643  dp = 0.0;
644  ge = 0.0;
645  }
646 
647  // otherwise if n(Pmin) >= 1/Beta -- photons generated
648 
649  else if (nMin > BetaInverse) {
650  dp = Pmax - Pmin;
651  ge = CAImax;
652  }
653 
654  // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then
655  // we need to find a P such that the value of n(P) == 1/Beta.
656  // Interpolation is performed by the GetEnergy() and
657  // Value() methods of the G4MaterialPropertiesTable and
658  // the GetValue() method of G4PhysicsVector.
659 
660  else {
661  Pmin = Rindex->GetEnergy(BetaInverse);
662  dp = Pmax - Pmin;
663 
664  // need boolean for current implementation of G4PhysicsVector
665  // ==> being phased out
666  G4bool isOutRange;
667  G4double CAImin = CerenkovAngleIntegrals->GetValue(Pmin, isOutRange);
668  ge = CAImax - CAImin;
669 
670  if (verboseLevel>0) {
671  G4cout << "CAImin = " << CAImin << G4endl;
672  G4cout << "ge = " << ge << G4endl;
673  }
674  }
675 
676  // Calculate number of photons
677  G4double NumPhotons = Rfact * charge/eplus * charge/eplus *
678  (dp - ge * BetaInverse*BetaInverse);
679 
680  return NumPhotons;
681 }
G4double G4ParticleHPJENDLHEData::G4double result
G4double condition(const G4ErrorSymMatrix &m)
const int N
Definition: mixmax.h:43
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:150
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
static G4LossTableManager * Instance()
G4double GetValue(G4double theEnergy, G4bool &isOutRange) const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetKineticEnergy() const
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetTotalEnergy() const
G4double GetStepLength() const
double x() const
const G4DynamicParticle * GetDynamicParticle() const
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:145
size_t GetIndex() const
Definition: G4Material.hh:262
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:160
G4TrackStatus GetTrackStatus() const
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
std::vector< G4Material * > G4MaterialTable
void SetTouchableHandle(const G4TouchableHandle &apValue)
void InsertValues(G4double energy, G4double value)
G4ParticleDefinition * GetDefinition() const
G4double GetVelocity() const
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:212
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double GetTotalMomentum() const
G4StepPoint * GetPreStepPoint() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:485
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
static constexpr double cm
Definition: G4SIunits.hh:119
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static constexpr double eplus
Definition: G4SIunits.hh:199
static constexpr double eV
Definition: G4SIunits.hh:215
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static constexpr double perCent
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:155
Definition: G4Step.hh:76
G4int GetTrackID() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void SetKineticEnergy(G4double aEnergy)
G4double GetEnergy(G4double aValue)
G4Material * GetMaterial() const
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition: G4Cerenkov.cc:492
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:133
double y() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4Cerenkov.cc:169
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
void insertAt(size_t, G4PhysicsVector *)
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition: G4Cerenkov.cc:96
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4double GetPDGCharge() const
G4ThreeVector GetDeltaPosition() const
G4bool IsFilledVectorExist() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4MaterialPropertyVector * GetProperty(const char *key)
const G4TouchableHandle & GetTouchableHandle() const
void clearAndDestroy()
G4ProcessType
G4double GetBeta() const