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