ParticleFilterPlacement Class Reference

Inheritance diagram for ParticleFilterPlacement:

BLManager::SteppingAction BLManager::TrackingAction BLCallback

List of all members.


Detailed Description

class ParticleFilterPlacement is one placement of a ParticleFilter element.

It inherits from G4Decay which does the heavy lifting, all we do here is decide whether or not to decay. We also kill tracks in BLManager::SteppingAction, because you cannot do that inside PostStepGetPhysicalInteractionLength().

Because geometry construction time is too early to manipulate physics processes, it registers a callback with BLManager to do it before the reference particle.

Public Member Functions

 ParticleFilterPlacement (BLCMDparticlefilter *_pf, G4VPhysicalVolume *_pv, G4String _name)
void callback (int type)
G4bool IsApplicable (const G4ParticleDefinition &)
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void UserSteppingAction (const G4Step *step)
void PreUserTrackingAction (const G4Track *track)
void PostUserTrackingAction (const G4Track *track)

Private Attributes

BLCMDparticlefilterfilter
G4VPhysicalVolume * vol
unsigned long nDecayed
unsigned long nKilled
int nWait
int referenceWait


Constructor & Destructor Documentation

ParticleFilterPlacement::ParticleFilterPlacement ( BLCMDparticlefilter _pf,
G4VPhysicalVolume *  _pv,
G4String  _name 
)

References filter, BLManager::getObject(), nDecayed, nKilled, BLManager::registerCallback(), and vol.

00434                                                         : 
00435                                         G4Decay(_name), BLCallback()
00436 {
00437         filter=_pf;
00438         vol=_pv; 
00439         nDecayed = 0;
00440         nKilled = 0;
00441 
00442         // now is too early to manage the physics processes
00443         BLManager::getObject()->registerCallback(this,0);
00444 
00445         // also register a callback after tracking, to display statistics
00446         BLManager::getObject()->registerCallback(this,2);
00447 }


Member Function Documentation

void ParticleFilterPlacement::callback ( int  type  )  [virtual]

Reimplemented from BLCallback.

References filter, BLElement::getName(), IsApplicable(), nDecayed, and nKilled.

00450 {
00451         if(type == 0) {
00452             // add ourselves to the process list for each applicable particle
00453             // and remove the Decay process if necessary.
00454             G4ParticleTable::G4PTblDicIterator *theParticleIterator =
00455                         G4ParticleTable::GetParticleTable()->GetIterator();
00456             theParticleIterator->reset();
00457             while((*theParticleIterator)()) {
00458                 G4ParticleDefinition *pd = theParticleIterator->value();
00459                 if(IsApplicable(*pd)) {
00460                         G4ProcessManager *pmgr = pd->GetProcessManager();
00461 again:                  G4ProcessVector *pv = pmgr->GetProcessList();
00462                         int n = pv->size();
00463                         for(int i=0; i<n; ++i) {
00464                                 G4VProcess *process = (*pv)[i];
00465                                 if(process->GetProcessName()=="Decay") {
00466                                         pmgr->RemoveProcess(process);
00467                                         goto again;
00468                                 }
00469                         }
00470                         pmgr->AddProcess(this);
00471                         pmgr->SetProcessOrdering(this,idxPostStep);
00472                         pmgr->SetProcessOrdering(this,idxAtRest);
00473 //for(int i=0; i<pv->size(); ++i) printf("ProcessName='%s'\n",(*pv)[i]->GetProcessName().c_str());
00474                 }
00475             }
00476         } else if(type == 2) {
00477                 printf("particlefilter %s: %lu decays  %lu killed\n",
00478                         filter->getName().c_str(),nDecayed,nKilled);
00479         }
00480 }

G4bool ParticleFilterPlacement::IsApplicable ( const G4ParticleDefinition &  pd  ) 

References BLCMDparticlefilter::decaySet, and filter.

Referenced by callback().

00483 {
00484         if(filter->decaySet.count(pd.GetParticleName()) > 0)
00485                 return G4Decay::IsApplicable(pd);
00486         return false;
00487 }

G4double ParticleFilterPlacement::AtRestGetPhysicalInteractionLength ( const G4Track &  track,
G4ForceCondition *  condition 
)

References filter, BLElement::getName(), nDecayed, BLCMDparticlefilter::steppingVerbose, and vol.

00491 {
00492         // we are only applied to particlies listed in decaySet.
00493         if(track.GetVolume() == vol || track.GetNextVolume() == vol) {
00494                 ++nDecayed;
00495                 *condition = Forced;
00496                 if(filter->steppingVerbose) 
00497                         printf("particlefilter '%s' forced Decay\n",
00498                                         filter->getName().c_str());
00499                 return DBL_MIN;
00500         }
00501 
00502         *condition = NotForced;
00503         return DBL_MAX;
00504 }

G4double ParticleFilterPlacement::PostStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4ForceCondition *  condition 
)

References filter, BLElement::getName(), nDecayed, BLCMDparticlefilter::steppingVerbose, and vol.

