Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
larg4::MCTruthEventActionService Class Reference

#include <MCTruthEventAction_service.h>

Inheritance diagram for larg4::MCTruthEventActionService:

Public Member Functions

 MCTruthEventActionService (fhicl::ParameterSet const &)
 
 ~MCTruthEventActionService ()
 
void setInputCollections (std::vector< art::Handle< std::vector< simb::MCTruth >>> const &mclists)
 

Private Member Functions

void generatePrimaries (G4Event *anEvent) override
 

Private Attributes

std::vector< art::Handle< std::vector< simb::MCTruth > > > const * fMCLists
 MCTruthCollection input lists. More...
 
std::map< G4int, G4int > fUnknownPDG
 map of unknown PDG codes to instances More...
 
std::map< G4int, G4int > fNon1StatusPDG
 PDG codes skipped because not status 1. More...
 
std::map< G4int, G4int > fProcessedPDG
 PDG codes processed. More...
 

Static Private Attributes

static G4ParticleTable * fParticleTable = nullptr
 Geant4's table of particle definitions. More...
 

Detailed Description

Definition at line 36 of file MCTruthEventAction_service.h.

Constructor & Destructor Documentation

larg4::MCTruthEventActionService::MCTruthEventActionService ( fhicl::ParameterSet const &  p)

Definition at line 36 of file MCTruthEventAction_service.cc.

37  : PrimaryGeneratorActionBase(p.get<string>("name", "MCTruthEventActionService"))
38 {}
p
Definition: test.py:223
larg4::MCTruthEventActionService::~MCTruthEventActionService ( )

Definition at line 40 of file MCTruthEventAction_service.cc.

41 {
42  // Print out a list of "unknown" PDG codes we saw in the input.
43  if (!fUnknownPDG.empty()) {
44  std::ostringstream badtxt;
45  for (auto const [pdg, n] : fUnknownPDG) {
46  badtxt << "\n Unknown PDG code = " << pdg << ", seen " << n << " times.";
47  if (genie_specific(pdg)) { badtxt << " (GENIE specific)"; }
48  }
49 
50  mf::LogWarning("MCTruthEventAction") << "The following unknown PDG codes were present in the "
51  << "simb::MCTruth input.\n"
52  << "They were not processed by Geant4." << badtxt.str();
53  }
54 
55  // Print out a list of PDG codes we saw in the input that weren't status=1
56  if (!fNon1StatusPDG.empty()) {
57  std::ostringstream non1txt;
58  for (auto const [pdg, n] : fNon1StatusPDG) {
59  non1txt << "\n PDG code = " << pdg << ", seen " << n << " times.";
60  if (genie_specific(pdg)) { non1txt << " (GENIE specific)"; }
61  }
62 
63  mf::LogDebug("MCTruthEventAction") << "The following PDG codes were present in the "
64  << "simb::MCTruth input with Status != 1.\n"
65  << "They were not processed by Geant4." << non1txt.str();
66  }
67 
68  // Print out a list of PDG codes that were processed
69  if (!fProcessedPDG.empty()) {
70  std::ostringstream goodtxt;
71  for (auto const [pdg, n] : fProcessedPDG) {
72  goodtxt << "\n PDG code = " << pdg << ", seen " << n << " times.";
73  }
74  mf::LogDebug("MCTruthEventAction") << "The following PDG codes were present in the "
75  << "simb::MCTruth input with Status = 1 "
76  << "and were processed by Geant4" << goodtxt.str();
77  }
78 }
std::void_t< T > n
std::map< G4int, G4int > fProcessedPDG
PDG codes processed.
std::map< G4int, G4int > fNon1StatusPDG
PDG codes skipped because not status 1.
std::map< G4int, G4int > fUnknownPDG
map of unknown PDG codes to instances
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning

Member Function Documentation

void larg4::MCTruthEventActionService::generatePrimaries ( G4Event *  anEvent)
overrideprivate

Definition at line 83 of file MCTruthEventAction_service.cc.

