ParticleListAction.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file ParticleListAction.cxx
3 /// \brief Use Geant4's user "hooks" to maintain a list of particles generated by Geant4.
4 ///
5 /// \author seligman@nevis.columbia.edu
6 ////////////////////////////////////////////////////////////////////////
7 
9 #include "lardataobj/Simulation/sim.h" // sim::NoParticleId
10 #include "nug4/G4Base/PrimaryParticleInformation.h"
11 #include "nug4/ParticleNavigation/ParticleList.h"
12 
13 #include "cetlib_except/exception.h"
15 
16 #include "Geant4/G4DynamicParticle.hh"
17 #include "Geant4/G4ParticleDefinition.hh"
18 #include "Geant4/G4PrimaryParticle.hh"
19 #include "Geant4/G4Step.hh"
20 #include "Geant4/G4StepPoint.hh"
21 #include "Geant4/G4String.hh"
22 #include "Geant4/G4ThreeVector.hh"
23 #include "Geant4/G4Track.hh"
24 #include "Geant4/G4VProcess.hh"
25 #include "Geant4/G4VUserPrimaryParticleInformation.hh"
26 
27 #include <TLorentzVector.h>
28 
29 #include <algorithm>
30 #include <cassert>
31 
32 //const G4bool debug = false; // unused
33 
34 // Photon variables defined at each step, for use
35 // in temporary velocity bug fix. -wforeman
37 bool entra = true;
38 
39 namespace larg4 {
40 
41  // Initialize static members.
45 
46  //----------------------------------------------------------------------------
47  // Dropped particle test
48 
49  bool
51  {
52  return !p || p->Trajectory().empty();
53  } // ParticleListAction::isDropped()
54 
55  //----------------------------------------------------------------------------
56  // Constructor.
58  bool storeTrajectories,
59  bool keepEMShowerDaughters,
60  bool keepMCParticleList /* = true */
61  )
62  : fenergyCut(energyCut * CLHEP::GeV)
63  , fparticleList(keepMCParticleList ? std::make_unique<sim::ParticleList>() : nullptr)
64  , fstoreTrajectories(storeTrajectories)
65  , fKeepEMShowerDaughters(keepEMShowerDaughters)
66  {}
67 
68  //----------------------------------------------------------------------------
69  // Begin the event
70  void
72  {
73  // Clear any previous particle information.
75  if (fparticleList) fparticleList->clear();
76  fParentIDMap.clear();
78  fCurrentPdgCode = 0;
79  }
80 
81  //-------------------------------------------------------------
82  // figure out the ultimate parentage of the particle with track ID
83  // trackid
84  // assume that the current track id has already been added to
85  // the fParentIDMap
86  int
88  {
89  int parentid = sim::NoParticleId;
90 
91  // search the fParentIDMap recursively until we have the parent id
92  // of the first EM particle that led to this one
94  while (itr != fParentIDMap.end()) {
95  MF_LOG_DEBUG("ParticleListAction") << "parentage for " << trackid << " " << (*itr).second;
96 
97  // set the parentid to the current parent ID, when the loop ends
98  // this id will be the first EM particle
99  parentid = (*itr).second;
100  itr = fParentIDMap.find(parentid);
101  }
102  MF_LOG_DEBUG("ParticleListAction") << "final parent ID " << parentid;
103 
104  return parentid;
105  }
106 
107  //----------------------------------------------------------------------------
108  // Create our initial simb::MCParticle object and add it to the sim::ParticleList.
109  void
111  {
112  // Particle type.
113  G4ParticleDefinition* particleDefinition = track->GetDefinition();
114  G4int pdgCode = particleDefinition->GetPDGEncoding();
115 
116  // Get Geant4's ID number for this track. This will be the same
117  // ID number that we'll use in the ParticleList.
118  // It is offset by the number of tracks accumulated from the previous Geant4
119  // runs (if any)
120  G4int trackID = track->GetTrackID() + fTrackIDOffset;
121  fCurrentTrackID = trackID;
122  fCurrentPdgCode = pdgCode;
123 
124  if (!fparticleList) {
125  // the rest is about adding a new particle to the list: skip
126  return; // note that fCurrentParticle is clear()'ed
127  }
128 
129  // And the particle's parent (same offset as above):
130  G4int parentID = track->GetParentID() + fTrackIDOffset;
131 
132  std::string process_name = "unknown";
133 
134  // Is there an MCTruth object associated with this G4Track? We
135  // have to go up a "chain" of information to find out:
136  const G4DynamicParticle* dynamicParticle = track->GetDynamicParticle();
137  const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
139  if (primaryParticle != 0) {
140  const G4VUserPrimaryParticleInformation* gppi = primaryParticle->GetUserInformation();
141  const g4b::PrimaryParticleInformation* ppi =
142  dynamic_cast<const g4b::PrimaryParticleInformation*>(gppi);
143  if (ppi != 0) {
144  primaryIndex = ppi->MCParticleIndex();
145 
146  // If we've made it this far, a PrimaryParticleInformation
147  // object exists and we are using a primary particle, set the
148  // process name accordingly
149  process_name = "primary";
150 
151  // primary particles should have parentID = 0, even if there
152  // are multiple MCTruths for this event
153  parentID = 0;
154  } // end else no primary particle information
155  } // Is there a G4PrimaryParticle?
156  // If this is not a primary particle...
157  else {
158  // check if this particle was made in an EM shower, don't put it in the particle
159  // list as we don't care about secondaries, tertiaries, etc for these showers
160  // figure out what process is making this track - skip it if it is
161  // one of pair production, compton scattering, photoelectric effect
162  // bremstrahlung, annihilation, any ionization - who wants to save
163  // a buttload of electrons that arent from a CC interaction?
164  process_name = track->GetCreatorProcess()->GetProcessName();
165  if (!fKeepEMShowerDaughters && (process_name.find("conv") != std::string::npos ||
166  process_name.find("LowEnConversion") != std::string::npos ||
167  process_name.find("Pair") != std::string::npos ||
168  process_name.find("compt") != std::string::npos ||
169  process_name.find("Compt") != std::string::npos ||
170  process_name.find("Brem") != std::string::npos ||
171  process_name.find("phot") != std::string::npos ||
172  process_name.find("Photo") != std::string::npos ||
173  process_name.find("Ion") != std::string::npos ||
174  process_name.find("annihil") != std::string::npos)) {
175 
176  // figure out the ultimate parentage of this particle
177  // first add this track id and its parent to the fParentIDMap
178  fParentIDMap[trackID] = parentID;
179 
180  fCurrentTrackID = -1 * this->GetParentage(trackID);
181 
182  // check that fCurrentTrackID is in the particle list - it is possible
183  // that this particle's parent is a particle that did not get tracked.
184  // An example is a partent that was made due to muMinusCaptureAtRest
185  // and the daughter was made by the phot process. The parent likely
186  // isn't saved in the particle list because it is below the energy cut
187  // which will put a bogus track id value into the sim::IDE object for
188  // the sim::SimChannel if we don't check it.
190 
191  // clear current particle as we are not stepping this particle and
192  // adding trajectory points to it
194  return;
195 
196  } // end if keeping EM shower daughters
197 
198  // Check the energy of the particle. If it falls below the energy
199  // cut, don't add it to our list.
200  G4double energy = track->GetKineticEnergy();
201  if (energy < fenergyCut) {
203 
204  // do add the particle to the parent id map though
205  // and set the current track id to be it's ultimate parent
206  fParentIDMap[trackID] = parentID;
207 
208  fCurrentTrackID = -1 * this->GetParentage(trackID);
209 
210  return;
211  }
212 
213  // check to see if the parent particle has been stored in the particle navigator
214  // if not, then see if it is possible to walk up the fParentIDMap to find the
215  // ultimate parent of this particle. Use that ID as the parent ID for this
216  // particle
217  if (!fparticleList->KnownParticle(parentID)) {
218  // do add the particle to the parent id map
219  // just in case it makes a daughter that we have to track as well
220  fParentIDMap[trackID] = parentID;
221  int pid = this->GetParentage(parentID);
222 
223  // if we still can't find the parent in the particle navigator,
224  // we have to give up
225  if (!fparticleList->KnownParticle(pid)) {
226  MF_LOG_WARNING("ParticleListAction")
227  << "can't find parent id: " << parentID << " in the particle list, or fParentIDMap."
228  << " Make " << parentID << " the mother ID for"
229  << " track ID " << fCurrentTrackID << " in the hope that it will aid debugging.";
230  }
231  else
232  parentID = pid;
233  }
234 
235  } // end if not a primary particle
236 
237  // This is probably the PDG mass, but just in case:
238  double mass = dynamicParticle->GetMass() / CLHEP::GeV;
239 
240  // Create the sim::Particle object.
243  new simb::MCParticle(trackID, pdgCode, process_name, parentID, mass);
244  fCurrentParticle.truthIndex = primaryIndex;
245 
246  // if we are not filtering, we have a decision already
247  if (!fFilter) fCurrentParticle.keep = true;
248 
249  // Polarization.
250  const G4ThreeVector& polarization = track->GetPolarization();
251  fCurrentParticle.particle->SetPolarization(
252  TVector3(polarization.x(), polarization.y(), polarization.z()));
253 
254  // Save the particle in the ParticleList.
256  }
257 
258  //----------------------------------------------------------------------------
259  void
261  {
262  if (!fCurrentParticle.hasParticle()) return;
263  assert(fparticleList);
264 
265  // if we have found no reason to keep it, drop it!
266  // (we might still need parentage information though)
267  if (!fCurrentParticle.keep) {
269  // after the particle is archived, it is deleted
271  return;
272  }
273 
274  if (aTrack) {
275  fCurrentParticle.particle->SetWeight(aTrack->GetWeight());
276  G4String process =
277  aTrack->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
278  fCurrentParticle.particle->SetEndProcess(process);
279  }
280 
281  // store truth record pointer, only if it is available
282  if (fCurrentParticle.isPrimary()) {
284  }
285 
286  return;
287  }
288 
289  //----------------------------------------------------------------------------
290  // With every step, add to the particle's trajectory.
291  void
293  {
294 
295  if (!fCurrentParticle.hasParticle()) { return; }
296 
297  // Temporary fix for problem where DeltaTime on the first step
298  // of optical photon propagation is calculated incorrectly. -wforeman
299  globalTime = step->GetTrack()->GetGlobalTime();
300  velocity_G4 = step->GetTrack()->GetVelocity();
301  velocity_step = step->GetStepLength() / step->GetDeltaTime();
302  if ((step->GetTrack()->GetDefinition()->GetPDGEncoding() == 0) &&
303  fabs(velocity_G4 - velocity_step) > 0.0001) {
304  // Subtract the faulty step time from the global time,
305  // and add the correct step time based on G4 velocity.
306  step->GetPostStepPoint()->SetGlobalTime(globalTime - step->GetDeltaTime() +
307  step->GetStepLength() / velocity_G4);
308  }
309 
310  // For the most part, we just want to add the post-step
311  // information to the particle's trajectory. There's one
312  // exception: In PreTrackingAction, the correct time information
313  // is not available. So add the correct vertex information here.
314 
315  if (fCurrentParticle.particle->NumberTrajectoryPoints() == 0) {
316 
317  // Get the pre/along-step information from the G4Step.
318  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
319 
320  const G4ThreeVector position = preStepPoint->GetPosition();
321  G4double time = preStepPoint->GetGlobalTime();
322 
323  // Remember that LArSoft uses cm, ns, GeV.
324  TLorentzVector fourPos(position.x() / CLHEP::cm,
325  position.y() / CLHEP::cm,
326  position.z() / CLHEP::cm,
327  time / CLHEP::ns);
328 
329  const G4ThreeVector momentum = preStepPoint->GetMomentum();
330  const G4double energy = preStepPoint->GetTotalEnergy();
331  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
332  momentum.y() / CLHEP::GeV,
333  momentum.z() / CLHEP::GeV,
334  energy / CLHEP::GeV);
335 
336  // Add the first point in the trajectory.
337  AddPointToCurrentParticle(fourPos, fourMom, "Start");
338 
339  } // end if this is the first step
340 
341  // At this point, the particle is being transported through the
342  // simulation. This method is being called for every voxel that
343  // the track passes through, but we don't want to update the
344  // trajectory information if we're just updating voxels. To check
345  // for this, look at the process name for the step, and compare it
346  // against the voxelization process name (set in PhysicsList.cxx).
347  G4String process = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
348  G4bool ignoreProcess = process.contains("LArVoxel") || process.contains("OpDetReadout");
349 
350  MF_LOG_DEBUG("ParticleListAction::SteppingAction")
351  << ": DEBUG - process='" << process << "'"
352  << " ignoreProcess=" << ignoreProcess << " fstoreTrajectories=" << fstoreTrajectories;
353 
354  // We store the initial creation point of the particle
355  // and its final position (ie where it has no more energy, or at least < 1 eV) no matter
356  // what, but whether we store the rest of the trajectory depends
357  // on the process, and on a user switch.
358  if (fstoreTrajectories && !ignoreProcess) {
359  // Get the post-step information from the G4Step.
360  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
361 
362  const G4ThreeVector position = postStepPoint->GetPosition();
363  G4double time = postStepPoint->GetGlobalTime();
364 
365  // Remember that LArSoft uses cm, ns, GeV.
366  TLorentzVector fourPos(position.x() / CLHEP::cm,
367  position.y() / CLHEP::cm,
368  position.z() / CLHEP::cm,
369  time / CLHEP::ns);
370 
371  const G4ThreeVector momentum = postStepPoint->GetMomentum();
372  const G4double energy = postStepPoint->GetTotalEnergy();
373  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
374  momentum.y() / CLHEP::GeV,
375  momentum.z() / CLHEP::GeV,
376  energy / CLHEP::GeV);
377 
378  // Add another point in the trajectory.
379  AddPointToCurrentParticle(fourPos, fourMom, std::string(process));
380  }
381  }
382 
383  //----------------------------------------------------------------------------
384  /// Utility class for the EndOfEventAction method: update the
385  /// daughter relationships in the particle list.
387  : public std::unary_function<sim::ParticleList::value_type, void> {
388  public:
389  UpdateDaughterInformation() : particleList(0) {}
390  void
391  SetParticleList(sim::ParticleList* p)
392  {
393  particleList = p;
394  }
395  void
396  operator()(sim::ParticleList::value_type& particleListEntry)
397  {
398  // We're looking at this Particle in the list.
399  int particleID = particleListEntry.first;
400 
401  // The parent ID of this particle;
402  // we ask the particle list since the particle itself might have been lost
403  // ("archived"), but the particle list still holds the information we need
404  int parentID = particleList->GetMotherOf(particleID);
405 
406  // If the parentID <= 0, this is a primary particle.
407  if (parentID <= 0) return;
408 
409  // If we get here, this particle is somebody's daughter. Add
410  // it to the list of daughter particles for that parent.
411 
412  // Get the parent particle from the list.
413  sim::ParticleList::iterator parentEntry = particleList->find(parentID);
414 
415  if (parentEntry == particleList->end()) {
416  // We have an "orphan": a particle whose parent isn't
417  // recorded in the particle list. This is not signficant;
418  // it's possible for a particle not to be saved in the list
419  // because it failed an energy cut, but for it to have a
420  // daughter that passed the cut (e.g., a nuclear decay).
421  return;
422  }
423  if (!parentEntry->second) return; // particle archived, nothing to update
424 
425  // Add the current particle to the daughter list of the
426  // parent.
427  simb::MCParticle* parent = (*parentEntry).second;
428  parent->AddDaughter(particleID);
429  }
430 
431  private:
432  sim::ParticleList* particleList;
433  };
434 
435  //----------------------------------------------------------------------------
436  // There's one last thing to do: All the particles have their
437  // parent IDs set (in PostTrackingAction), but we haven't set the
438  // daughters yet. That's done in this method.
439  void
441  {
442  if (!fparticleList) return;
443 
444  // Set up the utility class for the "for_each" algorithm. (We only
445  // need a separate set-up for the utility class because we need to
446  // give it the pointer to the particle list. We're using the STL
447  // "for_each" instead of the C++ "for loop" because it's supposed
448  // to be faster.
449  UpdateDaughterInformation updateDaughterInformation;
450  updateDaughterInformation.SetParticleList(fparticleList.get());
451 
452  // Update the daughter information for each particle in the list.
453  std::for_each(fparticleList->begin(), fparticleList->end(), updateDaughterInformation);
454  }
455 
456  //----------------------------------------------------------------------------
457  // Returns the ParticleList accumulated during the current event.
458  const sim::ParticleList*
460  {
461  if (!fparticleList) return nullptr;
462 
463  // check if the ParticleNavigator has entries, and if
464  // so grab the highest track id value from it to
465  // add to the fTrackIDOffset
466  int highestID = 0;
467  for (auto pn = fparticleList->begin(); pn != fparticleList->end(); pn++)
468  if ((*pn).first > highestID) highestID = (*pn).first;
469 
470  //Only change the fTrackIDOffset if there is in fact a particle to add to the event
471  if ((fparticleList->size()) != 0) { fTrackIDOffset = highestID + 1; }
472 
473  return fparticleList.get();
474  }
475 
476  //----------------------------------------------------------------------------
479  {
480  auto const iInfo = GetPrimaryTruthMap().find(trackId);
481  return (iInfo == GetPrimaryTruthMap().end()) ? simb::NoGeneratedParticleIndex : iInfo->second;
482  } // ParticleListAction::GetPrimaryTruthIndex()
483 
484  //----------------------------------------------------------------------------
485  // Yields the ParticleList accumulated during the current event.
486  sim::ParticleList&&
488  {
489  if (!fparticleList) {
490  // check with `hasList()` before calling this method
491  throw cet::exception("ParticleListAction")
492  << "ParticleListAction::YieldList(): particle list not build by user request.\n";
493  }
494  // check if the ParticleNavigator has entries, and if
495  // so grab the highest track id value from it to
496  // add to the fTrackIDOffset
497  int highestID = 0;
498  for (auto pn = fparticleList->begin(); pn != fparticleList->end(); pn++)
499  if ((*pn).first > highestID) highestID = (*pn).first;
500 
501  //Only change the fTrackIDOffset if there is in fact a particle to add to the event
502  if ((fparticleList->size()) != 0) { fTrackIDOffset = highestID + 1; }
503 
504  return std::move(*fparticleList);
505  } // ParticleList&& ParticleListAction::YieldList()
506 
507  //----------------------------------------------------------------------------
508  void
510  TLorentzVector const& mom,
511  std::string const& process)
512  {
513 
514  // Add the first point in the trajectory.
515  fCurrentParticle.particle->AddTrajectoryPoint(pos, mom, process);
516 
517  // also see if we can decide to keep the particle
518  if (!fCurrentParticle.keep) fCurrentParticle.keep = fFilter->mustKeep(pos);
519 
520  } // ParticleListAction::AddPointToCurrentParticle()
521 
522  //----------------------------------------------------------------------------
523 
524 } // namespace LArG4
static constexpr double cm
Definition: Units.h:68
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
intermediate_table::iterator iterator
int GetParentage(int trackid) const
G4bool fstoreTrajectories
Whether to store particle trajectories with each particle.
void AddDaughter(const int trackID)
Definition: MCParticle.h:268
void AddPointToCurrentParticle(TLorentzVector const &pos, TLorentzVector const &mom, std::string const &process)
Adds a trajectory point to the current particle, and runs the filter.
static int fCurrentPdgCode
pdg code of current particle
virtual void BeginOfEventAction(const G4Event *)
std::map< int, GeneratedParticleIndex_t > fPrimaryTruthMap
Map: particle track ID -> index of primary information in MC truth.
std::string string
Definition: nybbler.cc:12
const simb::MCTrajectory & Trajectory() const
Definition: MCParticle.h:253
GeneratedParticleIndex_t truthIndex
Index of the particle in the original generator truth record.
Geant4 interface.
double globalTime
STL namespace.
constexpr pointer get() const noexcept
Definition: exempt_ptr.h:148
const sim::ParticleList * GetList() const
intermediate_table::const_iterator const_iterator
void operator()(sim::ParticleList::value_type &particleListEntry)
GeneratedParticleIndex_t truthInfoIndex() const
Returns the index of the particle in the generator truth record.
bool empty() const
Definition: MCTrajectory.h:167
std::map< int, int > fParentIDMap
key is current track ID, value is parent ID
virtual void PreTrackingAction(const G4Track *)
def process(f, kind)
Definition: search.py:254
constexpr GeneratedParticleIndex_t NoGeneratedParticleIndex
Constant representing the absence of generator truth information.
Definition: simb.h:34
void clear()
Resets the information (does not release memory it does not own)
sim::ParticleList && YieldList()
static bool isDropped(simb::MCParticle const *p)
returns whether the specified particle has been marked as dropped
bool fKeepEMShowerDaughters
whether to keep EM shower secondaries, tertiaries, etc
static constexpr double GeV
Definition: Units.h:28
virtual void PostTrackingAction(const G4Track *)
static const int NoParticleId
Definition: sim.h:28
def move(depos, offset)
Definition: depos.py:107
ParticleListAction(double energyCut, bool storeTrajectories=false, bool keepEMShowerDaughters=false, bool keepMCParticleList=true)
p
Definition: test.py:223
void SetParticleList(sim::ParticleList *p)
double velocity_step
GeneratedParticleIndex_t GetPrimaryTruthIndex(int trackId) const
Returns the index of primary truth (sim::NoGeneratorIndex if none).
Code to link reconstructed objects back to the MC truth information.
std::unique_ptr< util::PositionInVolumeFilter > fFilter
filter for particles to be kept
std::unique_ptr< sim::ParticleList > fparticleList
bool isPrimary() const
Returns whether there is a particle.
cet::exempt_ptr< simb::MCParticle > particle
Object representing particle.
#define MF_LOG_DEBUG(id)
bool entra
double velocity_G4
virtual void EndOfEventAction(const G4Event *)
bool hasParticle() const
Returns whether there is a particle.
Tools and modules for checking out the basics of the Monte Carlo.
#define MF_LOG_WARNING(category)
def momentum(x1, x2, x3, scale=1.)
virtual void SteppingAction(const G4Step *)
QAsciiDict< Entry > ns
def parent(G, child, parent_type)
Definition: graph.py:67
std::map< int, GeneratedParticleIndex_t > const & GetPrimaryTruthMap() const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
std::size_t GeneratedParticleIndex_t
Type of particle index in the generator truth record (simb::MCTruth).
Definition: simb.h:30