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

Last change on this file since 1347 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

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