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

Last change on this file since 880 was 819, checked in by garnier, 17 years ago

import all except CVS

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