CMS 3D CMS Logo

CountProcessesAction Class Reference

#include <SimG4Core/CountProcesses/interface/CountProcessesAction.h>

Inheritance diagram for CountProcessesAction:

SimWatcher Observer< const BeginOfRun * > Observer< const EndOfRun * > Observer< const BeginOfTrack * > Observer< const G4Step * >

List of all members.

Public Member Functions

 CountProcessesAction (edm::ParameterSet const &p)
void DumpCreatorProcessList (bool printNsteps, std::ostream &out=std::cout)
void DumpParticleList (std::ostream &out=std::cout)
void DumpProcessList (bool printNsteps, std::ostream &out=std::cout)
void update (const G4Step *track)
 This routine will be called when the appropriate signal arrives.
void update (const EndOfRun *track)
 This routine will be called when the appropriate signal arrives.
void update (const BeginOfTrack *trk)
 This routine will be called when the appropriate signal arrives.
void update (const BeginOfRun *run)
 This routine will be called when the appropriate signal arrives.
 ~CountProcessesAction ()

Private Attributes

bool fDEBUG
mpssi theCreatorProcessList
psi theParticleList
mpssi theProcessList


Detailed Description

Definition at line 23 of file CountProcessesAction.h.


Constructor & Destructor Documentation

CountProcessesAction::CountProcessesAction ( edm::ParameterSet const &  p  ) 

Definition at line 14 of file CountProcessesAction.cc.

00015     : fDEBUG(p.getUntrackedParameter<bool>("DEBUG",false))
00016 {}

CountProcessesAction::~CountProcessesAction (  ) 

Definition at line 18 of file CountProcessesAction.cc.

00018 {}


Member Function Documentation

void CountProcessesAction::DumpCreatorProcessList ( bool  printNsteps,
std::ostream &  out = std::cout 
)

Definition at line 100 of file CountProcessesAction.cc.

References lat::endl(), and theCreatorProcessList.

Referenced by update().

00101 {    
00102     mpssi::iterator ite;
00103     for (ite = theCreatorProcessList.begin(); ite != theCreatorProcessList.end(); ite++) 
00104     {
00105         if (!printNsteps) 
00106             out << "PROC-CREATOR_LIST " << (*ite).first.first << " : " 
00107                 <<(*ite) .first.second << std::endl; 
00108         else if ((*ite).second != 0) 
00109             out << "PROC_CREATOR_COUNT " << (*ite).first.first << " : " 
00110                 <<(*ite) .first.second << " = " << (*ite).second << std::endl; 
00111     }
00112 }

void CountProcessesAction::DumpParticleList ( std::ostream &  out = std::cout  ) 

Definition at line 114 of file CountProcessesAction.cc.

References lat::endl(), and theParticleList.

Referenced by update().

00115 {    
00116     psi::iterator ite;
00117     for (ite = theParticleList.begin(); ite != theParticleList.end(); ite++) 
00118     {
00119         if ((*ite).second != 0) 
00120             out << "PART_LIST: " << (*ite).first << " = " << (*ite).second << std::endl; 
00121     }
00122 }

void CountProcessesAction::DumpProcessList ( bool  printNsteps,
std::ostream &  out = std::cout 
)

Definition at line 86 of file CountProcessesAction.cc.

References lat::endl(), and theProcessList.

Referenced by update().

00087 {    
00088     mpssi::iterator ite;
00089     for (ite = theProcessList.begin(); ite != theProcessList.end(); ite++) 
00090     {
00091         if (!printNsteps) 
00092             out << "PROC_LIST " << (*ite).first.first << " : " 
00093                 << (*ite) .first.second << std::endl; 
00094         else if ((*ite).second != 0)
00095             out << "PROC_COUNT " << (*ite).first.first << " : " 
00096                 << (*ite) .first.second << " = " << (*ite).second << std::endl; 
00097     }
00098 }

