CMS 3D CMS Logo

CountProcessesAction.cc
Go to the documentation of this file.
2 
6 
7 #include "G4Event.hh"
8 #include "G4ParticleTable.hh"
9 #include "G4ProcessManager.hh"
10 #include "G4Run.hh"
11 #include "G4Step.hh"
12 #include "G4Track.hh"
13 
15  : fDEBUG(p.getUntrackedParameter<bool>("DEBUG", false)) {}
16 
18 
20  G4ParticleTable *partTable = G4ParticleTable::GetParticleTable();
21  int siz = partTable->size();
22  for (int ii = 0; ii < siz; ii++) {
23  G4ParticleDefinition *particle = partTable->GetParticle(ii);
24  std::string particleName = particle->GetParticleName();
25  if (fDEBUG)
26  std::cout << ii << " PCA " << particleName << " " << particle->GetPDGStable() << " " << particle->IsShortLived()
27  << std::endl;
29 
30  //--- All processes of this particle
31  G4ProcessManager *pmanager = particle->GetProcessManager();
32  G4ProcessVector *pvect = pmanager->GetProcessList();
33  int sizproc = pvect->size();
34  for (int jj = 0; jj < sizproc; jj++) {
35  std::string processName = (*pvect)[jj]->GetProcessName();
36  if (fDEBUG)
37  std::cout << jj << " PCR " << processName << std::endl;
39  }
40  }
41  DumpProcessList(false);
42 }
43 
45  //----- Fill counter of particles
46  const G4Track *aTrack = (*trk)();
47  std::string particleName = aTrack->GetDefinition()->GetParticleName();
49 
50  //----- Fill counter of Creator Processes
51  const G4VProcess *proc = aTrack->GetCreatorProcess();
53  if (proc != nullptr)
54  processName = proc->GetProcessName();
55  else
56  processName = "Primary";
57  pss parproc(particleName, processName);
58  mpssi::iterator ite = theCreatorProcessList.find(parproc);
59  if (ite == theCreatorProcessList.end())
60  theCreatorProcessList[parproc] = 1;
61  else
62  (*ite).second = (*ite).second + 1;
63  if (fDEBUG)
64  std::cout << " creator " << particleName << " " << processName << theCreatorProcessList.size() << std::endl;
65 }
66 
67 void CountProcessesAction::update(const G4Step *aStep) {
69  if (aStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr)
70  processName = aStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
71  else
72  processName = "User Limit";
73  std::string particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
75 }
76 
78  DumpProcessList(true);
81 }
82 
83 void CountProcessesAction::DumpProcessList(bool printNsteps, std::ostream &out) {
84  mpssi::iterator ite;
85  for (ite = theProcessList.begin(); ite != theProcessList.end(); ite++) {
86  if (!printNsteps)
87  out << "PROC_LIST " << (*ite).first.first << " : " << (*ite).first.second << std::endl;
88  else if ((*ite).second != 0)
89  out << "PROC_COUNT " << (*ite).first.first << " : " << (*ite).first.second << " = " << (*ite).second << std::endl;
90  }
91 }
92 
93 void CountProcessesAction::DumpCreatorProcessList(bool printNsteps, std::ostream &out) {
94  mpssi::iterator ite;
95  for (ite = theCreatorProcessList.begin(); ite != theCreatorProcessList.end(); ite++) {
96  if (!printNsteps)
97  out << "PROC-CREATOR_LIST " << (*ite).first.first << " : " << (*ite).first.second << std::endl;
98  else if ((*ite).second != 0)
99  out << "PROC_CREATOR_COUNT " << (*ite).first.first << " : " << (*ite).first.second << " = " << (*ite).second
100  << std::endl;
101  }
102 }
103 
105  psi::iterator ite;
106  for (ite = theParticleList.begin(); ite != theParticleList.end(); ite++) {
107  if ((*ite).second != 0)
108  out << "PART_LIST: " << (*ite).first << " = " << (*ite).second << std::endl;
109  }
110 }
void DumpProcessList(bool printNsteps, std::ostream &out=std::cout)
std::pair< ALIstring, ALIstring > pss
Definition: Fit.h:25
void DumpParticleList(std::ostream &out=std::cout)
void DumpCreatorProcessList(bool printNsteps, std::ostream &out=std::cout)
ii
Definition: cuy.py:589
CountProcessesAction(edm::ParameterSet const &p)
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.