#include <G4Cerenkov.hh>
Inheritance diagram for G4Cerenkov:
Definition at line 81 of file G4Cerenkov.hh.
G4Cerenkov::G4Cerenkov | ( | const G4String & | processName = "Cerenkov" , |
|
G4ProcessType | type = fElectromagnetic | |||
) |
Definition at line 92 of file G4Cerenkov.cc.
References fCerenkov, G4cout, G4endl, G4VProcess::GetProcessName(), G4VProcess::SetProcessSubType(), thePhysicsTable, and G4VProcess::verboseLevel.
00093 : G4VProcess(processName, type) 00094 { 00095 SetProcessSubType(fCerenkov); 00096 00097 fTrackSecondariesFirst = false; 00098 fMaxBetaChange = 0.; 00099 fMaxPhotons = 0; 00100 00101 thePhysicsTable = NULL; 00102 00103 if (verboseLevel>0) { 00104 G4cout << GetProcessName() << " is created " << G4endl; 00105 } 00106 00107 BuildThePhysicsTable(); 00108 }
G4Cerenkov::~G4Cerenkov | ( | ) |
Definition at line 118 of file G4Cerenkov.cc.
References G4PhysicsTable::clearAndDestroy(), and thePhysicsTable.
00119 { 00120 if (thePhysicsTable != NULL) { 00121 thePhysicsTable->clearAndDestroy(); 00122 delete thePhysicsTable; 00123 } 00124 }
G4Cerenkov::G4Cerenkov | ( | const G4Cerenkov & | right | ) |
virtual G4VParticleChange* G4Cerenkov::AlongStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
virtual G4VParticleChange* G4Cerenkov::AtRestDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
virtual G4double G4Cerenkov::AtRestGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [inline, virtual] |
void G4Cerenkov::DumpPhysicsTable | ( | ) | const [inline] |
Definition at line 241 of file G4Cerenkov.hh.
References G4PhysicsOrderedFreeVector::DumpValues(), G4PhysicsTable::entries(), and thePhysicsTable.
00242 { 00243 G4int PhysicsTableSize = thePhysicsTable->entries(); 00244 G4PhysicsOrderedFreeVector *v; 00245 00246 for (G4int i = 0 ; i < PhysicsTableSize ; i++ ) 00247 { 00248 v = (G4PhysicsOrderedFreeVector*)(*thePhysicsTable)[i]; 00249 v->DumpValues(); 00250 } 00251 }
G4double G4Cerenkov::GetMeanFreePath | ( | const G4Track & | aTrack, | |
G4double | , | |||
G4ForceCondition * | ||||
) |
G4PhysicsTable * G4Cerenkov::GetPhysicsTable | ( | ) | const [inline] |
Definition at line 254 of file G4Cerenkov.hh.
References thePhysicsTable.
00255 { 00256 return thePhysicsTable; 00257 }
G4bool G4Cerenkov::IsApplicable | ( | const G4ParticleDefinition & | aParticleType | ) | [inline, virtual] |
Reimplemented from G4VProcess.
Definition at line 214 of file G4Cerenkov.hh.
References G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), and G4ParticleDefinition::IsShortLived().
Referenced by G4OpticalPhysics::ConstructProcess().
00215 { 00216 if (aParticleType.GetParticleName() == "chargedgeantino") return false; 00217 if (aParticleType.IsShortLived()) return false; 00218 00219 return (aParticleType.GetPDGCharge() != 0); 00220 }
G4VParticleChange * G4Cerenkov::PostStepDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) | [virtual] |
Implements G4VProcess.
Definition at line 134 of file G4Cerenkov.cc.
References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, fAlive, fSuspend, G4cout, G4endl, G4Poisson(), G4UniformRand, G4StepPoint::GetBeta(), G4DynamicParticle::GetDefinition(), G4Step::GetDeltaPosition(), G4Track::GetDynamicParticle(), G4StepPoint::GetGlobalTime(), G4Track::GetMaterial(), G4Material::GetMaterialPropertiesTable(), G4PhysicsOrderedFreeVector::GetMaxLowEdgeEnergy(), G4PhysicsOrderedFreeVector::GetMaxValue(), G4PhysicsOrderedFreeVector::GetMinLowEdgeEnergy(), G4VParticleChange::GetNumberOfSecondaries(), G4ParticleDefinition::GetPDGCharge(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4MaterialPropertiesTable::GetProperty(), G4Step::GetStepLength(), G4StepPoint::GetTouchableHandle(), G4Track::GetTrackID(), G4Track::GetTrackStatus(), G4StepPoint::GetVelocity(), G4ParticleChange::Initialize(), G4OpticalPhoton::OpticalPhoton(), G4VProcess::pParticleChange, G4VParticleChange::ProposeTrackStatus(), G4DynamicParticle::SetKineticEnergy(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetParentID(), G4DynamicParticle::SetPolarization(), G4Track::SetTouchableHandle(), G4PhysicsVector::Value(), and G4VProcess::verboseLevel.
00143 { 00145 // Should we ensure that the material is dispersive? 00147 00148 aParticleChange.Initialize(aTrack); 00149 00150 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle(); 00151 const G4Material* aMaterial = aTrack.GetMaterial(); 00152 00153 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint(); 00154 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint(); 00155 00156 G4ThreeVector x0 = pPreStepPoint->GetPosition(); 00157 G4ThreeVector p0 = aStep.GetDeltaPosition().unit(); 00158 G4double t0 = pPreStepPoint->GetGlobalTime(); 00159 00160 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00161 aMaterial->GetMaterialPropertiesTable(); 00162 if (!aMaterialPropertiesTable) return pParticleChange; 00163 00164 G4MaterialPropertyVector* Rindex = 00165 aMaterialPropertiesTable->GetProperty("RINDEX"); 00166 if (!Rindex) return pParticleChange; 00167 00168 // particle charge 00169 const G4double charge = aParticle->GetDefinition()->GetPDGCharge(); 00170 00171 // particle beta 00172 const G4double beta = (pPreStepPoint ->GetBeta() + 00173 pPostStepPoint->GetBeta())/2.; 00174 00175 G4double MeanNumberOfPhotons = 00176 GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex); 00177 00178 if (MeanNumberOfPhotons <= 0.0) { 00179 00180 // return unchanged particle and no secondaries 00181 00182 aParticleChange.SetNumberOfSecondaries(0); 00183 00184 return pParticleChange; 00185 00186 } 00187 00188 G4double step_length; 00189 step_length = aStep.GetStepLength(); 00190 00191 MeanNumberOfPhotons = MeanNumberOfPhotons * step_length; 00192 00193 G4int NumPhotons = (G4int) G4Poisson(MeanNumberOfPhotons); 00194 00195 if (NumPhotons <= 0) { 00196 00197 // return unchanged particle and no secondaries 00198 00199 aParticleChange.SetNumberOfSecondaries(0); 00200 00201 return pParticleChange; 00202 } 00203 00205 00206 aParticleChange.SetNumberOfSecondaries(NumPhotons); 00207 00208 if (fTrackSecondariesFirst) { 00209 if (aTrack.GetTrackStatus() == fAlive ) 00210 aParticleChange.ProposeTrackStatus(fSuspend); 00211 } 00212 00214 00215 G4double Pmin = Rindex->GetMinLowEdgeEnergy(); 00216 G4double Pmax = Rindex->GetMaxLowEdgeEnergy(); 00217 G4double dp = Pmax - Pmin; 00218 00219 G4double nMax = Rindex->GetMaxValue(); 00220 00221 G4double BetaInverse = 1./beta; 00222 00223 G4double maxCos = BetaInverse / nMax; 00224 G4double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos); 00225 00226 const G4double beta1 = pPreStepPoint ->GetBeta(); 00227 const G4double beta2 = pPostStepPoint->GetBeta(); 00228 00229 G4double MeanNumberOfPhotons1 = 00230 GetAverageNumberOfPhotons(charge,beta1,aMaterial,Rindex); 00231 G4double MeanNumberOfPhotons2 = 00232 GetAverageNumberOfPhotons(charge,beta2,aMaterial,Rindex); 00233 00234 for (G4int i = 0; i < NumPhotons; i++) { 00235 00236 // Determine photon energy 00237 00238 G4double rand; 00239 G4double sampledEnergy, sampledRI; 00240 G4double cosTheta, sin2Theta; 00241 00242 // sample an energy 00243 00244 do { 00245 rand = G4UniformRand(); 00246 sampledEnergy = Pmin + rand * dp; 00247 sampledRI = Rindex->Value(sampledEnergy); 00248 cosTheta = BetaInverse / sampledRI; 00249 00250 sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta); 00251 rand = G4UniformRand(); 00252 00253 } while (rand*maxSin2 > sin2Theta); 00254 00255 // Generate random position of photon on cone surface 00256 // defined by Theta 00257 00258 rand = G4UniformRand(); 00259 00260 G4double phi = twopi*rand; 00261 G4double sinPhi = std::sin(phi); 00262 G4double cosPhi = std::cos(phi); 00263 00264 // calculate x,y, and z components of photon energy 00265 // (in coord system with primary particle direction 00266 // aligned with the z axis) 00267 00268 G4double sinTheta = std::sqrt(sin2Theta); 00269 G4double px = sinTheta*cosPhi; 00270 G4double py = sinTheta*sinPhi; 00271 G4double pz = cosTheta; 00272 00273 // Create photon momentum direction vector 00274 // The momentum direction is still with respect 00275 // to the coordinate system where the primary 00276 // particle direction is aligned with the z axis 00277 00278 G4ParticleMomentum photonMomentum(px, py, pz); 00279 00280 // Rotate momentum direction back to global reference 00281 // system 00282 00283 photonMomentum.rotateUz(p0); 00284 00285 // Determine polarization of new photon 00286 00287 G4double sx = cosTheta*cosPhi; 00288 G4double sy = cosTheta*sinPhi; 00289 G4double sz = -sinTheta; 00290 00291 G4ThreeVector photonPolarization(sx, sy, sz); 00292 00293 // Rotate back to original coord system 00294 00295 photonPolarization.rotateUz(p0); 00296 00297 // Generate a new photon: 00298 00299 G4DynamicParticle* aCerenkovPhoton = 00300 new G4DynamicParticle(G4OpticalPhoton::OpticalPhoton(), 00301 photonMomentum); 00302 aCerenkovPhoton->SetPolarization 00303 (photonPolarization.x(), 00304 photonPolarization.y(), 00305 photonPolarization.z()); 00306 00307 aCerenkovPhoton->SetKineticEnergy(sampledEnergy); 00308 00309 // Generate new G4Track object: 00310 00311 G4double delta, NumberOfPhotons, N; 00312 00313 do { 00314 rand = G4UniformRand(); 00315 delta = rand * aStep.GetStepLength(); 00316 NumberOfPhotons = MeanNumberOfPhotons1 - delta * 00317 (MeanNumberOfPhotons1-MeanNumberOfPhotons2)/ 00318 aStep.GetStepLength(); 00319 N = G4UniformRand() * 00320 std::max(MeanNumberOfPhotons1,MeanNumberOfPhotons2); 00321 } while (N > NumberOfPhotons); 00322 00323 G4double deltaTime = delta / 00324 ((pPreStepPoint->GetVelocity()+ 00325 pPostStepPoint->GetVelocity())/2.); 00326 00327 G4double aSecondaryTime = t0 + deltaTime; 00328 00329 G4ThreeVector aSecondaryPosition = 00330 x0 + rand * aStep.GetDeltaPosition(); 00331 00332 G4Track* aSecondaryTrack = 00333 new G4Track(aCerenkovPhoton,aSecondaryTime,aSecondaryPosition); 00334 00335 aSecondaryTrack->SetTouchableHandle( 00336 aStep.GetPreStepPoint()->GetTouchableHandle()); 00337 00338 aSecondaryTrack->SetParentID(aTrack.GetTrackID()); 00339 00340 aParticleChange.AddSecondary(aSecondaryTrack); 00341 } 00342 00343 if (verboseLevel>0) { 00344 G4cout <<"\n Exiting from G4Cerenkov::DoIt -- NumberOfSecondaries = " 00345 << aParticleChange.GetNumberOfSecondaries() << G4endl; 00346 } 00347 00348 return pParticleChange; 00349 }
G4double G4Cerenkov::PostStepGetPhysicalInteractionLength | ( | const G4Track & | aTrack, | |
G4double | , | |||
G4ForceCondition * | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 461 of file G4Cerenkov.cc.
References DBL_MAX, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetMaterial(), G4Track::GetMaterialCutsCouple(), G4Material::GetMaterialPropertiesTable(), G4PhysicsOrderedFreeVector::GetMaxValue(), G4ParticleDefinition::GetPDGMass(), G4MaterialPropertiesTable::GetProperty(), G4DynamicParticle::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), G4LossTableManager::Instance(), NotForced, and StronglyForced.
00465 { 00466 *condition = NotForced; 00467 G4double StepLimit = DBL_MAX; 00468 00469 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle(); 00470 const G4Material* aMaterial = aTrack.GetMaterial(); 00471 const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple(); 00472 00473 G4double kineticEnergy = aParticle->GetKineticEnergy(); 00474 const G4ParticleDefinition* particleType = aParticle->GetDefinition(); 00475 G4double mass = particleType->GetPDGMass(); 00476 00477 // particle beta 00478 G4double beta = aParticle->GetTotalMomentum() / 00479 aParticle->GetTotalEnergy(); 00480 // particle gamma 00481 G4double gamma = aParticle->GetTotalEnergy()/mass; 00482 00483 G4MaterialPropertiesTable* aMaterialPropertiesTable = 00484 aMaterial->GetMaterialPropertiesTable(); 00485 00486 G4MaterialPropertyVector* Rindex = NULL; 00487 00488 if (aMaterialPropertiesTable) 00489 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX"); 00490 00491 G4double nMax; 00492 if (Rindex) { 00493 nMax = Rindex->GetMaxValue(); 00494 } else { 00495 return StepLimit; 00496 } 00497 00498 G4double BetaMin = 1./nMax; 00499 if ( BetaMin >= 1. ) return StepLimit; 00500 00501 G4double GammaMin = 1./std::sqrt(1.-BetaMin*BetaMin); 00502 00503 if (gamma < GammaMin ) return StepLimit; 00504 00505 G4double kinEmin = mass*(GammaMin-1.); 00506 00507 G4double RangeMin = G4LossTableManager::Instance()-> 00508 GetRange(particleType, 00509 kinEmin, 00510 couple); 00511 G4double Range = G4LossTableManager::Instance()-> 00512 GetRange(particleType, 00513 kineticEnergy, 00514 couple); 00515 00516 G4double Step = Range - RangeMin; 00517 if (Step < 1.*um ) return StepLimit; 00518 00519 if (Step > 0. && Step < StepLimit) StepLimit = Step; 00520 00521 // If user has defined an average maximum number of photons to 00522 // be generated in a Step, then calculate the Step length for 00523 // that number of photons. 00524 00525 if (fMaxPhotons > 0) { 00526 00527 // particle charge 00528 const G4double charge = aParticle-> 00529 GetDefinition()->GetPDGCharge(); 00530 00531 G4double MeanNumberOfPhotons = 00532 GetAverageNumberOfPhotons(charge,beta,aMaterial,Rindex); 00533 00534 Step = 0.; 00535 if (MeanNumberOfPhotons > 0.0) Step = fMaxPhotons / 00536 MeanNumberOfPhotons; 00537 00538 if (Step > 0. && Step < StepLimit) StepLimit = Step; 00539 } 00540 00541 // If user has defined an maximum allowed change in beta per step 00542 if (fMaxBetaChange > 0.) { 00543 00544 G4double dedx = G4LossTableManager::Instance()-> 00545 GetDEDX(particleType, 00546 kineticEnergy, 00547 couple); 00548 00549 G4double deltaGamma = gamma - 00550 1./std::sqrt(1.-beta*beta* 00551 (1.-fMaxBetaChange)* 00552 (1.-fMaxBetaChange)); 00553 00554 Step = mass * deltaGamma / dedx; 00555 00556 if (Step > 0. && Step < StepLimit) StepLimit = Step; 00557 00558 } 00559 00560 *condition = StronglyForced; 00561 return StepLimit; 00562 }
void G4Cerenkov::SetMaxBetaChangePerStep | ( | const G4double | d | ) | [inline] |
Definition at line 229 of file G4Cerenkov.hh.
Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetMaxBetaChangePerStep().
void G4Cerenkov::SetMaxNumPhotonsPerStep | ( | const G4int | NumPhotons | ) | [inline] |
Definition at line 235 of file G4Cerenkov.hh.
Referenced by G4OpticalPhysics::ConstructProcess(), and G4OpticalPhysics::SetMaxNumPhotonsPerStep().
void G4Cerenkov::SetTrackSecondariesFirst | ( | const G4bool | state | ) | [inline] |
Definition at line 223 of file G4Cerenkov.hh.
Referenced by G4OpticalPhysics::SetTrackSecondariesFirst().
G4PhysicsTable* G4Cerenkov::thePhysicsTable [protected] |
Definition at line 197 of file G4Cerenkov.hh.
Referenced by DumpPhysicsTable(), G4Cerenkov(), GetPhysicsTable(), and ~G4Cerenkov().