Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:19

0001 #include "SimG4Core/CountProcesses/interface/CountProcessesAction.h"
0002 
0003 #include "SimG4Core/Notification/interface/BeginOfRun.h"
0004 #include "SimG4Core/Notification/interface/BeginOfTrack.h"
0005 #include "SimG4Core/Notification/interface/EndOfRun.h"
0006 
0007 #include "G4Event.hh"
0008 #include "G4ParticleTable.hh"
0009 #include "G4ProcessManager.hh"
0010 #include "G4Run.hh"
0011 #include "G4Step.hh"
0012 #include "G4Track.hh"
0013 
0014 CountProcessesAction::CountProcessesAction(edm::ParameterSet const &p)
0015     : fDEBUG(p.getUntrackedParameter<bool>("DEBUG", false)) {}
0016 
0017 CountProcessesAction::~CountProcessesAction() {}
0018 
0019 void CountProcessesAction::update(const BeginOfRun *run) {
0020   G4ParticleTable *partTable = G4ParticleTable::GetParticleTable();
0021   int siz = partTable->size();
0022   for (int ii = 0; ii < siz; ii++) {
0023     G4ParticleDefinition *particle = partTable->GetParticle(ii);
0024     std::string particleName = particle->GetParticleName();
0025     if (fDEBUG)
0026       std::cout << ii << " PCA " << particleName << " " << particle->GetPDGStable() << " " << particle->IsShortLived()
0027                 << std::endl;
0028     theParticleList[particleName] = 0;
0029 
0030     //--- All processes of this particle
0031     G4ProcessManager *pmanager = particle->GetProcessManager();
0032     G4ProcessVector *pvect = pmanager->GetProcessList();
0033     int sizproc = pvect->size();
0034     for (int jj = 0; jj < sizproc; jj++) {
0035       std::string processName = (*pvect)[jj]->GetProcessName();
0036       if (fDEBUG)
0037         std::cout << jj << " PCR " << processName << std::endl;
0038       theProcessList[pss(particleName, processName)] = 0;
0039     }
0040   }
0041   DumpProcessList(false);
0042 }
0043 
0044 void CountProcessesAction::update(const BeginOfTrack *trk) {
0045   //----- Fill counter of particles
0046   const G4Track *aTrack = (*trk)();
0047   std::string particleName = aTrack->GetDefinition()->GetParticleName();
0048   theParticleList[particleName]++;
0049 
0050   //----- Fill counter of Creator Processes
0051   const G4VProcess *proc = aTrack->GetCreatorProcess();
0052   std::string processName;
0053   if (proc != nullptr)
0054     processName = proc->GetProcessName();
0055   else
0056     processName = "Primary";
0057   pss parproc(particleName, processName);
0058   mpssi::iterator ite = theCreatorProcessList.find(parproc);
0059   if (ite == theCreatorProcessList.end())
0060     theCreatorProcessList[parproc] = 1;
0061   else
0062     (*ite).second = (*ite).second + 1;
0063   if (fDEBUG)
0064     std::cout << " creator " << particleName << " " << processName << theCreatorProcessList.size() << std::endl;
0065 }
0066 
0067 void CountProcessesAction::update(const G4Step *aStep) {
0068   std::string processName;
0069   if (aStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr)
0070     processName = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
0071   else
0072     processName = "User Limit";
0073   std::string particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
0074   theProcessList[pss(particleName, processName)] = theProcessList[pss(particleName, processName)] + 1;
0075 }
0076 
0077 void CountProcessesAction::update(const EndOfRun *run) {
0078   DumpProcessList(true);
0079   DumpCreatorProcessList(true);
0080   DumpParticleList();
0081 }
0082 
0083 void CountProcessesAction::DumpProcessList(bool printNsteps, std::ostream &out) {
0084   mpssi::iterator ite;
0085   for (ite = theProcessList.begin(); ite != theProcessList.end(); ite++) {
0086     if (!printNsteps)
0087       out << "PROC_LIST " << (*ite).first.first << " : " << (*ite).first.second << std::endl;
0088     else if ((*ite).second != 0)
0089       out << "PROC_COUNT " << (*ite).first.first << " : " << (*ite).first.second << " = " << (*ite).second << std::endl;
0090   }
0091 }
0092 
0093 void CountProcessesAction::DumpCreatorProcessList(bool printNsteps, std::ostream &out) {
0094   mpssi::iterator ite;
0095   for (ite = theCreatorProcessList.begin(); ite != theCreatorProcessList.end(); ite++) {
0096     if (!printNsteps)
0097       out << "PROC-CREATOR_LIST " << (*ite).first.first << " : " << (*ite).first.second << std::endl;
0098     else if ((*ite).second != 0)
0099       out << "PROC_CREATOR_COUNT " << (*ite).first.first << " : " << (*ite).first.second << " = " << (*ite).second
0100           << std::endl;
0101   }
0102 }
0103 
0104 void CountProcessesAction::DumpParticleList(std::ostream &out) {
0105   psi::iterator ite;
0106   for (ite = theParticleList.begin(); ite != theParticleList.end(); ite++) {
0107     if ((*ite).second != 0)
0108       out << "PART_LIST: " << (*ite).first << " = " << (*ite).second << std::endl;
0109   }
0110 }