84 {
85  // For each MCTruth (probably only one, but you never know):
86  // index keeps track of which MCTruth object you are using
87  auto const& mclistHandles = *fMCLists;
88 
89  size_t mclSize = mclistHandles.size(); // -- should match the number of generators
90  mf::LogDebug("generatePrimaries") << "MCTruth Handles Size: " << mclSize;
91 
92  // -- Loop over MCTruth Handle List
93  size_t index = 0;
94  std::map<CLHEP::HepLorentzVector, G4PrimaryVertex*> vertexMap;
95  for (size_t mcl = 0; mcl < mclSize; ++mcl) {
96  mf::LogDebug("generatePrimaries")
97  << "MCTruth Handle Number: " << (mcl + 1) << " of " << mclSize;
98  art::Handle<std::vector<simb::MCTruth>> mclistHandle = mclistHandles[mcl];
99  // Loop over all MCTruth handle entries for a given generator,
100  // usually only one, but you never know
101  for (size_t i = 0; i < mclistHandle->size(); ++i) {
102  art::Ptr<simb::MCTruth> mclist(mclistHandle, i);
103 
104  mf::LogDebug("generatePrimaries")
105  << "art::Ptr number " << (i + 1) << " of " << mclistHandle->size() << ", Ptr: " << mclist;
106  int nPart = mclist->NParticles();
107  MF_LOG_INFO("generatePrimaries") << "Generating " << nPart << " particles";
108 
109  // -- Loop over all particles in MCTruth Object
110  for (int m = 0; m != nPart; ++m) {
111  simb::MCParticle const& particle = mclist->GetParticle(m);
112  G4int const pdgCode = particle.PdgCode();
113 
114  if (particle.StatusCode() != 1) {
115  MF_LOG_DEBUG("generatePrimaries")
116  << "Status code != 1, skipping particle number with MCTruth index = "
117  << "(" << mcl << "," << i << ")"
118  << " and particle index = " << m;
119  fNon1StatusPDG[pdgCode] += 1;
120  continue;
121  }
122 
123  if (((m + 1) % nPart) < 2) // -- only first and last will satisfy this
124  {
125  mf::LogDebug("generatePrimaries") << "Particle Number: " << (m + 1) << " of " << nPart;
126  mf::LogDebug("generatePrimaries") << "TrackID: " << particle.TrackId();
127  mf::LogDebug("generatePrimaries") << "index: " << index;
128  }
129 
130  G4double x = particle.Vx() * CLHEP::cm;
131  G4double y = particle.Vy() * CLHEP::cm;
132  G4double z = particle.Vz() * CLHEP::cm;
133  G4double t = particle.T() * CLHEP::ns;
134 
135  // Create a CLHEP four-vector from the particle's vertex.
136  CLHEP::HepLorentzVector fourpos(x, y, z, t);
137 
138  // Is this vertex already in our map?
139  G4PrimaryVertex* vertex = nullptr;
140  auto result = vertexMap.find(fourpos);
141  if (result == vertexMap.end()) {
142  // No, it's not, so create a new vertex and add it to the map.
143  vertex = new G4PrimaryVertex(x, y, z, t);
144  vertexMap[fourpos] = vertex;
145 
146  // Add the vertex to the G4Event.
147  anEvent->AddPrimaryVertex(vertex);
148  }
149  else {
150  // Yes, it is, so use the existing vertex.
151  vertex = result->second;
152  }
153 
154  // Get additional particle information.
155  TLorentzVector const& momentum = particle.Momentum(); // (px,py,pz,E)
156  TVector3 const& polarization = particle.Polarization();
157 
158  // Get the particle table if necessary. (Note: we're doing
159  // this "late" because I'm not sure at what point the G4
160  // particle table is initialized in the loading process.)
161  if (fParticleTable == nullptr) { fParticleTable = G4ParticleTable::GetParticleTable(); }
162 
163  // Get Geant4's definition of the particle.
164  G4ParticleDefinition* particleDefinition = nullptr;
165 
166  if (pdgCode == 0) { particleDefinition = fParticleTable->FindParticle("opticalphoton"); }
167  else {
168  particleDefinition = fParticleTable->FindParticle(pdgCode);
169  }
170 
171  if (pdgCode >= 2'000'000'000) { // If the particle is generator-specific
172  mf::LogDebug("ConvertPrimaryToGeant4")
173  << ": %%% Will skip particle with generator-specific PDG code = " << pdgCode;
174  fUnknownPDG[pdgCode] += 1;
175  continue;
176  }
177  if (pdgCode > 1'000'000'000) { // If the particle is a nucleus
178  mf::LogDebug("ConvertPrimaryToGeant4")
179  << ": %%% Nuclear PDG code = " << pdgCode << " (x,y,z,t)=(" << x << "," << y << "," << z
180  << "," << t << ")"
181  << " P=" << momentum.P() << ", E=" << momentum.E();
182  // If the particle table doesn't have a definition yet, ask
183  // the ion table for one. This will create a new ion
184  // definition as needed.
185  if (!particleDefinition) {
186  int const Z = (pdgCode % 10'000'000) / 10'000; // atomic number
187  int const A = (pdgCode % 10'000) / 10; // mass number
188  particleDefinition = fParticleTable->GetIonTable()->GetIon(Z, A, 0.);
189  }
190  }
191 
192  // What if the PDG code is unknown? This has been a known
193  // issue with GENIE.
194  if (particleDefinition == nullptr) {
195  mf::LogDebug("ConvertPrimaryToGeant4") << ": %%% Code not found = " << pdgCode;
196  fUnknownPDG[pdgCode] += 1;
197  continue;
198  }
199 
200  fProcessedPDG[pdgCode] += 1;
201 
202  // Create a Geant4 particle to add to the vertex.
203  auto* g4particle = new G4PrimaryParticle(particleDefinition,
204  momentum.Px() * CLHEP::GeV,
205  momentum.Py() * CLHEP::GeV,
206  momentum.Pz() * CLHEP::GeV);
207 
208  // Add more particle information the Geant4 particle.
209  g4particle->SetCharge(particleDefinition->GetPDGCharge());
210  g4particle->SetPolarization(polarization.x(), polarization.y(), polarization.z());
211 
212  // Add the particle to the vertex.
213  vertex->SetPrimary(g4particle);
214 
215  // Create a PrimaryParticleInformation object, and save the
216  // MCTruth pointer in it. This will allow the
217  // ParticleActionList class to access MCTruth information
218  // during Geant4's tracking.
219 
220  auto* primaryParticleInfo = new g4b::PrimaryParticleInformation;
221  primaryParticleInfo->SetMCTruth(mclist.get(), index, m);
222 
223  // Save the PrimaryParticleInformation in the
224  // G4PrimaryParticle for access during tracking.
225  g4particle->SetUserInformation(primaryParticleInfo);
226  } // -- for each particle in MCTruth
227  ++index;
228  } // -- for each MCTruth entry
229  } // -- for each MCTruth handle
230 } // -- generatePrimaries()
static constexpr double cm
Definition: Units.h:68
const TVector3 & Polarization() const
Definition: MCParticle.h:214
int PdgCode() const
Definition: MCParticle.h:212
static QCString result
static G4ParticleTable * fParticleTable
Geant4&#39;s table of particle definitions.
int StatusCode() const
Definition: MCParticle.h:211
int TrackId() const
Definition: MCParticle.h:210
std::vector< art::Handle< std::vector< simb::MCTruth > > > const * fMCLists
MCTruthCollection input lists.
static constexpr double GeV
Definition: Units.h:28
double T(const int i=0) const
Definition: MCParticle.h:224
std::map< G4int, G4int > fProcessedPDG
PDG codes processed.
#define MF_LOG_INFO(category)
std::map< G4int, G4int > fNon1StatusPDG
PDG codes skipped because not status 1.
double Vx(const int i=0) const
Definition: MCParticle.h:221
std::map< G4int, G4int > fUnknownPDG
map of unknown PDG codes to instances
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
#define MF_LOG_DEBUG(id)
double Vz(const int i=0) const
Definition: MCParticle.h:223
list x
Definition: train.py:276
def momentum(x1, x2, x3, scale=1.)
QAsciiDict< Entry > ns
double Vy(const int i=0) const
Definition: MCParticle.h:222
vertex reconstruction
void larg4::MCTruthEventActionService::setInputCollections ( std::vector< art::Handle< std::vector< simb::MCTruth >>> const &  mclists)
inline

Definition at line 42 of file MCTruthEventAction_service.h.

43  {
44  fMCLists = &mclists;
45  }
std::vector< art::Handle< std::vector< simb::MCTruth > > > const * fMCLists
MCTruthCollection input lists.

Member Data Documentation

std::vector<art::Handle<std::vector<simb::MCTruth> > > const* larg4::MCTruthEventActionService::fMCLists
private

MCTruthCollection input lists.

Definition at line 55 of file MCTruthEventAction_service.h.

std::map<G4int, G4int> larg4::MCTruthEventActionService::fNon1StatusPDG
private

PDG codes skipped because not status 1.

Definition at line 57 of file MCTruthEventAction_service.h.

G4ParticleTable * larg4::MCTruthEventActionService::fParticleTable = nullptr
staticprivate

Geant4's table of particle definitions.

Definition at line 53 of file MCTruthEventAction_service.h.

std::map<G4int, G4int> larg4::MCTruthEventActionService::fProcessedPDG
private

PDG codes processed.

Definition at line 58 of file MCTruthEventAction_service.h.

std::map<G4int, G4int> larg4::MCTruthEventActionService::fUnknownPDG
private

map of unknown PDG codes to instances

Definition at line 56 of file MCTruthEventAction_service.h.


The documentation for this class was generated from the following files: