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