MaterialFilter Class Reference

Inheritance diagram for MaterialFilter:

BLManager::SteppingAction BLCallback

List of all members.

Public Member Functions

 MaterialFilter (G4String keep, G4String kill, G4String _require, G4Material *_material)
void UserSteppingAction (const G4Step *step)
void callback (int type)

Private Attributes

std::set< G4ParticleDefinition * > keepSet
std::set< G4ParticleDefinition * > killSet
BLEvaluatoreval
G4String require
long nKilled
bool steppingVerbose
G4Material * material
G4ParticleTable * table


Constructor & Destructor Documentation

MaterialFilter::MaterialFilter ( G4String  keep,
G4String  kill,
G4String  _require,
G4Material *  _material 
)

References eval, BLParam::getInt(), BLManager::getObject(), keepSet, killSet, material, nKilled, Param, BLCommand::printError(), BLManager::registerCallback(), BLManager::registerSteppingAction(), require, BLCommand::splitString(), steppingVerbose, and table.

00271                                                :  BLManager::SteppingAction(),
00272                                         BLCallback(), keepSet(), killSet()
00273 {
00274         eval = 0;
00275         require = _require;
00276         nKilled = 0;
00277         steppingVerbose = Param.getInt("steppingVerbose");
00278         material = _material;
00279         table = G4ParticleTable::GetParticleTable();
00280 
00281         if(keep.size() != 0) {
00282                 printf("                keep=%s\n",keep.c_str());
00283                 std::vector<G4String> v=BLCommand::splitString(keep,',',true);
00284                 for(unsigned i=0; i<v.size(); ++i) {
00285                         G4ParticleDefinition *p = table->FindParticle(v[i]);
00286                         if(!p) {
00287                                 BLCommand::printError("material: particle '%s' not found",
00288                                         v[i].c_str());
00289                         }
00290                         keepSet.insert(p);
00291                 }
00292         }
00293 
00294         if(kill.size() != 0) {
00295                 printf("                kill=%s\n",kill.c_str());
00296                 std::vector<G4String> v=BLCommand::splitString(kill,',',true);
00297                 for(unsigned i=0; i<v.size(); ++i) {
00298                         G4ParticleDefinition *p = table->FindParticle(v[i]);
00299                         if(!p) {
00300                                 BLCommand::printError("material: particle '%s' not found",
00301                                         v[i].c_str());
00302                         }
00303                         killSet.insert(p);
00304                 }
00305         }
00306 
00307         if(require.size() != 0) {
00308                 printf("                require=%s\n",require.c_str());
00309                 eval = new BLEvaluator();
00310         }
00311 
00312         BLManager::getObject()->registerSteppingAction(this);
00313         BLManager::getObject()->registerCallback(this,2);
00314 }


Member Function Documentation

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

Implements BLManager::SteppingAction.

References BLCOORD_GLOBAL, eval, BLEvaluator::evaluate(), keepSet, killSet, material, nKilled, require, BLEvaluator::setTrackVariables(), and steppingVerbose.

00317 {
00318         // active in all states
00319 
00320         // only handle steps that enter a new volume
00321         G4VPhysicalVolume *preVol = 
00322                                 step->GetPreStepPoint()->GetPhysicalVolume();
00323         G4VPhysicalVolume *postVol = 
00324                                 step->GetPostStepPoint()->GetPhysicalVolume();
00325         if(preVol == postVol) return;
00326 
00327         // only handle this material (in the new volume)
00328         if(postVol->GetLogicalVolume()->GetMaterial() != material) return;
00329 
00330         G4Track *track = step->GetTrack();
00331         G4ParticleDefinition *def = track->GetDefinition();
00332         if(killSet.count(def) > 0 ||
00333                              (keepSet.size() > 0 && keepSet.count(def) == 0)) {
00334                 track->SetTrackStatus(fStopAndKill);
00335                 ++nKilled;
00336                 if(steppingVerbose) 
00337                         printf("material %s killed track\n",
00338                                   postVol->GetLogicalVolume()->GetMaterial()->
00339                                                         GetName().c_str());
00340         } else if(eval != 0) {
00341                 eval->setTrackVariables(track,BLCOORD_GLOBAL);
00342                 if(eval->evaluate(require) == 0.0) {
00343                         track->SetTrackStatus(fStopAndKill);
00344                         ++nKilled;
00345                         if(steppingVerbose) 
00346                                 printf("material %s killed track\n",
00347                                   postVol->GetLogicalVolume()->GetMaterial()->
00348                                                         GetName().c_str());
00349                 }
00350         }
00351 }

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

Reimplemented from BLCallback.

References material, and nKilled.

00354 {
00355         printf("material %s: %lu killed\n",material->GetName().c_str(),nKilled);
00356 }


Member Data Documentation

std::set<G4ParticleDefinition*> MaterialFilter::keepSet [private]

std::set<G4ParticleDefinition*> MaterialFilter::killSet [private]

G4String MaterialFilter::require [private]

long MaterialFilter::nKilled [private]

G4Material* MaterialFilter::material [private]

G4ParticleTable* MaterialFilter::table [private]

Referenced by MaterialFilter().


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