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