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