Public Member Functions | |
ReweightProcess (G4ProcessManager *_pmgr, G4VProcess *proc, double _ratio) | |
virtual void | StartTracking (G4Track *track) |
virtual void | EndTracking () |
virtual G4double | AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) |
virtual G4VParticleChange * | AlongStepDoIt (const G4Track &track, const G4Step &stepData) |
virtual G4double | PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) |
virtual G4VParticleChange * | PostStepDoIt (const G4Track &track, const G4Step &stepData) |
virtual G4double | AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) |
virtual G4VParticleChange * | AtRestDoIt (const G4Track &track, const G4Step &stepData) |
Private Attributes | |
double | ratio |
G4ProcessManager * | pmgr |
ReweightProcess::ReweightProcess | ( | G4ProcessManager * | _pmgr, | |
G4VProcess * | proc, | |||
double | _ratio | |||
) |
virtual void ReweightProcess::StartTracking | ( | G4Track * | track | ) | [inline, virtual] |
virtual void ReweightProcess::EndTracking | ( | ) | [inline, virtual] |
G4double ReweightProcess::AlongStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4double | currentMinimumStep, | |||
G4double & | proposedSafety, | |||
G4GPILSelection * | selection | |||
) | [virtual] |
G4VParticleChange * ReweightProcess::AlongStepDoIt | ( | const G4Track & | track, | |
const G4Step & | stepData | |||
) | [virtual] |
References BLAssert, and pmgr.
00290 { 00291 G4VParticleChange *pc = &aParticleChange; 00292 pc->Initialize(track); 00293 00294 // Re-weight the track because it survived the artificial interaction 00295 // length. Do it here, because all GPIL routines have been called, 00296 // but no PostStepDoIt() or AtRestDoIt() have been called. So all 00297 // processes have computed their currentInteractionLength, and this 00298 // re-weighting will apply to all discrete processes. 00299 // Note that only the first ReweightProcess in the list will do the 00300 // re-weight; others return with no change. 00301 G4double effectiveIL=0.0, realIL=0.0; 00302 int n=pmgr->GetProcessListLength(); 00303 G4ProcessVector *pv=pmgr->GetProcessList(); 00304 int nReWeight = 0; 00305 for(int i=0; i<n; ++i) { 00306 G4VProcess *p=(*pv)[i]; 00307 G4double currentIL=p->GetCurrentInteractionLength(); 00308 if(currentIL <= 0.0) continue; 00309 effectiveIL += 1.0/currentIL; 00310 ReweightProcess *rw = dynamic_cast<ReweightProcess*>(p); 00311 if(rw != 0) { 00312 if(nReWeight++ != 0 && rw == this) return pc; 00313 currentIL = rw->pRegProcess->GetCurrentInteractionLength(); 00314 BLAssert(currentIL > DBL_MIN); 00315 } 00316 realIL += 1.0/currentIL; 00317 } 00318 effectiveIL = 1.0/effectiveIL; 00319 realIL = 1.0/realIL; 00320 G4double stepLength = stepData.GetStepLength(); 00321 G4double factor = exp(stepLength/effectiveIL - stepLength/realIL); 00322 pc->ProposeParentWeight(track.GetWeight()*factor); 00323 00324 return pc; 00325 }
G4double ReweightProcess::PostStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4ForceCondition * | condition | |||
) | [virtual] |
References ratio.
00329 { 00330 if(previousStepSize > 0.0) { 00331 SubtractNumberOfInteractionLengthLeft(previousStepSize); 00332 if(theNumberOfInteractionLengthLeft < 0.0) { 00333 theNumberOfInteractionLengthLeft=perMillion; 00334 } 00335 } 00336 00337 pRegProcess->PostStepGetPhysicalInteractionLength(track, 00338 previousStepSize, condition); 00339 currentInteractionLength = pRegProcess->GetCurrentInteractionLength(); 00340 if(currentInteractionLength < 0.0) 00341 G4Exception("reweightprocess","Process cannot be re-weighted", 00342 FatalException, pRegProcess->GetProcessName().c_str()); 00343 currentInteractionLength /= ratio; 00344 00345 G4double value = theNumberOfInteractionLengthLeft * 00346 currentInteractionLength; 00347 return value; 00348 }
G4VParticleChange * ReweightProcess::PostStepDoIt | ( | const G4Track & | track, | |
const G4Step & | stepData | |||
) | [virtual] |
References ratio.
00352 { 00353 G4VParticleChange *pc = pRegProcess->PostStepDoIt(track,stepData); 00354 00355 // re-weight track 00356 pc->ProposeParentWeight(track.GetWeight()/ratio); 00357 // (Geant4 bug requires the following update to the weight, sometimes.) 00358 stepData.GetPostStepPoint()->SetWeight(track.GetWeight()/ratio); 00359 00360 // re-weight secondaries 00361 int n=pc->GetNumberOfSecondaries(); 00362 for(int i=0; i<n; ++i) { 00363 G4Track *t=pc->GetSecondary(i); 00364 t->SetWeight(t->GetWeight()/ratio); 00365 } 00366 00367 ClearNumberOfInteractionLengthLeft(); 00368 return pc; 00369 }
G4double ReweightProcess::AtRestGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4ForceCondition * | condition | |||
) | [virtual] |
References ratio.
00373 { 00374 G4double value = pRegProcess-> 00375 AtRestGetPhysicalInteractionLength(track,condition); 00376 00377 value /= ratio; 00378 00379 return value; 00380 }
G4VParticleChange * ReweightProcess::AtRestDoIt | ( | const G4Track & | track, | |
const G4Step & | stepData | |||
) | [virtual] |
References ratio.
00384 { 00385 G4VParticleChange *pc = pRegProcess->AtRestDoIt(track,stepData); 00386 00387 // re-weight track 00388 pc->ProposeParentWeight(track.GetWeight()/ratio); 00389 // (Geant4 bug requires the following update to the weight, sometimes.) 00390 stepData.GetPostStepPoint()->SetWeight(track.GetWeight()/ratio); 00391 00392 // re-weight secondaries 00393 int n=pc->GetNumberOfSecondaries(); 00394 for(int i=0; i<n; ++i) { 00395 G4Track *t=pc->GetSecondary(i); 00396 t->SetWeight(t->GetWeight()/ratio); 00397 } 00398 00399 ClearNumberOfInteractionLengthLeft(); 00400 return pc; 00401 }
double ReweightProcess::ratio [private] |
G4ProcessManager* ReweightProcess::pmgr [private] |
Referenced by AlongStepDoIt(), and ReweightProcess().