Geant4  10.01.p03
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 85355 2014-10-28 09:58:59Z 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 
81  // static data members
83 
84 //G4bool G4Cerenkov::fTrackSecondariesFirst = false;
85 //G4double G4Cerenkov::fMaxBetaChange = 0.;
86 //G4int G4Cerenkov::fMaxPhotons = 0;
87 
89  // Operators
91 
92 // G4Cerenkov::operator=(const G4Cerenkov &right)
93 // {
94 // }
95 
97  // Constructors
99 
101  : G4VProcess(processName, type) ,
102  fTrackSecondariesFirst(false),
103  fMaxBetaChange(0),
104  fMaxPhotons(0)
105 {
107 
108  thePhysicsTable = NULL;
109 
110  if (verboseLevel>0) {
111  G4cout << GetProcessName() << " is created " << G4endl;
112  }
113 }
114 
115 // G4Cerenkov::G4Cerenkov(const G4Cerenkov &right)
116 // {
117 // }
118 
120  // Destructors
122 
124 {
125  if (thePhysicsTable != NULL) {
127  delete thePhysicsTable;
128  }
129 }
130 
132  // Methods
134 
136 {
137  G4bool result = false;
138  if (aParticleType.GetPDGCharge() != 0.0 &&
139  aParticleType.GetPDGMass() != 0.0 &&
140  aParticleType.GetParticleName() != "chargedgeantino" &&
141  !aParticleType.IsShortLived() ) { result = true; }
142 
143  return result;
144 }
145 
147 {
148  fTrackSecondariesFirst = state;
149 }
150 
152 {
154 }
155 
157 {
158  fMaxPhotons = NumPhotons;
159 }
160 
162 {
164 }
165 
166 // PostStepDoIt
167 // -------------
168 //
170 G4Cerenkov::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
171 
172 // This routine is called for each tracking Step of a charged particle
173 // in a radiator. A Poisson-distributed number of photons is generated
174 // according to the Cerenkov formula, distributed evenly along the track
175 // segment and uniformly azimuth w.r.t. the particle direction. The
176 // parameters are then transformed into the Master Reference System, and
177 // they are added to the particle change.
178 
179 {
181  // Should we ensure that the material is dispersive?
183 
184  aParticleChange.Initialize(aTrack);
185 
186  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
187  const G4Material* aMaterial = aTrack.GetMaterial();
188 
189  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
190  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
191 
192  G4ThreeVector x0 = pPreStepPoint->GetPosition();
193  G4ThreeVector p0 = aStep.GetDeltaPosition().unit();
194  G4double t0 = pPreStepPoint->GetGlobalTime();
195 
196  G4MaterialPropertiesTable* aMaterialPropertiesTable =
197  aMaterial->GetMaterialPropertiesTable();
198  if (!aMaterialPropertiesTable) return pParticleChange;
199 
200  G4MaterialPropertyVector* Rindex =
201  aMaterialPropertiesTable->GetProperty("RINDEX");
202  if (!Rindex) return pParticleChange;
203 
204  // particle charge
205  const G4double charge = aParticle->GetDefinition()->GetPDGCharge();
206 
207  // particle beta
208  const G4double beta = (pPreStepPoint ->GetBeta() +
209  pPostStepPoint->GetBeta())/2.;
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;
225  step_length = aStep.GetStepLength();
226 
227  MeanNumberOfPhotons = MeanNumberOfPhotons * step_length;
228 
229  G4int NumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons);
230 
231  if (NumPhotons <= 0) {
232 
233  // return unchanged particle and no secondaries
234 
236 
237  return pParticleChange;
238  }
239 
241 
243 
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  const G4double beta1 = pPreStepPoint ->GetBeta();
263  const 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 < NumPhotons; 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  } while (rand*maxSin2 > sin2Theta);
290 
291  // Generate random position of photon on cone surface
292  // defined by Theta
293 
294  rand = G4UniformRand();
295 
296  G4double phi = twopi*rand;
297  G4double sinPhi = std::sin(phi);
298  G4double cosPhi = std::cos(phi);
299 
300  // calculate x,y, and z components of photon energy
301  // (in coord system with primary particle direction
302  // aligned with the z axis)
303 
304  G4double sinTheta = std::sqrt(sin2Theta);
305  G4double px = sinTheta*cosPhi;
306  G4double py = sinTheta*sinPhi;
307  G4double pz = cosTheta;
308 
309  // Create photon momentum direction vector
310  // The momentum direction is still with respect
311  // to the coordinate system where the primary
312  // particle direction is aligned with the z axis
313 
314  G4ParticleMomentum photonMomentum(px, py, pz);
315 
316  // Rotate momentum direction back to global reference
317  // system
318 
319  photonMomentum.rotateUz(p0);
320 
321  // Determine polarization of new photon
322 
323  G4double sx = cosTheta*cosPhi;
324  G4double sy = cosTheta*sinPhi;
325  G4double sz = -sinTheta;
326 
327  G4ThreeVector photonPolarization(sx, sy, sz);
328 
329  // Rotate back to original coord system
330 
331  photonPolarization.rotateUz(p0);
332 
333  // Generate a new photon:
334 
335  G4DynamicParticle* aCerenkovPhoton =
337  photonMomentum);
338  aCerenkovPhoton->SetPolarization
339  (photonPolarization.x(),
340  photonPolarization.y(),
341  photonPolarization.z());
342 
343  aCerenkovPhoton->SetKineticEnergy(sampledEnergy);
344 
345  // Generate new G4Track object:
346 
347  G4double delta, NumberOfPhotons, N;
348 
349  do {
350  rand = G4UniformRand();
351  delta = rand * aStep.GetStepLength();
352  NumberOfPhotons = MeanNumberOfPhotons1 - delta *
353  (MeanNumberOfPhotons1-MeanNumberOfPhotons2)/
354  aStep.GetStepLength();
355  N = G4UniformRand() *
356  std::max(MeanNumberOfPhotons1,MeanNumberOfPhotons2);
357  } while (N > NumberOfPhotons);
358 
359  G4double deltaTime = delta /
360  ((pPreStepPoint->GetVelocity()+
361  pPostStepPoint->GetVelocity())/2.);
362 
363  G4double aSecondaryTime = t0 + deltaTime;
364 
365  G4ThreeVector aSecondaryPosition =
366  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 
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 
419  aPhysicsOrderedFreeVector = new G4PhysicsOrderedFreeVector();
420  G4MaterialPropertyVector* theRefractionIndexVector =
421  aMaterialPropertiesTable->GetProperty("RINDEX");
422 
423  if (theRefractionIndexVector) {
424 
425  // Retrieve the first refraction index in vector
426  // of (photon energy, refraction index) pairs
427 
428  G4double currentRI = (*theRefractionIndexVector)[0];
429 
430  if (currentRI > 1.0) {
431 
432  // Create first (photon energy, Cerenkov Integral)
433  // pair
434 
435  G4double currentPM = theRefractionIndexVector->
436  Energy(0);
437  G4double currentCAI = 0.0;
438 
439  aPhysicsOrderedFreeVector->
440  InsertValues(currentPM , currentCAI);
441 
442  // Set previous values to current ones prior to loop
443 
444  G4double prevPM = currentPM;
445  G4double prevCAI = currentCAI;
446  G4double prevRI = currentRI;
447 
448  // loop over all (photon energy, refraction index)
449  // pairs stored for this material
450 
451  for (size_t ii = 1;
452  ii < theRefractionIndexVector->GetVectorLength();
453  ++ii)
454  {
455  currentRI = (*theRefractionIndexVector)[ii];
456  currentPM = theRefractionIndexVector->Energy(ii);
457 
458  currentCAI = 0.5*(1.0/(prevRI*prevRI) +
459  1.0/(currentRI*currentRI));
460 
461  currentCAI = prevCAI +
462  (currentPM - prevPM) * currentCAI;
463 
464  aPhysicsOrderedFreeVector->
465  InsertValues(currentPM, currentCAI);
466 
467  prevPM = currentPM;
468  prevCAI = currentCAI;
469  prevRI = currentRI;
470  }
471 
472  }
473  }
474  }
475 
476  // The Cerenkov integral for a given material
477  // will be inserted in thePhysicsTable
478  // according to the position of the material in
479  // the material table.
480 
481  thePhysicsTable->insertAt(i,aPhysicsOrderedFreeVector);
482 
483  }
484 }
485 
486 // GetMeanFreePath
487 // ---------------
488 //
489 
491  G4double,
493 {
494  return 1.;
495 }
496 
498  const G4Track& aTrack,
499  G4double,
501 {
502  *condition = NotForced;
503  G4double StepLimit = DBL_MAX;
504 
505  const G4Material* aMaterial = aTrack.GetMaterial();
506  G4int materialIndex = aMaterial->GetIndex();
507 
508  // If Physics Vector is not defined no Cerenkov photons
509  // this check avoid string comparison below
510  if(!(*thePhysicsTable)[materialIndex]) { return StepLimit; }
511 
512  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
513  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
514 
515  G4double kineticEnergy = aParticle->GetKineticEnergy();
516  const G4ParticleDefinition* particleType = aParticle->GetDefinition();
517  G4double mass = particleType->GetPDGMass();
518 
519  // particle beta
520  G4double beta = aParticle->GetTotalMomentum() /
521  aParticle->GetTotalEnergy();
522  // particle gamma
523  G4double gamma = aParticle->GetTotalEnergy()/mass;
524 
525  G4MaterialPropertiesTable* aMaterialPropertiesTable =
526  aMaterial->GetMaterialPropertiesTable();
527 
528  G4MaterialPropertyVector* Rindex = NULL;
529 
530  if (aMaterialPropertiesTable)
531  Rindex = aMaterialPropertiesTable->GetProperty("RINDEX");
532 
533  G4double nMax;
534  if (Rindex) {
535  nMax = Rindex->GetMaxValue();
536  } else {
537  return StepLimit;
538  }
539 
540  G4double BetaMin = 1./nMax;
541  if ( BetaMin >= 1. ) return StepLimit;
542 
543  G4double GammaMin = 1./std::sqrt(1.-BetaMin*BetaMin);
544 
545  if (gamma < GammaMin ) return StepLimit;
546 
547  G4double kinEmin = mass*(GammaMin-1.);
548 
550  GetRange(particleType,
551  kinEmin,
552  couple);
554  GetRange(particleType,
555  kineticEnergy,
556  couple);
557 
558  G4double Step = Range - RangeMin;
559  if (Step < 1.*um ) return StepLimit;
560 
561  if (Step > 0. && Step < StepLimit) StepLimit = Step;
562 
563  // If user has defined an average maximum number of photons to
564  // be generated in a Step, then calculate the Step length for
565  // that number of photons.
566 
567  if (fMaxPhotons > 0) {
568 
569  // particle charge
570  const G4double charge = aParticle->
571  GetDefinition()->GetPDGCharge();
572 
573  G4double MeanNumberOfPhotons =
574  GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex);
575 
576  Step = 0.;
577  if (MeanNumberOfPhotons > 0.0) Step = fMaxPhotons /
578  MeanNumberOfPhotons;
579 
580  if (Step > 0. && Step < StepLimit) StepLimit = Step;
581  }
582 
583  // If user has defined an maximum allowed change in beta per step
584  if (fMaxBetaChange > 0.) {
585 
587  GetDEDX(particleType,
588  kineticEnergy,
589  couple);
590 
591  G4double deltaGamma = gamma -
592  1./std::sqrt(1.-beta*beta*
593  (1.-fMaxBetaChange)*
594  (1.-fMaxBetaChange));
595 
596  Step = mass * deltaGamma / dedx;
597 
598  if (Step > 0. && Step < StepLimit) StepLimit = Step;
599 
600  }
601 
602  *condition = StronglyForced;
603  return StepLimit;
604 }
605 
606 // GetAverageNumberOfPhotons
607 // -------------------------
608 // This routine computes the number of Cerenkov photons produced per
609 // GEANT-unit (millimeter) in the current medium.
610 // ^^^^^^^^^^
611 
612 G4double
614  const G4double beta,
615  const G4Material* aMaterial,
616  G4MaterialPropertyVector* Rindex) const
617 {
618  const G4double Rfact = 369.81/(eV * cm);
619 
620  if(beta <= 0.0)return 0.0;
621 
622  G4double BetaInverse = 1./beta;
623 
624  // Vectors used in computation of Cerenkov Angle Integral:
625  // - Refraction Indices for the current material
626  // - new G4PhysicsOrderedFreeVector allocated to hold CAI's
627 
628  G4int materialIndex = aMaterial->GetIndex();
629 
630  // Retrieve the Cerenkov Angle Integrals for this material
631 
632  G4PhysicsOrderedFreeVector* CerenkovAngleIntegrals =
633  (G4PhysicsOrderedFreeVector*)((*thePhysicsTable)(materialIndex));
634 
635  if(!(CerenkovAngleIntegrals->IsFilledVectorExist()))return 0.0;
636 
637  // Min and Max photon energies
638  G4double Pmin = Rindex->GetMinLowEdgeEnergy();
639  G4double Pmax = Rindex->GetMaxLowEdgeEnergy();
640 
641  // Min and Max Refraction Indices
642  G4double nMin = Rindex->GetMinValue();
643  G4double nMax = Rindex->GetMaxValue();
644 
645  // Max Cerenkov Angle Integral
646  G4double CAImax = CerenkovAngleIntegrals->GetMaxValue();
647 
648  G4double dp, ge;
649 
650  // If n(Pmax) < 1/Beta -- no photons generated
651 
652  if (nMax < BetaInverse) {
653  dp = 0;
654  ge = 0;
655  }
656 
657  // otherwise if n(Pmin) >= 1/Beta -- photons generated
658 
659  else if (nMin > BetaInverse) {
660  dp = Pmax - Pmin;
661  ge = CAImax;
662  }
663 
664  // If n(Pmin) < 1/Beta, and n(Pmax) >= 1/Beta, then
665  // we need to find a P such that the value of n(P) == 1/Beta.
666  // Interpolation is performed by the GetEnergy() and
667  // Value() methods of the G4MaterialPropertiesTable and
668  // the GetValue() method of G4PhysicsVector.
669 
670  else {
671  Pmin = Rindex->GetEnergy(BetaInverse);
672  dp = Pmax - Pmin;
673 
674  // need boolean for current implementation of G4PhysicsVector
675  // ==> being phased out
676  G4bool isOutRange;
677  G4double CAImin = CerenkovAngleIntegrals->
678  GetValue(Pmin, isOutRange);
679  ge = CAImax - CAImin;
680 
681  if (verboseLevel>0) {
682  G4cout << "CAImin = " << CAImin << G4endl;
683  G4cout << "ge = " << ge << G4endl;
684  }
685  }
686 
687  // Calculate number of photons
688  G4double NumPhotons = Rfact * charge/eplus * charge/eplus *
689  (dp - ge * BetaInverse*BetaInverse);
690 
691  return NumPhotons;
692 }
static const double cm
Definition: G4SIunits.hh:106
G4double condition(const G4ErrorSymMatrix &m)
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:151
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
static G4LossTableManager * Instance()
G4MaterialPropertyVector * GetProperty(const char *key)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:146
size_t GetIndex() const
Definition: G4Material.hh:262
G4double GetAverageNumberOfPhotons(const G4double charge, const G4double beta, const G4Material *aMaterial, G4MaterialPropertyVector *Rindex) const
Definition: G4Cerenkov.cc:613
G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:497
G4TrackStatus GetTrackStatus() const
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:588
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
std::vector< G4Material * > G4MaterialTable
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4ParticleDefinition * GetDefinition() const
G4double GetVelocity() const
G4PhysicsTable * thePhysicsTable
Definition: G4Cerenkov.hh:210
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int fMaxPhotons
Definition: G4Cerenkov.hh:219
G4double GetTotalMomentum() const
G4StepPoint * GetPreStepPoint() const
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4Cerenkov.cc:490
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
G4double fMaxBetaChange
Definition: G4Cerenkov.hh:218
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double perCent
Definition: G4SIunits.hh:296
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:156
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:595
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void SetKineticEnergy(G4double aEnergy)
G4double GetEnergy(G4double aValue)
G4Material * GetMaterial() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Cerenkov.cc:170
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
static const double eV
Definition: G4SIunits.hh:194
const G4double p0
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)
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
Definition: Step.hh:41
void AddSecondary(G4Track *aSecondary)
void insertAt(size_t, G4PhysicsVector *)
G4Cerenkov(const G4String &processName="Cerenkov", G4ProcessType type=fElectromagnetic)
Definition: G4Cerenkov.cc:100
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:161
static const double eplus
Definition: G4SIunits.hh:178
G4ForceCondition
G4double GetPDGCharge() const
G4ThreeVector G4ParticleMomentum
G4ThreeVector GetDeltaPosition() const
G4bool IsFilledVectorExist() const
#define DBL_MAX
Definition: templates.hh:83
const G4TouchableHandle & GetTouchableHandle() const
void clearAndDestroy()
void BuildThePhysicsTable()
Definition: G4Cerenkov.cc:391
G4ProcessType
G4double GetBeta() const
G4bool fTrackSecondariesFirst
Definition: G4Cerenkov.hh:217
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:135