void CountProcessesAction::update ( const G4Step *   )  [virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step * >.

Definition at line 69 of file CountProcessesAction.cc.

References theProcessList.

00070 {
00071     std::string processName;
00072     if(aStep->GetPostStepPoint()->GetProcessDefinedStep() != 0)
00073         processName = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
00074     else processName = "User Limit";
00075     std::string particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
00076     theProcessList[pss(particleName,processName)] = theProcessList[pss(particleName,processName)] + 1;
00077 }

void CountProcessesAction::update ( const EndOfRun  )  [virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfRun * >.

Definition at line 79 of file CountProcessesAction.cc.

References DumpCreatorProcessList(), DumpParticleList(), and DumpProcessList().

00080 {
00081     DumpProcessList(1);
00082     DumpCreatorProcessList(1);
00083     DumpParticleList();
00084 }

void CountProcessesAction::update ( const BeginOfTrack  )  [virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack * >.

Definition at line 48 of file CountProcessesAction.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), fDEBUG, if(), proc, theCreatorProcessList, and theParticleList.

00049 {
00050     //----- Fill counter of particles
00051     const G4Track * aTrack = (*trk)();
00052     std::string particleName = aTrack->GetDefinition()->GetParticleName();
00053     theParticleList[particleName]++;
00054 
00055     //----- Fill counter of Creator Processes
00056     const G4VProcess * proc = aTrack->GetCreatorProcess();
00057     std::string processName;
00058     if (proc != 0) processName = proc->GetProcessName();
00059     else processName = "Primary";
00060     pss parproc(particleName,processName);
00061     mpssi::iterator ite = theCreatorProcessList.find(parproc);
00062     if (ite == theCreatorProcessList.end()) theCreatorProcessList[ parproc ] = 1;
00063     else (*ite).second = (*ite).second +1; 
00064     if (fDEBUG) 
00065         std::cout << " creator " << particleName << " " << processName 
00066                   << theCreatorProcessList.size() << std::endl;
00067 }

void CountProcessesAction::update ( const BeginOfRun  )  [virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun * >.

Definition at line 20 of file CountProcessesAction.cc.

References GenMuonPlsPt100GeV_cfg::cout, DumpProcessList(), lat::endl(), fDEBUG, siz, theParticleList, and theProcessList.

00021 {
00022     G4ParticleTable * partTable = G4ParticleTable::GetParticleTable();
00023     int siz = partTable->size();
00024     for (int ii= 0; ii < siz; ii++)
00025     {
00026         G4ParticleDefinition * particle = partTable->GetParticle(ii);
00027         std::string particleName = particle->GetParticleName();
00028         if (fDEBUG)
00029             std::cout << ii << " PCA " << particleName<< " " << particle->GetPDGStable() 
00030                       << " " << particle->IsShortLived() << std::endl;
00031         theParticleList[particleName] = 0;
00032 
00033         //--- All processes of this particle 
00034         G4ProcessManager * pmanager = particle->GetProcessManager();
00035         G4ProcessVector * pvect = pmanager->GetProcessList();
00036         int sizproc = pvect->size();
00037         for (int jj = 0; jj < sizproc; jj++) 
00038         {
00039             std::string processName = (*pvect)[jj]->GetProcessName();
00040             if (fDEBUG)
00041                 std::cout << jj << " PCR " << processName<< std::endl;
00042             theProcessList[pss(particleName,processName)] = 0;
00043         }
00044     }
00045     DumpProcessList(0);
00046 }


Member Data Documentation

bool CountProcessesAction::fDEBUG [private]

Definition at line 43 of file CountProcessesAction.h.

Referenced by update().

mpssi CountProcessesAction::theCreatorProcessList [private]

Definition at line 45 of file CountProcessesAction.h.

Referenced by DumpCreatorProcessList(), and update().

psi CountProcessesAction::theParticleList [private]

Definition at line 46 of file CountProcessesAction.h.

Referenced by DumpParticleList(), and update().

mpssi CountProcessesAction::theProcessList [private]

Definition at line 44 of file CountProcessesAction.h.

Referenced by DumpProcessList(), and update().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:16:58 2009 for CMSSW by  doxygen 1.5.4