00509 {
00510         // we are only applied to particlies listed in decaySet.
00511         if(track.GetVolume() == vol || track.GetNextVolume() == vol) {
00512                 ++nDecayed;
00513                 *condition = Forced;
00514                 if(filter->steppingVerbose) 
00515                         printf("particlefilter '%s' forced Decay\n",
00516                                         filter->getName().c_str());
00517                 return DBL_MIN;
00518         }
00519 
00520         *condition = NotForced;
00521         return DBL_MAX;
00522 }

void ParticleFilterPlacement::UserSteppingAction ( const G4Step *  step  )  [virtual]

Implements BLManager::SteppingAction.

References BEAM, BLEvaluator::evaluate(), filter, BLManager::getObject(), BLCMDparticlefilter::keepSet, BLCMDparticlefilter::killSet, nKilled, nWait, referenceWait, BLCMDparticlefilter::require, BLCMDparticlefilter::steppingVerbose, VISUAL, and vol.

00524 {
00525         // handle all states
00526 
00527         G4Track *track = step->GetTrack();
00528         G4VPhysicalVolume *preVol = 
00529                                 step->GetPreStepPoint()->GetPhysicalVolume();
00530         G4VPhysicalVolume *postVol = 
00531                                 step->GetPostStepPoint()->GetPhysicalVolume();
00532 
00533         // only handle the step that enters the volume
00534         if(preVol == postVol || postVol != vol) return;
00535 
00536 //printf("ParticleFilterPlacement nWait=%d referenceWait=%d\n",nWait,referenceWait);
00537         // wait until track has hit us nWait times.
00538         if(BLManager::getObject()->getState() == BEAM) {
00539                 if(--nWait > 0) return;
00540         } else if(BLManager::getObject()->getState() == VISUAL) {
00541                 if(--nWait > 0) return;
00542         } else {
00543                 if(--referenceWait > 0) return;
00544         }
00545 //printf("ParticleFilterPlacement-2 nWait=%d referenceWait=%d\n",nWait,referenceWait);
00546 
00547         // we get called for every particle transiting our volume, so we
00548         // need to check for presence in killSet and keepSet, then require.
00549         G4ParticleDefinition *def = track->GetDefinition();
00550         if(filter->killSet.count(def->GetParticleName()) > 0 ||
00551                         (filter->keepSet.size() > 0 &&
00552                          filter->keepSet.count(def->GetParticleName()) == 0)) {
00553                 track->SetTrackStatus(fStopAndKill);
00554                 ++nKilled;
00555                 if(filter->steppingVerbose) 
00556                         printf("ParticleFilterPlacement killed track\n");
00557         } else if(filter->require != "") {
00558                 static BLEvaluator *e=0;
00559                 if(!e)
00560                         e = new BLEvaluator();
00561                 G4ThreeVector position = track->GetPosition();
00562                 e->setVariable("x",position[0]);
00563                 e->setVariable("y",position[1]);
00564                 e->setVariable("z",position[2]);
00565                 e->setVariable("t",track->GetGlobalTime());
00566                 G4ThreeVector mom = track->GetMomentum();
00567                 e->setVariable("Px",mom[0]);
00568                 e->setVariable("Py",mom[1]);
00569                 e->setVariable("Pz",mom[2]);
00570                 e->setVariable("Ptot",sqrt(mom[0]*mom[0]+mom[1]*mom[1]+
00571                                                         mom[2]*mom[2]));
00572                 e->setVariable("PDGid",def->GetPDGEncoding());
00573                 e->setVariable("EventID",BLManager::getObject()->
00574                                 getEventID());
00575                 e->setVariable("TrackID",BLManager::getObject()->
00576                                 getExternalTrackID(track));
00577                 e->setVariable("ParentID",BLManager::getObject()->
00578                                 getExternalParentID(track));
00579                 if(e->evaluate(filter->require) == 0.0) {
00580                         track->SetTrackStatus(fStopAndKill);
00581                         ++nKilled;
00582                         if(filter->steppingVerbose) 
00583                                 printf("ParticleFilterPlacement killed track\n");
00584                 }
00585         }
00586 }

void ParticleFilterPlacement::PreUserTrackingAction ( const G4Track *  track  )  [virtual]

Implements BLManager::TrackingAction.

References filter, BLCMDparticlefilter::nWait, nWait, BLCMDparticlefilter::referenceWait, and referenceWait.

00589 {
00590         nWait = filter->nWait;
00591         referenceWait = filter->referenceWait;
00592 //printf("ParticleFilterPlacement-PreTracking nWait=%d referenceWait=%d\n",nWait,referenceWait);
00593 }

void ParticleFilterPlacement::PostUserTrackingAction ( const G4Track *  track  )  [inline, virtual]

Implements BLManager::TrackingAction.

00160 { }


Member Data Documentation

G4VPhysicalVolume* ParticleFilterPlacement::vol [private]

unsigned long ParticleFilterPlacement::nDecayed [private]

unsigned long ParticleFilterPlacement::nKilled [private]


The documentation for this class was generated from the following file:
g4beamline