LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ParticleListAction.cc
Go to the documentation of this file.
1 
30 
32 
34 
35 // Framework includes
36 #include "fhiclcpp/ParameterSet.h"
37 
38 // ROOT includes
39 #include "TLorentzVector.h"
40 
41 #include "CLHEP/Units/SystemOfUnits.h"
42 
43 // Geant4 includes
44 #include "Geant4/G4DynamicParticle.hh"
45 #include "Geant4/G4ParticleDefinition.hh"
46 #include "Geant4/G4PrimaryParticle.hh"
47 #include "Geant4/G4StepPoint.hh"
48 #include "Geant4/G4ThreeVector.hh"
49 #include "Geant4/G4Track.hh"
50 #include "Geant4/G4VProcess.hh"
51 #include "Geant4/G4VUserPrimaryParticleInformation.hh"
52 
53 #include "range/v3/view.hpp"
54 
55 // STL includes
56 #include <algorithm>
57 #include <cassert>
58 #include <sstream>
59 #include <string>
60 #include <vector>
61 
62 namespace larg4 {
63 
64  //----------------------------------------------------------------------------
65  // Constructor.
67  : artg4tk::EventActionBase("PLASEventActionBase")
68  , artg4tk::TrackingActionBase("PLASTrackingActionBase")
69  , artg4tk::SteppingActionBase("PLASSteppingActionBase")
70  , fenergyCut(p.get<double>("EnergyCut", 0.0 * CLHEP::GeV))
71  , fstoreTrajectories(p.get<bool>("storeTrajectories", true))
72  , fkeepGenTrajectories(p.get<std::vector<std::string>>("keepGenTrajectories", {}))
73  , fKeepEMShowerDaughters(p.get<bool>("keepEMShowerDaughters", true))
74  , fNotStoredPhysics(p.get<std::vector<std::string>>("NotStoredPhysics", {}))
75  , fkeepOnlyPrimaryFullTraj(p.get<bool>("keepOnlyPrimaryFullTrajectories", false))
76  , fSparsifyTrajectories(p.get<bool>("SparsifyTrajectories", false))
77  , fSparsifyMargin(p.get<double>("SparsifyMargin", 0.015))
78  , fKeepTransportation(p.get<bool>("KeepTransportation", false))
79  , fKeepSecondToLast(p.get<bool>("KeepSecondToLast", false))
80  , fKeepParticlesInVolumes(p.get<std::vector<std::string>>("KeepParticlesInVolumes", {}))
82  p.get<std::vector<std::string>>("KeepDroppedParticlesInVolumes", {}))
85  !fKeepDroppedParticlesInVolumes.empty() ? std::make_unique<sim::ParticleList>() : nullptr)
86  {
87  //Assert that keepEMShowerDaughters and storeDroppedMCParticles are not both true
90  << "ParticleListActionService: keepEMShowerDaughters and storeDroppedMCParticles cannot "
91  "both be true.\n";
92  }
93  // -- D.R. If a custom list of not storable physics is provided, use it, otherwise
94  // use the default list. This preserves the behavior of the keepEmShowerDaughters
95  // parameter
96  bool customNotStored = not fNotStoredPhysics.empty();
98  // -- Don't keep all processes
99  if (!customNotStored) // -- Don't keep but haven't provided a list
100  { // -- default list of not stored physics
101  fNotStoredPhysics = {"conv",
102  "LowEnConversion",
103  "Pair",
104  "compt",
105  "Compt",
106  "Brem",
107  "phot",
108  "Photo",
109  "Ion",
110  "annihil"};
111  }
112 
113  std::stringstream sstored;
114  sstored << "The full tracking information will not be stored for particles"
115  << " resulting from the following processes: \n{ ";
116  for (auto const& i : fNotStoredPhysics) {
117  sstored << "\"" << i << "\" ";
118  fNotStoredCounterUMap.emplace(i, 0); // -- initialize counter
119  }
120  mf::LogInfo("ParticleListActionService") << sstored.str() << "}\n";
121  }
122  else { // -- Keep all processes
123  mf::LogInfo("ParticleListActionService")
124  << "Storing full tracking information for all processes. \n";
125  if (customNotStored) // -- custom list will be ignored
126  {
127  mf::LogWarning("StoredPhysics")
128  << "NotStoredPhysics provided, but will be ignored."
129  << " To use NotStoredPhysics, set keepEMShowerDaughters to false";
130  }
131  }
132 
133  // -- sparsify info
135  mf::LogInfo("ParticleListActionService")
136  << "Trajectory sparsification enabled with SparsifyMargin : " << fSparsifyMargin << "\n";
137  } // end constructor
138 
139  //----------------------------------------------------------------------------
140  // Begin the event
142  {
143  // Clear any previous particle information.
146  fParentIDMap.clear();
147  fTargetIDMap.clear();
148  fMCTIndexMap.clear();
149  fMCTPrimProcessKeepMap.clear();
151  fTrackIDOffset = 0;
152  fPrimaryTruthMap.clear();
153  fMCTIndexToGeneratorMap.clear();
154  fNotStoredCounterUMap.clear();
155  fdroppedTracksMap.clear();
157  // -- D.R. If a custom list of keepGenTrajectories is provided, use it, otherwise
158  // keep or drop decision made based storeTrajectories parameter. This preserves
159  // the behavior of the storeTrajectories fhicl param
160  bool customKeepTraj = not fkeepGenTrajectories.empty();
161  if (!fstoreTrajectories) { // -- fstoreTrajectories : false
162  mf::LogDebug("beginOfEventAction::Generator") << "Trajectory points will not be stored.";
163  }
164  else if (!customKeepTraj) { // -- fstoretrajectories : true and empty keepGenTrajectories list
165  mf::LogDebug("beginOfEventAction::Generator")
166  << "keepGenTrajectories list is empty. Will"
167  << " store trajectory points for all generators";
168  }
169 
170  // -- D.R. determine mapping between MCTruthIndex(s) and generator(s) for later reference
171  size_t nKeep = 0;
172  std::string generator_name = "unknown";
173  for (size_t mcti = 0; mcti < fMCLists->size(); ++mcti) {
174 
175  std::stringstream sskeepgen;
176  sskeepgen << "MCTruth object summary :";
177  sskeepgen << "\n\tPrimary MCTIndex : " << mcti;
178 
179  // -- Obtain the generator (provenance) corresponding to the mctruth index:
180  auto const& mclistHandle = (*fMCLists)[mcti];
181  generator_name = mclistHandle.provenance()->inputTag().label();
182  sskeepgen << "\n\tProvenance/Generator : " << generator_name;
183 
184  G4bool keepGen = false;
185  if (fstoreTrajectories) // -- storeTrajectories set to true; check which
186  {
187  if (!customKeepTraj) { // -- no custom list, so keep all
188  keepGen = true;
189  nKeep++;
190  }
191  else { // -- custom list, so check the ones in the event against provided keep list
192  for (auto keepableGen : fkeepGenTrajectories) {
193  if (generator_name == keepableGen) { // -- exit upon finding match; false by default
194  keepGen = true;
195  nKeep++;
196  break;
197  }
198  }
199  }
200  }
201  fMCTIndexToGeneratorMap.emplace(mcti, std::make_pair(generator_name, keepGen));
202  sskeepgen << "\n\tTrajectory points storable : " << (keepGen ? "true" : "false") << "\n";
203  mf::LogDebug("beginOfEventAction::Generator") << sskeepgen.str();
204  }
205 
206  if (nKeep == 0 && customKeepTraj && fstoreTrajectories) {
207  mf::LogWarning("beginOfEventAction::keepableGenerators")
208  << "storeTrajectories"
209  << " set to true and a non-empty keepGenTrajectories list provided in configuration file, "
210  "but"
211  << " none of the generators in this list are present in the event! Double check list or "
212  "don't"
213  << " provide keepGenTrajectories in the configuration to keep all trajectories from all"
214  << " generator labels. This may be expected for generators that have a nonzero probability "
215  "of"
216  << " producing no particles (e.g. some radiologicals)";
217  }
218  }
219 
220  //-------------------------------------------------------------
221  // figure out the ultimate parentage of the particle with track ID
222  // trackid
223  // assume that the current track id has already been added to
224  // the fParentIDMap
226  {
227  int parentid = sim::NoParticleId;
228 
229  // search the fParentIDMap recursively until we have the parent id
230  // of the first EM particle that led to this one
231  auto itr = fParentIDMap.find(trackid);
232  while (itr != fParentIDMap.end()) {
233  // set the parentid to the current parent ID, when the loop ends
234  // this id will be the first EM particle
235  parentid = (*itr).second;
236  itr = fParentIDMap.find(parentid);
237  }
238 
239  return parentid;
240  }
241 
242  //----------------------------------------------------------------------------
243  // Create our initial simb::MCParticle object and add it to the sim::ParticleList.
245  {
246  // Particle type.
247  G4ParticleDefinition* particleDefinition = track->GetDefinition();
248  G4int pdgCode = particleDefinition->GetPDGEncoding();
249 
250  // Get Geant4's ID number for this track. This will be the same
251  // ID number that we'll use in the ParticleList.
252  // It is offset by the number of tracks accumulated from the previous Geant4
253  // runs (if any)
254  int const trackID = track->GetTrackID() + fTrackIDOffset;
255  fCurrentTrackID = trackID;
256  fTargetIDMap[trackID] = fCurrentTrackID;
257  // And the particle's parent (same offset as above):
258  int parentID = track->GetParentID() + fTrackIDOffset;
259 
260  std::string process_name = "unknown";
261  std::string mct_primary_process = "unknown";
262  bool isFromMCTProcessPrimary = false;
263  bool notstore = false;
264 
265  // Is there an MCTruth object associated with this G4Track? We
266  // have to go up a "chain" of information to find out:
267  const G4DynamicParticle* dynamicParticle = track->GetDynamicParticle();
268  const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
270  size_t primarymctIndex = 0;
271  if (primaryParticle != nullptr) {
272  const G4VUserPrimaryParticleInformation* gppi = primaryParticle->GetUserInformation();
274  dynamic_cast<const g4b::PrimaryParticleInformation*>(gppi);
275  if (ppi != nullptr) {
276  primaryIndex = ppi->MCParticleIndex();
277  primarymctIndex = ppi->MCTruthIndex();
278  mct_primary_process = ppi->GetMCParticle()->Process();
279 
280  // If we've made it this far, a PrimaryParticleInformation
281  // object exists and we are using a primary particle, set the
282  // process name accordingly
283 
284  // -- D.R. : if process == "primary" exactly (this will most likely be the case), mark it as
285  // a primary with keepable trajectory for itself and its descendants.
286  // -- elsif it simply starts with "primary", accept it but mark it as non-keep
287  // -- else force it to be primary because we have determined that it is a primary particle
288  // This was the original behavior of the code i.e. process was _set_ to "primary"
289  // regardless of what the process name was in the generator MCTruth object.
290  //
291  // -- NOTE: This enforces a convention for process names assigned in the gen stage.
292  if (mct_primary_process.compare("primary") == 0) {
293  process_name = "primary";
294  isFromMCTProcessPrimary = true;
295  }
296  else if (mct_primary_process.find("primary") == 0) {
297  process_name = mct_primary_process;
298  isFromMCTProcessPrimary = false;
299  mf::LogDebug("PrimaryParticle") << "MCTruth primary process name contains \"primary\" "
300  << " but is not solely \"primary\" : " << process_name
301  << ".\nWill not store full set of trajectory points.";
302  }
303  else { // -- override it
304  process_name = "primary";
305  isFromMCTProcessPrimary = true;
306  mf::LogWarning("PrimaryParticle")
307  << "MCTruth primary process does not beging with string"
308  << " literal \"primary\" : " << process_name << "\nOVERRIDING it to \"primary\"";
309  }
310  // -- The process_name check is simply to allow an additional way to reduce memory usage,
311  // namely, by creating MCTruth input particles with multiple process labels (e.g. "primary"
312  // and "primaryBackground") that can be used to veto storage of trajectory points.
313 
314  // primary particles should have parentID = 0, even if there
315  // are multiple MCTruths for this event
316  parentID = 0;
317  } // end else no primary particle information
318  } // Is there a G4PrimaryParticle?
319  // If this is not a primary particle...
320  else {
321  // check if this particle was made in an undesirable process. For example:
322  // if one is not interested in EM shower particles, don't put it in the particle
323  // list as one wouldn't care about secondaries, tertiaries, etc. For these showers
324  // figure out what process is making this track - skip it if it is
325  // one of pair production, compton scattering, photoelectric effect
326  // bremstrahlung, annihilation, or ionization
327  process_name = track->GetCreatorProcess()->GetProcessName();
328  if (!fKeepEMShowerDaughters) {
329  for (auto const& p : fNotStoredPhysics) {
330  if (process_name.find(p) != std::string::npos) {
331  notstore = true;
332  mf::LogDebug("NotStoredPhysics") << "Found process : " << process_name;
333 
334  int old = 0;
335  auto search = fNotStoredCounterUMap.find(p);
336  if (search != fNotStoredCounterUMap.end()) { old = search->second; }
337  fNotStoredCounterUMap.insert_or_assign(p, (old + 1));
338  break;
339  }
340  }
341  if (notstore) {
342  // figure out the ultimate parentage of this particle
343  // first add this track id and its parent to the fParentIDMap
344  fParentIDMap[trackID] = parentID;
345  fCurrentTrackID = -1 * this->GetParentage(trackID);
346  // check that fCurrentTrackID is in the particle list - it is possible
347  // that this particle's parent is a particle that did not get tracked.
348  // An example is a parent that was made due to muMinusCaptureAtRest
349  // and the daughter was made by the phot process. The parent likely
350  // isn't saved in the particle list because it is below the energy cut
351  // which will put a bogus track id value into the sim::IDE object for
352  // the sim::SimChannel if we don't check it.
354  fTargetIDMap[trackID] = fCurrentTrackID;
355  // clear current particle as we are not stepping this particle and
356  // adding trajectory points to it
357  fdroppedTracksMap[this->GetParentage(trackID)].insert(trackID);
358  // keep track of this particle in the fMCTIndexMap as well, as we may keep a daughter
359  if (auto it = fMCTIndexMap.find(parentID); it != cend(fMCTIndexMap)) {
360  fMCTIndexMap[trackID] = it->second;
361  }
362  // keep track of this particle in the fMCTIndexMap as well, as we may keep a daughter
363  if (auto it = fMCTIndexMap.find(parentID); it != cend(fMCTIndexMap)) {
364  fMCTIndexMap[trackID] = it->second;
365  }
366  if (!fStoreDroppedMCParticles) { //Only clear if not storing dropped particles
368  return;
369  }
370  } // end if process matches an undesired process
371  } // end if not keeping EM shower daughters
372 
373  // Check the energy of the particle. If it falls below the energy
374  // cut, don't add it to our list.
375  G4double energy = track->GetKineticEnergy();
376  if (energy < fenergyCut && pdgCode != 0) {
377  fdroppedTracksMap[this->GetParentage(trackID)].insert(trackID);
379  // do add the particle to the parent id map though
380  // and set the current track id to be it's ultimate parent
381  fParentIDMap[trackID] = parentID;
382  fCurrentTrackID = -1 * this->GetParentage(trackID);
383  fTargetIDMap[trackID] = fCurrentTrackID;
384  // keep track of this particle in the fMCTIndexMap as well, as we may keep a daughter
385  if (auto it = fMCTIndexMap.find(parentID); it != cend(fMCTIndexMap)) {
386  fMCTIndexMap[trackID] = it->second;
387  }
388  return;
389  }
390 
391  // check to see if the parent particle has been stored in the particle navigator
392  // if not, then see if it is possible to walk up the fParentIDMap to find the
393  // ultimate parent of this particle. Use that ID as the parent ID for this
394  // particle
395  if (!fParticleList.KnownParticle(parentID) &&
396  (fMCTIndexMap.count(parentID) == 0 ||
397  !(fdroppedParticleList && fdroppedParticleList->KnownParticle(parentID)))) {
398  // do add the particle to the parent id map
399  // just in case it makes a daughter that we have to track as well
400  fParentIDMap[trackID] = parentID;
401  int pid = this->GetParentage(parentID);
402 
403  // if we still can't find the parent in the particle navigator,
404  // we have to give up
405  if (!fParticleList.KnownParticle(pid) &&
406  (fMCTIndexMap.count(pid) == 0 ||
407  !(fdroppedParticleList && fdroppedParticleList->KnownParticle(parentID)))) {
408  MF_LOG_WARNING("ParticleListActionService")
409  << "can't find parent id: " << parentID << " in the particle list, or fParentIDMap."
410  << " Make " << parentID << " the mother ID for"
411  << " track ID " << fCurrentTrackID << " in the hope that it will aid debugging.";
412  }
413  else
414  parentID = pid;
415  }
416 
417  // Once the parentID is secured, inherit the MCTruth Index
418  // which should have been set already
419  if (auto it = fMCTIndexMap.find(parentID); it != cend(fMCTIndexMap)) {
420  primarymctIndex = it->second;
421  }
422  else {
424  << "Could not locate MCT Index for parent trackID of " << parentID;
425  }
426 
427  // Inherit whether the parent is from a primary with MCTruth process_name == "primary"
428  isFromMCTProcessPrimary = fMCTPrimProcessKeepMap[parentID];
429  } // end if not a primary particle
430 
431  // This is probably the PDG mass, but just in case:
432  double mass = dynamicParticle->GetMass() / CLHEP::GeV;
433 
434  // Create the sim::Particle object.
436  fCurrentParticle.isDropped = notstore; //mark if the particle is dropped
438  new simb::MCParticle{trackID, pdgCode, process_name, parentID, mass};
439  fCurrentParticle.truthIndex = primaryIndex;
440 
441  fMCTIndexMap[trackID] = primarymctIndex;
442 
443  fMCTPrimProcessKeepMap[trackID] = isFromMCTProcessPrimary;
444 
445  // -- determine whether full set of trajectorie points should be stored or only the start and end points
447  (!fstoreTrajectories ||
448  (notstore && fStoreDroppedMCParticles)) ? //also don't store if dropped particle
449  false : /*don't want trajectory points at all, bail*/
450  (!(fMCTIndexToGeneratorMap[primarymctIndex].second)) ?
451  false : /*particle is not from a storable generator*/
453  true : /*want all primaries tracked for a storable generator*/
454  (isFromMCTProcessPrimary) ?
455  true : /*only descendants from primaries with MCTruth process == "primary"*/
456  false; /*not from MCTruth process "primary"*/
457  // Polarization.
458  const G4ThreeVector& polarization = track->GetPolarization();
460  TVector3{polarization.x(), polarization.y(), polarization.z()});
461 
462  if (track->GetProperTime() != 0) { return; }
463 
464  // if KeepEMShowerDaughters = False and we decided to drop this particle,
465  if (notstore) { // this bool checks if particle is eliminated by NotStoredPhysics
467  return;
468  }
469  // Save the particle in the ParticleList.
470 
471  // if we are not filtering, we have a decision already. The extra check is to see if we are dropping
472  // particle from a process that is not stored. We don't do the same for dropped particles
473  // since if it doesn't have a filter, we don't produce a separate list anyways
476  }
477 
478  //----------------------------------------------------------------------------
480  {
481  if (!fCurrentParticle.hasParticle()) { return; }
482 
483  if (aTrack) {
484  fCurrentParticle.particle->SetWeight(aTrack->GetWeight());
485 
486  // Get the post-step information from the G4Step.
487  const G4StepPoint* postStepPoint = aTrack->GetStep()->GetPostStepPoint();
488  if (!postStepPoint->GetProcessDefinedStep()) {
489  // Now we get to do some awkward cleanup because the
490  // fParticleList was augmented during the
491  // preUserTrackingAction. We cannot call 'Archive' because
492  // that only sets the mapped type of the entry to
493  // nullptr...which is really bad whenever we iterate through
494  // entries in the endOfEventAction and dereference the mapped
495  // type. We have to entirely erase the entry.
496  auto key_to_erase = fParticleList.key(fCurrentParticle.particle);
497  fParticleList.erase(key_to_erase);
499  //Check if particle is in dropped list - edge case
500  if (fdroppedParticleList->KnownParticle(fCurrentParticle.particle->TrackId())) {
501  //If it is, archive it
503  }
504  }
505  // after the particle is archived, it is deleted
507  return;
508  }
509  G4String process = postStepPoint->GetProcessDefinedStep()->GetProcessName();
511 
512  // -- D.R. Store the final point only for particles that have not had intermediate trajectory
513  // points saved. This avoids double counting the final trajectory point for particles from
514  // generators with storable trajectory points.
515 
517  const G4ThreeVector position = postStepPoint->GetPosition();
518  G4double time = postStepPoint->GetGlobalTime();
519 
520  // Remember that LArSoft uses cm, ns, GeV.
521  TLorentzVector fourPos(position.x() / CLHEP::cm,
522  position.y() / CLHEP::cm,
523  position.z() / CLHEP::cm,
524  time / CLHEP::ns);
525 
526  const G4ThreeVector momentum = postStepPoint->GetMomentum();
527  const G4double energy = postStepPoint->GetTotalEnergy();
528  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
529  momentum.y() / CLHEP::GeV,
530  momentum.z() / CLHEP::GeV,
531  energy / CLHEP::GeV);
532 
533  // Add another point in the trajectory.
534  AddPointToCurrentParticle(fourPos, fourMom, std::string(process));
535  }
536  // -- particle has a full trajectory, apply SparsifyTrajectory method if enabled
537  else if (fSparsifyTrajectories) {
539  }
540  }
541 
543  auto key_to_erase = fParticleList.key(fCurrentParticle.particle);
544  //Erase primaries
545  if (!fCurrentParticle.isDropped) fParticleList.erase(key_to_erase);
546  //Erase dropped particles
549  }
550  //
551  int const trackID = aTrack->GetTrackID() + fTrackIDOffset;
552  int parentID = aTrack->GetParentID() + fTrackIDOffset;
553  fdroppedTracksMap[parentID].insert(trackID);
555  // do add the particle to the parent id map though
556  // and set the current track id to be it's ultimate parent
557  fParentIDMap[trackID] = parentID;
558  fCurrentTrackID = -1 * this->GetParentage(trackID);
559  fTargetIDMap[trackID] = fCurrentTrackID;
560  }
561 
562  // store truth record pointer, only if it is available
563  if (fCurrentParticle.isPrimary()) {
565  }
566  return;
567  }
568 
569  //----------------------------------------------------------------------------
570  // With every step, add to the particle's trajectory.
572  {
573  // N.B. G4 guarantees that following are non-null:
574  // - step
575  // - step->GetPostStepPoint()
576  if (!fCurrentParticle.hasParticle() || !step->GetPostStepPoint()->GetProcessDefinedStep()) {
577  return;
578  }
579  // Temporary fix for problem where DeltaTime on the first step
580  // of optical photon propagation is calculated incorrectly. -wforeman
581  double const globalTime = step->GetTrack()->GetGlobalTime();
582  double const velocity_G4 = step->GetTrack()->GetVelocity();
583  double const velocity_step = step->GetStepLength() / step->GetDeltaTime();
584  if ((step->GetTrack()->GetDefinition()->GetPDGEncoding() == 0) &&
585  fabs(velocity_G4 - velocity_step) > 0.0001) {
586  // Subtract the faulty step time from the global time,
587  // and add the correct step time based on G4 velocity.
588  step->GetPostStepPoint()->SetGlobalTime(globalTime - step->GetDeltaTime() +
589  step->GetStepLength() / velocity_G4);
590  }
591 
592  // For the most part, we just want to add the post-step
593  // information to the particle's trajectory. There's one
594  // exception: In PreTrackingAction, the correct time information
595  // is not available. So add the correct vertex information here.
596 
598 
599  // Get the pre/along-step information from the G4Step.
600  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
601 
602  const G4ThreeVector position = preStepPoint->GetPosition();
603  G4double time = preStepPoint->GetGlobalTime();
604 
605  // Remember that LArSoft uses cm, ns, GeV.
606  TLorentzVector fourPos(position.x() / CLHEP::cm,
607  position.y() / CLHEP::cm,
608  position.z() / CLHEP::cm,
609  time / CLHEP::ns);
610 
611  const G4ThreeVector momentum = preStepPoint->GetMomentum();
612  const G4double energy = preStepPoint->GetTotalEnergy();
613  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
614  momentum.y() / CLHEP::GeV,
615  momentum.z() / CLHEP::GeV,
616  energy / CLHEP::GeV);
617 
618  // Add the first point in the trajectory.
619  AddPointToCurrentParticle(fourPos, fourMom, "Start");
620  } // end if this is the first step
621 
622  // At this point, the particle is being transported through the
623  // simulation.
624  // change below:
625  // This method is being called for every step that
626  // the track passes through, but we don't want to update the
627  // trajectory information if the step was defined by the StepLimiter.
628  G4String process = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
629  G4bool ignoreProcess = process.contains("StepLimiter");
630 
631  // We store the initial creation point of the particle
632  // and its final position (ie where it has no more energy, or at least < 1 eV) no matter
633  // what, but whether we store the rest of the trajectory depends
634  // on the process, and on a user switch.
635  // -- D.R. Store additional trajectory points only for desired generators and processes
636  if (!ignoreProcess && fCurrentParticle.keepFullTrajectory) {
637 
638  // Get the post-step information from the G4Step.
639  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
640 
641  const G4ThreeVector position = postStepPoint->GetPosition();
642  G4double time = postStepPoint->GetGlobalTime();
643 
644  // Remember that LArSoft uses cm, ns, GeV.
645  TLorentzVector fourPos(position.x() / CLHEP::cm,
646  position.y() / CLHEP::cm,
647  position.z() / CLHEP::cm,
648  time / CLHEP::ns);
649 
650  const G4ThreeVector momentum = postStepPoint->GetMomentum();
651  const G4double energy = postStepPoint->GetTotalEnergy();
652  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
653  momentum.y() / CLHEP::GeV,
654  momentum.z() / CLHEP::GeV,
655  energy / CLHEP::GeV);
656 
657  // Add another point in the trajectory.
658  AddPointToCurrentParticle(fourPos, fourMom, std::string(process));
659  }
660  }
661 
662  //----------------------------------------------------------------------------
666  public:
667  explicit UpdateDaughterInformation(sim::ParticleList& p) : particleList{&p} {}
668  void operator()(sim::ParticleList::value_type& particleListEntry) const
669  {
670  // We're looking at this Particle in the list.
671  int particleID = particleListEntry.first;
672 
673  // The parent ID of this particle;
674  // we ask the particle list since the particle itself might have been lost
675  // ("archived"), but the particle list still holds the information we need
676  int parentID = particleList->GetMotherOf(particleID);
677 
678  // If the parentID <= 0, this is a primary particle.
679  if (parentID <= 0) return;
680 
681  // If we get here, this particle is somebody's daughter. Add
682  // it to the list of daughter particles for that parent.
683 
684  // Get the parent particle from the list.
685  sim::ParticleList::iterator parentEntry = particleList->find(parentID);
686 
687  if (parentEntry == particleList->end()) {
688  // We have an "orphan": a particle whose parent isn't
689  // recorded in the particle list. This is not signficant;
690  // it's possible for a particle not to be saved in the list
691  // because it failed an energy cut, but for it to have a
692  // daughter that passed the cut (e.g., a nuclear decay).
693  return;
694  }
695  if (!parentEntry->second) return; // particle archived, nothing to update
696 
697  // Add the current particle to the daughter list of the parent.
698  simb::MCParticle* parent = parentEntry->second;
699  parent->AddDaughter(particleID);
700  }
701 
702  private:
704  };
705 
706  //----------------------------------------------------------------------------
708  {
709  auto const it = fPrimaryTruthMap.find(trackId);
710  return it == fPrimaryTruthMap.end() ? simb::NoGeneratedParticleIndex : it->second;
711  }
712 
713  //----------------------------------------------------------------------------
714  // Yields the ParticleList accumulated during the current event.
716  {
717  // check if the ParticleNavigator has entries, and if
718  // so grab the highest track id value from it to
719  // add to the fTrackIDOffset
720  int highestID = 0;
721  for (auto pn = fParticleList.begin(); pn != fParticleList.end(); pn++)
722  if ((*pn).first > highestID) highestID = (*pn).first;
723 
724  // If we have stored dropped particles,
725  // include them in the offset.
726  if (fdroppedParticleList) {
727  for (auto pn = fdroppedParticleList->begin(); pn != fdroppedParticleList->end(); pn++)
728  if ((*pn).first > highestID) highestID = (*pn).first;
729  }
730 
731  //Only change the fTrackIDOffset if there is in fact a particle to add to the event
732  if ((fParticleList.size()) != 0) {
733  fTrackIDOffset = highestID + 1;
734  mf::LogDebug("YieldList:fTrackIDOffset")
735  << "highestID = " << highestID << "\nfTrackIDOffset= " << fTrackIDOffset;
736  }
737  return std::move(fParticleList);
738  } // ParticleList&& ParticleListActionService::YieldList()
739  //----------------------------------------------------------------------------
740  // Yields the (dropped) ParticleList accumulated during the current event.
742  {
743  if (!fdroppedParticleList) {
744  throw cet::exception("ParticleListAction")
745  << "ParticleListAction::YieldDroppedList(): dropped particle list not build by user "
746  "request.\n";
747  }
748  // check if the ParticleNavigator has entries, and if
749  // so grab the highest track id value from it to
750  // add to the fTrackIDOffset
751  int highestID = 0;
752  for (auto pn = fParticleList.begin(); pn != fParticleList.end(); pn++)
753  if ((*pn).first > highestID) highestID = (*pn).first;
754 
755  // If we have stored dropped particles,
756  // include them in the offset.
757  for (auto pn = fdroppedParticleList->begin(); pn != fdroppedParticleList->end(); pn++)
758  if ((*pn).first > highestID) highestID = (*pn).first;
759 
760  //Only change the fTrackIDOffset if there is in fact a particle to add to the event
761  if ((fParticleList.size()) != 0) {
762  fTrackIDOffset = highestID + 1;
763  mf::LogDebug("YieldList:fTrackIDOffset")
764  << "highestID = " << highestID << "\nfTrackIDOffset= " << fTrackIDOffset;
765  }
766  return std::move(*fdroppedParticleList);
767  } // ParticleList&& ParticleListAction::YieldDroppedList()
768  //----------------------------------------------------------------------------
769  // Dropped particle test
770 
772  {
773  return !p || p->Trajectory().empty();
774  } // ParticleListAction::isDropped()
775  //----------------------------------------------------------------------------
777  TLorentzVector const& mom,
778  std::string const& process)
779  {
780  // Add the first point in the trajectory.
782 
783  // also see if we can decide to keep the particle
785  fCurrentParticle.isInVolume = fFilter->mustKeep(pos);
787  fCurrentParticle.isInVolume = fDroppedFilter->mustKeep(pos);
788  }
789 
790  // Called at the end of each event. Call detectors to convert hits for the
791  // event and pass the call on to the action objects.
793  {
794  // -- End of Run Report
795  if (!fNotStoredCounterUMap.empty()) { // -- Only if there is something to report
796  std::stringstream sscounter;
797  sscounter << "Not Stored Process summary:";
798  for (auto const& [process, count] : fNotStoredCounterUMap) {
799  sscounter << "\n\t" << process << " : " << count;
800  }
801  mf::LogInfo("ParticleListActionService") << sscounter.str();
802  }
803 
804  partCol_ = std::make_unique<std::vector<simb::MCParticle>>();
805  droppedCol_ = std::make_unique<sim::ParticleAncestryMap>();
806  droppedPartCol_ = std::make_unique<std::vector<simb::MCParticle>>();
807  tpassn_ =
808  std::make_unique<art::Assns<simb::MCTruth, simb::MCParticle, sim::GeneratedParticleInfo>>();
809  // Set up the utility class for the "for_each" algorithm. (We only
810  // need a separate set-up for the utility class because we need to
811  // give it the pointer to the particle list. We're using the STL
812  // "for_each" instead of the C++ "for loop" because it's supposed
813  // to be faster.
814  std::for_each(
816 
817  MF_LOG_INFO("endOfEventAction") << "MCTruth Handles Size: " << fMCLists->size();
818 
819  unsigned int nGeneratedParticles = 0;
820  unsigned int nMCTruths = 0;
821  sim::ParticleList particleList = YieldList();
822  // Request a list of dropped particles
823  sim::ParticleList droppedParticleList;
824  if (fdroppedParticleList) { droppedParticleList = YieldDroppedList(); }
825  for (size_t mcl = 0; mcl < fMCLists->size(); ++mcl) {
826  auto const& mclistHandle = (*fMCLists)[mcl];
827  MF_LOG_INFO("endOfEventAction") << "mclistHandle Size: " << mclistHandle->size();
828  for (size_t m = 0; m < mclistHandle->size(); ++m) {
829  art::Ptr<simb::MCTruth> mct(mclistHandle, m);
830  MF_LOG_INFO("endOfEventAction") << "Found " << mct->NParticles() << " particles";
831  for (simb::MCParticle* p : particleList | ranges::views::values) {
832  auto gen_index = fMCTIndexMap[p->TrackId()];
833  if (gen_index != nMCTruths) continue;
834  // if the particle has been marked as dropped, we don't save it
835  // (as of LArSoft ~v5.6 this does not ever happen because
836  // ParticleListAction has already taken care of deleting them)
837  // if (isDropped(p)) continue;
838  assert(p->NumberTrajectoryPoints() != 0ull);
839  ++nGeneratedParticles;
840  sim::GeneratedParticleInfo const truthInfo{GetPrimaryTruthIndex(p->TrackId())};
841  if (!truthInfo.hasGeneratedParticleIndex() && p->Mother() == 0) {
842  MF_LOG_WARNING("endOfEventAction") << "No GeneratedParticleIndex()!";
843  // this means it's primary but with no information; logic error!!
845  << "Failed to match primary particle:\n"
846  << "\nwith particles from the truth record '" << mclistHandle.provenance()->inputTag()
847  << "':\n\n";
848  }
849  partCol_->push_back(std::move(*p));
851  tpassn_->addSingle(mct, mcp_ptr, truthInfo);
852 
853  } // endfor p in particleList
854  if (fStoreDroppedMCParticles && droppedPartCol_) {
855  for (simb::MCParticle* p : droppedParticleList | ranges::views::values) {
856  if (isDropped(p)) continue; //Is it dropped??
857  if (p->StatusCode() != 1) continue; //Is it a primary particle??
858 
859  droppedPartCol_->push_back(std::move(*p));
860  } // for(droppedParticleList)
861  } // if (fStoreDroppedMCParticles && droppedPartCol_)
862  mf::LogDebug("Offset") << "nGeneratedParticles = " << nGeneratedParticles;
863  droppedCol_->SetMap(fdroppedTracksMap);
864  ++nMCTruths;
865  }
866  }
867  fTrackIDOffset = 0;
868  }
869 } // namespace LArG4
unsigned int NumberTrajectoryPoints() const
Definition: MCParticle.h:219
simb::MCParticle * particle
simple structure representing particle
bool isDropped(simb::MCParticle const *p)
void AddDaughter(const int trackID)
Definition: MCParticle.h:269
void userSteppingAction(const G4Step *) override
ParticleListActionService(fhicl::ParameterSet const &)
void beginOfEventAction(const G4Event *) override
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:258
decltype(auto) constexpr cend(T &&obj)
ADL-aware version of std::cend.
Definition: StdUtils.h:93
std::map< int, std::set< int > > fdroppedTracksMap
map <ParentID, set: list of track ids for which no MCParticle was created>
void preUserTrackingAction(const G4Track *) override
std::unique_ptr< sim::ParticleList > fdroppedParticleList
bool fKeepEMShowerDaughters
whether to keep EM shower secondaries, tertiaries, etc
key_type key(mapped_type const &part) const
Extracts the key from the specified value.
Definition: ParticleList.h:338
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
const simb::MCTrajectory & Trajectory() const
Definition: MCParticle.h:254
bool fKeepSecondToLast
tell whether or not to force keeping the second to last point
list_type::value_type value_type
Definition: ParticleList.h:130
std::unique_ptr< art::Assns< simb::MCTruth, simb::MCParticle, sim::GeneratedParticleInfo > > tpassn_
list_type::iterator iterator
Definition: ParticleList.h:131
std::unique_ptr< util::PositionInVolumeFilter > fFilter
filter for particles to be kept
void Add(simb::MCParticle *value)
Definition: ParticleList.h:315
Geant4 interface.
int GetParentage(int trackid) const
STL namespace.
double fSparsifyMargin
set the sparsification margin
void AddPointToCurrentParticle(TLorentzVector const &pos, TLorentzVector const &mom, std::string const &process)
Adds a trajectory point to the current particle, and runs the filter.
int NParticles() const
Definition: MCTruth.h:75
void postUserTrackingAction(const G4Track *) override
double velocity_step
std::string Process() const
Definition: MCParticle.h:216
simb::GeneratedParticleIndex_t GetPrimaryTruthIndex(int trackId) const
Returns the index of primary truth (sim::NoGeneratorIndex if none).
bool empty() const
Definition: MCTrajectory.h:167
int TrackId() const
Definition: MCParticle.h:211
bool fKeepTransportation
tell whether or not to keep the transportation process
std::unique_ptr< std::vector< simb::MCParticle > > partCol_
std::map< int, size_t > fMCTIndexMap
Map: particle track ID -> index of primary parent in std::vector<simb::MCTruth> object.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
std::vector< std::string > fKeepParticlesInVolumes
Only write particles that have trajectories through these volumes.
constexpr GeneratedParticleIndex_t NoGeneratedParticleIndex
Constant representing the absence of generator truth information.
Definition: simb.h:34
bool fSparsifyTrajectories
help reduce the number of trajectory points.
void SetPolarization(const TVector3 &p)
Definition: MCParticle.h:270
sim::ParticleList && YieldDroppedList()
static const int NoParticleId
Definition: sim.h:21
decltype(auto) values(Coll &&coll)
Range-for loop helper iterating across the values of the specified collection.
Use Geant4&#39;s user "hooks" to maintain a list of particles generated by Geant4.
std::map< size_t, std::pair< std::string, G4bool > > fMCTIndexToGeneratorMap
Map: MCTruthIndex -> generator, input label of generator and keepGenerator decision.
double energy
Definition: plottest35.C:25
simb::MCParticle const * GetMCParticle() const
std::unordered_map< std::string, int > fNotStoredCounterUMap
Map: not stored process and counter.
iterator begin()
Definition: ParticleList.h:305
bool fStoreDroppedMCParticles
Whether to keep a sim::MCParticle list of dropped particles.
std::unique_ptr< sim::ParticleAncestryMap > droppedCol_
std::map< int, int > fTargetIDMap
key is original track ID, value is ID to assign for downstream objs (e.g. SimEdeps) ...
std::unique_ptr< std::vector< simb::MCParticle > > droppedPartCol_
This collection will hold the MCParticleLite objects created from dropped particles.
void endOfEventAction(const G4Event *) override
double velocity_G4
#define MF_LOG_INFO(category)
std::vector< std::string > fKeepDroppedParticlesInVolumes
std::unique_ptr< util::PositionInVolumeFilter > fDroppedFilter
bool isPrimary() const
Returns whether there is a particle.
void SetWeight(double wt)
Definition: MCParticle.h:272
void SetEndProcess(std::string s)
Definition: MCParticle.cxx:111
std::vector< std::string > fNotStoredPhysics
Physics processes that will not be stored.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::map< int, int > fParentIDMap
key is current track ID, value is parent ID
simb::GeneratedParticleIndex_t truthInfoIndex() const
Returns the index of the particle in the generator truth record.
double globalTime
art::EDProductGetter const * productGetter_
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
void clear()
Resets the information (does not release memory it does not own)
std::vector< std::string > fkeepGenTrajectories
std::map< int, bool > fMCTPrimProcessKeepMap
Map: particle trakc ID -> boolean decision to keep or not full trajectory points. ...
simb::GeneratedParticleIndex_t truthIndex
Index of the particle in the original generator truth record.
bool hasParticle() const
Returns whether there is a particle.
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
Contains information about a generated particle.
G4bool fstoreTrajectories
Whether to store particle trajectories with each particle.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::map< int, simb::GeneratedParticleIndex_t > fPrimaryTruthMap
Map: particle track ID -> index of primary information in MC truth.
UpdateDaughterInformation(sim::ParticleList &p)
size_type erase(const key_type &key)
bool KnownParticle(int trackID) const
Returns whether we have had this particle, archived or live.
Definition: ParticleList.h:215
#define MF_LOG_WARNING(category)
void operator()(sim::ParticleList::value_type &particleListEntry) const
void SparsifyTrajectory(double margin=0.1, bool keep_second_to_last=false)
Definition: MCParticle.h:266
Float_t track
Definition: plot.C:35
std::vector< art::Handle< std::vector< simb::MCTruth > > > const * fMCLists
MCTruthCollection input lists.
size_type size() const
Definition: ParticleList.h:313
Tools and modules for checking out the basics of the Monte Carlo.
GeneratedParticleIndex_t MCParticleIndex() const
Returns the index of the corresponding particle in truth record.
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