source: trunk/source/processes/electromagnetic/xrays/src/G4Cerenkov.cc@ 1315

Last change on this file since 1315 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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