All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ParticleListAction.cc
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 /// Design considerations
8 /// ---------------------
9 ///
10 /// This class relies on the MCTruth index from
11 /// g4b::PrimaryParticleInformation to operate correctly. This index
12 /// is an integer value that corresponds to an MCTruth object, as
13 /// accessed through art::Handle<std::vector<simb::MCTruth>> objects.
14 /// However, the order in which MCTruth objects are processed must be
15 /// consistent between this service and the MCTruthEventAction
16 /// service, which creates the PrimaryParticleInformation object,
17 /// otherwise the Assns objects created here will be incorrect.
18 ///
19 /// Through art 3.09, one can rely on the order returned by a given
20 /// Event::getMany call to be predictable and consistent within the
21 /// same program. However, this behavior should not necessarily be
22 /// relied upon, and a different implementation of this class would
23 /// insulate users from such details, making the implementation
24 /// simpler. One should determine whether storing an art::ProductID
25 /// object along with an MCTruthIndex might be more helpful.
26 ///
27 ////////////////////////////////////////////////////////////////////////
28 
30 
32 
33 #include "nug4/G4Base/PrimaryParticleInformation.h"
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/G4ParticleDefinition.hh"
45 #include "Geant4/G4DynamicParticle.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  {
81  // -- D.R. If a custom list of not storable physics is provided, use it, otherwise
82  // use the default list. This preserves the behavior of the keepEmShowerDaughters
83  // parameter
84  bool customNotStored = not fNotStoredPhysics.empty();
86  // -- Don't keep all processes
87  if (!customNotStored) // -- Don't keep but haven't provided a list
88  { // -- default list of not stored physics
89  fNotStoredPhysics = {"conv",
90  "LowEnConversion",
91  "Pair",
92  "compt",
93  "Compt",
94  "Brem",
95  "phot",
96  "Photo",
97  "Ion",
98  "annihil"};
99  }
100 
101  std::stringstream sstored;
102  sstored << "The full tracking information will not be stored for particles"
103  << " resulting from the following processes: \n{ ";
104  for (auto const& i : fNotStoredPhysics) {
105  sstored << "\"" << i << "\" ";
106  fNotStoredCounterUMap.emplace(i, 0); // -- initialize counter
107  }
108  mf::LogInfo("ParticleListActionService") << sstored.str() << "}\n";
109  }
110  else { // -- Keep all processes
111  mf::LogInfo("ParticleListActionService")
112  << "Storing full tracking information for all processes. \n";
113  if (customNotStored) // -- custom list will be ignored
114  {
115  mf::LogWarning("StoredPhysics")
116  << "NotStoredPhysics provided, but will be ignored."
117  << " To use NotStoredPhysics, set keepEMShowerDaughters to false";
118  }
119  }
120 
121  // -- sparsify info
123  mf::LogInfo("ParticleListActionService")
124  << "Trajectory sparsification enabled with SparsifyMargin : " << fSparsifyMargin << "\n";
125  }
126 
127  //----------------------------------------------------------------------------
128  // Begin the event
129  void
130  ParticleListActionService::beginOfEventAction(const G4Event*)
131  {
132  // Clear any previous particle information.
133  fCurrentParticle.clear();
134  fParticleList.clear();
135  fParentIDMap.clear();
136  fTargetIDMap.clear();
137  fMCTIndexMap.clear();
138  fMCTPrimProcessKeepMap.clear();
139  fCurrentTrackID = sim::NoParticleId;
140  fTrackIDOffset = 0;
141  fPrimaryTruthMap.clear();
142  fMCTIndexToGeneratorMap.clear();
143  fNotStoredCounterUMap.clear();
144  fdroppedTracksMap.clear();
145  // -- D.R. If a custom list of keepGenTrajectories is provided, use it, otherwise
146  // keep or drop decision made based storeTrajectories parameter. This preserves
147  // the behavior of the storeTrajectories fhicl param
148  bool customKeepTraj = not fkeepGenTrajectories.empty();
149  if (!fstoreTrajectories) { // -- fstoreTrajectories : false
150  mf::LogDebug("beginOfEventAction::Generator") << "Trajectory points will not be stored.";
151  }
152  else if (!customKeepTraj) { // -- fstoretrajectories : true and empty keepGenTrajectories list
153  mf::LogDebug("beginOfEventAction::Generator")
154  << "keepGenTrajectories list is empty. Will"
155  << " store trajectory points for all generators";
156  }
157 
158  // -- D.R. determine mapping between MCTruthIndex(s) and generator(s) for later reference
159  size_t nKeep = 0;
160  std::string generator_name = "unknown";
161  for (size_t mcti = 0; mcti < fMCLists->size(); ++mcti) {
162 
163  std::stringstream sskeepgen;
164  sskeepgen << "MCTruth object summary :";
165  sskeepgen << "\n\tPrimary MCTIndex : " << mcti;
166 
167  // -- Obtain the generator (provenance) corresponding to the mctruth index:
168  auto const& mclistHandle = (*fMCLists)[mcti];
169  generator_name = mclistHandle.provenance()->inputTag().label();
170  sskeepgen << "\n\tProvenance/Generator : " << generator_name;
171 
172  G4bool keepGen = false;
173  if (fstoreTrajectories) // -- storeTrajectories set to true; check which
174  {
175  if (!customKeepTraj) { // -- no custom list, so keep all
176  keepGen = true;
177  nKeep++;
178  }
179  else { // -- custom list, so check the ones in the event against provided keep list
180  for (auto keepableGen : fkeepGenTrajectories) {
181  if (generator_name == keepableGen) { // -- exit upon finding match; false by default
182  keepGen = true;
183  nKeep++;
184  break;
185  }
186  }
187  }
188  }
189  fMCTIndexToGeneratorMap.emplace(mcti, std::make_pair(generator_name, keepGen));
190  sskeepgen << "\n\tTrajectory points storable : " << (keepGen ? "true" : "false") << "\n";
191  mf::LogDebug("beginOfEventAction::Generator") << sskeepgen.str();
192  }
193 
194  if (nKeep == 0 && customKeepTraj && fstoreTrajectories) {
195  mf::LogWarning("beginOfEventAction::keepableGenerators")
196  << "storeTrajectories"
197  << " set to true and a non-empty keepGenTrajectories list provided in configuration file, "
198  "but"
199  << " none of the generators in this list are present in the event! Double check list or "
200  "don't"
201  << " provide keepGenTrajectories in the configuration to keep all trajectories from all"
202  << " generator labels. This may be expected for generators that have a nonzero probability "
203  "of"
204  << " producing no particles (e.g. some radiologicals)";
205  }
206  }
207 
208 
209  //-------------------------------------------------------------
210  // figure out the ultimate parentage of the particle with track ID
211  // trackid
212  // assume that the current track id has already been added to
213  // the fParentIDMap
214  int
215  ParticleListActionService::GetParentage(int trackid) const
216  {
217  int parentid = sim::NoParticleId;
218 
219  // search the fParentIDMap recursively until we have the parent id
220  // of the first EM particle that led to this one
221  auto itr = fParentIDMap.find(trackid);
222  while (itr != fParentIDMap.end()) {
223  // set the parentid to the current parent ID, when the loop ends
224  // this id will be the first EM particle
225  parentid = (*itr).second;
226  itr = fParentIDMap.find(parentid);
227  }
228 
229  return parentid;
230  }
231 
232  //----------------------------------------------------------------------------
233  // Create our initial simb::MCParticle object and add it to the sim::ParticleList.
234  void
235  ParticleListActionService::preUserTrackingAction(const G4Track* track)
236  {
237  // Particle type.
238  G4ParticleDefinition* particleDefinition = track->GetDefinition();
239  G4int pdgCode = particleDefinition->GetPDGEncoding();
240 
241  // Get Geant4's ID number for this track. This will be the same
242  // ID number that we'll use in the ParticleList.
243  // It is offset by the number of tracks accumulated from the previous Geant4
244  // runs (if any)
245  int const trackID = track->GetTrackID() + fTrackIDOffset;
246  fCurrentTrackID = trackID;
247  fTargetIDMap[trackID] = fCurrentTrackID;
248  // And the particle's parent (same offset as above):
249  int parentID = track->GetParentID() + fTrackIDOffset;
250 
251  std::string process_name = "unknown";
252  std::string mct_primary_process = "unknown";
253  bool isFromMCTProcessPrimary = false;
254 
255  // Is there an MCTruth object associated with this G4Track? We
256  // have to go up a "chain" of information to find out:
257  const G4DynamicParticle* dynamicParticle = track->GetDynamicParticle();
258  const G4PrimaryParticle* primaryParticle = dynamicParticle->GetPrimaryParticle();
259  simb::GeneratedParticleIndex_t primaryIndex = simb::NoGeneratedParticleIndex;
260  size_t primarymctIndex = 0;
261  if (primaryParticle != nullptr) {
262  const G4VUserPrimaryParticleInformation* gppi = primaryParticle->GetUserInformation();
263  const g4b::PrimaryParticleInformation* ppi =
264  dynamic_cast<const g4b::PrimaryParticleInformation*>(gppi);
265  if (ppi != nullptr) {
266  primaryIndex = ppi->MCParticleIndex();
267  primarymctIndex = ppi->MCTruthIndex();
268  mct_primary_process = ppi->GetMCParticle()->Process();
269 
270  // If we've made it this far, a PrimaryParticleInformation
271  // object exists and we are using a primary particle, set the
272  // process name accordingly
273 
274  // -- D.R. : if process == "primary" exactly (this will most likely be the case), mark it as
275  // a primary with keepable trajectory for itself and its descendants.
276  // -- elsif it simply starts with "primary", accept it but mark it as non-keep
277  // -- else force it to be primary because we have determined that it is a primary particle
278  // This was the original behavior of the code i.e. process was _set_ to "primary"
279  // regardless of what the process name was in the generator MCTruth object.
280  //
281  // -- NOTE: This enforces a convention for process names assigned in the gen stage.
282  if (mct_primary_process.compare("primary") == 0) {
283  process_name = "primary";
284  isFromMCTProcessPrimary = true;
285  }
286  else if (mct_primary_process.find("primary") == 0) {
287  process_name = mct_primary_process;
288  isFromMCTProcessPrimary = false;
289  mf::LogDebug("PrimaryParticle") << "MCTruth primary process name contains \"primary\" "
290  << " but is not solely \"primary\" : " << process_name
291  << ".\nWill not store full set of trajectory points.";
292  }
293  else { // -- override it
294  process_name = "primary";
295  isFromMCTProcessPrimary = true;
296  mf::LogWarning("PrimaryParticle")
297  << "MCTruth primary process does not beging with string"
298  << " literal \"primary\" : " << process_name << "\nOVERRIDING it to \"primary\"";
299  }
300  // -- The process_name check is simply to allow an additional way to reduce memory usage,
301  // namely, by creating MCTruth input particles with multiple process labels (e.g. "primary"
302  // and "primaryBackground") that can be used to veto storage of trajectory points.
303 
304  // primary particles should have parentID = 0, even if there
305  // are multiple MCTruths for this event
306  parentID = 0;
307  } // end else no primary particle information
308  } // Is there a G4PrimaryParticle?
309  // If this is not a primary particle...
310  else {
311  // check if this particle was made in an undesirable process. For example:
312  // if one is not interested in EM shower particles, don't put it in the particle
313  // list as one wouldn't care about secondaries, tertiaries, etc. For these showers
314  // figure out what process is making this track - skip it if it is
315  // one of pair production, compton scattering, photoelectric effect
316  // bremstrahlung, annihilation, or ionization
317  process_name = track->GetCreatorProcess()->GetProcessName();
318  if (!fKeepEMShowerDaughters) {
319  bool notstore = false;
320  for (auto const& p : fNotStoredPhysics) {
321  if (process_name.find(p) != std::string::npos) {
322  notstore = true;
323  mf::LogDebug("NotStoredPhysics") << "Found process : " << process_name;
324 
325  int old = 0;
326  auto search = fNotStoredCounterUMap.find(p);
327  if (search != fNotStoredCounterUMap.end()) { old = search->second; }
328  fNotStoredCounterUMap.insert_or_assign(p, (old + 1));
329  break;
330  }
331  }
332  if (notstore) {
333  // figure out the ultimate parentage of this particle
334  // first add this track id and its parent to the fParentIDMap
335  fParentIDMap[trackID] = parentID;
336  fCurrentTrackID =-1 * this->GetParentage(trackID);
337  // check that fCurrentTrackID is in the particle list - it is possible
338  // that this particle's parent is a particle that did not get tracked.
339  // An example is a parent that was made due to muMinusCaptureAtRest
340  // and the daughter was made by the phot process. The parent likely
341  // isn't saved in the particle list because it is below the energy cut
342  // which will put a bogus track id value into the sim::IDE object for
343  // the sim::SimChannel if we don't check it.
344  if (!fParticleList.KnownParticle(fCurrentTrackID)) fCurrentTrackID = sim::NoParticleId;
345  fTargetIDMap[trackID] = fCurrentTrackID;
346  // clear current particle as we are not stepping this particle and
347  // adding trajectory points to it
348  fdroppedTracksMap[this->GetParentage(trackID)].insert(trackID);
349  fCurrentParticle.clear();
350  return;
351  } // end if process matches an undesired process
352  } // end if keeping EM shower daughters
353 
354  // Check the energy of the particle. If it falls below the energy
355  // cut, don't add it to our list.
356  G4double energy = track->GetKineticEnergy();
357  if (energy < fenergyCut) {
358  fdroppedTracksMap[this->GetParentage(trackID)].insert(trackID);
359  fCurrentParticle.clear();
360  // do add the particle to the parent id map though
361  // and set the current track id to be it's ultimate parent
362  fParentIDMap[trackID] = parentID;
363  fCurrentTrackID = -1 * this->GetParentage(trackID);
364  fTargetIDMap[trackID] = fCurrentTrackID;
365  return;
366  }
367 
368  // check to see if the parent particle has been stored in the particle navigator
369  // if not, then see if it is possible to walk up the fParentIDMap to find the
370  // ultimate parent of this particle. Use that ID as the parent ID for this
371  // particle
372  if (!fParticleList.KnownParticle(parentID)) {
373  // do add the particle to the parent id map
374  // just in case it makes a daughter that we have to track as well
375  fParentIDMap[trackID] = parentID;
376  int pid = this->GetParentage(parentID);
377 
378  // if we still can't find the parent in the particle navigator,
379  // we have to give up
380  if (!fParticleList.KnownParticle(pid)) {
381  MF_LOG_WARNING("ParticleListActionService")
382  << "can't find parent id: " << parentID << " in the particle list, or fParentIDMap."
383  << " Make " << parentID << " the mother ID for"
384  << " track ID " << fCurrentTrackID << " in the hope that it will aid debugging.";
385  }
386  else
387  parentID = pid;
388  }
389 
390  // Once the parentID is secured, inherit the MCTruth Index
391  // which should have been set already
392  if (auto it = fMCTIndexMap.find(parentID); it != cend(fMCTIndexMap)) {
393  primarymctIndex = it->second;
394  }
395  else {
396  throw art::Exception(art::errors::LogicError)
397  << "Could not locate MCT Index for parent trackID of " << parentID;
398  }
399 
400  // Inherit whether the parent is from a primary with MCTruth process_name == "primary"
401  isFromMCTProcessPrimary = fMCTPrimProcessKeepMap[parentID];
402  } // end if not a primary particle
403 
404  // This is probably the PDG mass, but just in case:
405  double mass = dynamicParticle->GetMass() / CLHEP::GeV;
406 
407  // Create the sim::Particle object.
408  fCurrentParticle.clear();
409  fCurrentParticle.particle =
410  new simb::MCParticle{trackID, pdgCode, process_name, parentID, mass};
411  fCurrentParticle.truthIndex = primaryIndex;
412 
413  fMCTIndexMap[trackID] = primarymctIndex;
414 
415  fMCTPrimProcessKeepMap[trackID] = isFromMCTProcessPrimary;
416 
417  // -- determine whether full set of trajectorie points should be stored or only the start and end points
418  fCurrentParticle.keepFullTrajectory =
419  (!fstoreTrajectories) ?
420  false : /*don't want trajectory points at all, bail*/
421  (!(fMCTIndexToGeneratorMap[primarymctIndex].second)) ?
422  false : /*particle is not from a storable generator*/
424  true : /*want all primaries tracked for a storable generator*/
425  (isFromMCTProcessPrimary) ?
426  true : /*only descendants from primaries with MCTruth process == "primary"*/
427  false; /*not from MCTruth process "primary"*/
428  // Polarization.
429  const G4ThreeVector& polarization = track->GetPolarization();
430  fCurrentParticle.particle->SetPolarization(
431  TVector3{polarization.x(), polarization.y(), polarization.z()});
432  // Save the particle in the ParticleList.
433  fParticleList.Add(fCurrentParticle.particle);
434  }
435 
436  //----------------------------------------------------------------------------
437  void
438  ParticleListActionService::postUserTrackingAction(const G4Track* aTrack)
439  {
440  if (!fCurrentParticle.hasParticle()) return;
441 
442  if (aTrack) {
443  fCurrentParticle.particle->SetWeight(aTrack->GetWeight());
444 
445  // Get the post-step information from the G4Step.
446  const G4StepPoint* postStepPoint = aTrack->GetStep()->GetPostStepPoint();
447  if (!postStepPoint->GetProcessDefinedStep()) {
448  // Now we get to do some awkward cleanup because the
449  // fparticleList was augmented during the
450  // preUserTrackingAction. We cannot call 'Archive' because
451  // that only sets the mapped type of the entry to
452  // nullptr...which is really bad whenever we iterate through
453  // entries in the endOfEventAction and dereference the mapped
454  // type. We have to entirely erase the entry.
455  auto key_to_erase = fParticleList.key(fCurrentParticle.particle);
456  fParticleList.erase(key_to_erase);
457  // after the particle is archived, it is deleted
458  fCurrentParticle.clear();
459  return;
460  }
461 
462  G4String process = postStepPoint->GetProcessDefinedStep()->GetProcessName();
463  fCurrentParticle.particle->SetEndProcess(process);
464 
465  // -- D.R. Store the final point only for particles that have not had intermediate trajectory
466  // points saved. This avoids double counting the final trajectory point for particles from
467  // generators with storable trajectory points.
468 
469  if (!fCurrentParticle.keepFullTrajectory) {
470  const G4ThreeVector position = postStepPoint->GetPosition();
471  G4double time = postStepPoint->GetGlobalTime();
472 
473  // Remember that LArSoft uses cm, ns, GeV.
474  TLorentzVector fourPos(position.x() / CLHEP::cm,
475  position.y() / CLHEP::cm,
476  position.z() / CLHEP::cm,
477  time / CLHEP::ns);
478 
479  const G4ThreeVector momentum = postStepPoint->GetMomentum();
480  const G4double energy = postStepPoint->GetTotalEnergy();
481  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
482  momentum.y() / CLHEP::GeV,
483  momentum.z() / CLHEP::GeV,
484  energy / CLHEP::GeV);
485 
486  // Add another point in the trajectory.
487  AddPointToCurrentParticle(fourPos, fourMom, std::string(process));
488  }
489  // -- particle has a full trajectory, apply SparsifyTrajectory method if enabled
490  else if (fSparsifyTrajectories) {
491  fCurrentParticle.particle->SparsifyTrajectory(fSparsifyMargin, fKeepSecondToLast);
492  }
493  }
494 
495  // store truth record pointer, only if it is available
496  if (fCurrentParticle.isPrimary()) {
497  fPrimaryTruthMap[fCurrentParticle.particle->TrackId()] = fCurrentParticle.truthInfoIndex();
498  }
499 
500  return;
501  }
502 
503  //----------------------------------------------------------------------------
504  // With every step, add to the particle's trajectory.
505  void
506  ParticleListActionService::userSteppingAction(const G4Step* step)
507  {
508  // N.B. G4 guarantees that following are non-null:
509  // - step
510  // - step->GetPostStepPoint()
511  if (!fCurrentParticle.hasParticle() || !step->GetPostStepPoint()->GetProcessDefinedStep()) {
512  return;
513  }
514  // Temporary fix for problem where DeltaTime on the first step
515  // of optical photon propagation is calculated incorrectly. -wforeman
516  double const globalTime = step->GetTrack()->GetGlobalTime();
517  double const velocity_G4 = step->GetTrack()->GetVelocity();
518  double const velocity_step = step->GetStepLength() / step->GetDeltaTime();
519  if ((step->GetTrack()->GetDefinition()->GetPDGEncoding() == 0) &&
520  fabs(velocity_G4 - velocity_step) > 0.0001) {
521  // Subtract the faulty step time from the global time,
522  // and add the correct step time based on G4 velocity.
523  step->GetPostStepPoint()->SetGlobalTime(globalTime - step->GetDeltaTime() +
524  step->GetStepLength() / velocity_G4);
525  }
526 
527  // For the most part, we just want to add the post-step
528  // information to the particle's trajectory. There's one
529  // exception: In PreTrackingAction, the correct time information
530  // is not available. So add the correct vertex information here.
531 
532  if (fCurrentParticle.particle->NumberTrajectoryPoints() == 0) {
533 
534  // Get the pre/along-step information from the G4Step.
535  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
536 
537  const G4ThreeVector position = preStepPoint->GetPosition();
538  G4double time = preStepPoint->GetGlobalTime();
539 
540  // Remember that LArSoft uses cm, ns, GeV.
541  TLorentzVector fourPos(position.x() / CLHEP::cm,
542  position.y() / CLHEP::cm,
543  position.z() / CLHEP::cm,
544  time / CLHEP::ns);
545 
546  const G4ThreeVector momentum = preStepPoint->GetMomentum();
547  const G4double energy = preStepPoint->GetTotalEnergy();
548  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
549  momentum.y() / CLHEP::GeV,
550  momentum.z() / CLHEP::GeV,
551  energy / CLHEP::GeV);
552 
553  // Add the first point in the trajectory.
554  AddPointToCurrentParticle(fourPos, fourMom, "Start");
555  } // end if this is the first step
556 
557  // At this point, the particle is being transported through the
558  // simulation.
559  // change below:
560  // This method is being called for every step that
561  // the track passes through, but we don't want to update the
562  // trajectory information if the step was defined by the StepLimiter.
563  G4String process = step->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName();
564  G4bool ignoreProcess = process.contains("StepLimiter");
565 
566  // We store the initial creation point of the particle
567  // and its final position (ie where it has no more energy, or at least < 1 eV) no matter
568  // what, but whether we store the rest of the trajectory depends
569  // on the process, and on a user switch.
570  // -- D.R. Store additional trajectory points only for desired generators and processes
571  if (!ignoreProcess && fCurrentParticle.keepFullTrajectory) {
572 
573  // Get the post-step information from the G4Step.
574  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
575 
576  const G4ThreeVector position = postStepPoint->GetPosition();
577  G4double time = postStepPoint->GetGlobalTime();
578 
579  // Remember that LArSoft uses cm, ns, GeV.
580  TLorentzVector fourPos(position.x() / CLHEP::cm,
581  position.y() / CLHEP::cm,
582  position.z() / CLHEP::cm,
583  time / CLHEP::ns);
584 
585  const G4ThreeVector momentum = postStepPoint->GetMomentum();
586  const G4double energy = postStepPoint->GetTotalEnergy();
587  TLorentzVector fourMom(momentum.x() / CLHEP::GeV,
588  momentum.y() / CLHEP::GeV,
589  momentum.z() / CLHEP::GeV,
590  energy / CLHEP::GeV);
591 
592  // Add another point in the trajectory.
593  AddPointToCurrentParticle(fourPos, fourMom, std::string(process));
594  }
595  }
596 
597  //----------------------------------------------------------------------------
598  /// Utility class for the EndOfEventAction method: update the
599  /// daughter relationships in the particle list.
601  public:
602  explicit UpdateDaughterInformation(sim::ParticleList& p) : particleList{&p} {}
603  void
604  operator()(sim::ParticleList::value_type& particleListEntry) const
605  {
606  // We're looking at this Particle in the list.
607  int particleID = particleListEntry.first;
608 
609  // The parent ID of this particle;
610  // we ask the particle list since the particle itself might have been lost
611  // ("archived"), but the particle list still holds the information we need
612  int parentID = particleList->GetMotherOf(particleID);
613 
614  // If the parentID <= 0, this is a primary particle.
615  if (parentID <= 0) return;
616 
617  // If we get here, this particle is somebody's daughter. Add
618  // it to the list of daughter particles for that parent.
619 
620  // Get the parent particle from the list.
621  sim::ParticleList::iterator parentEntry = particleList->find(parentID);
622 
623  if (parentEntry == particleList->end()) {
624  // We have an "orphan": a particle whose parent isn't
625  // recorded in the particle list. This is not signficant;
626  // it's possible for a particle not to be saved in the list
627  // because it failed an energy cut, but for it to have a
628  // daughter that passed the cut (e.g., a nuclear decay).
629  return;
630  }
631  if (!parentEntry->second) return; // particle archived, nothing to update
632 
633  // Add the current particle to the daughter list of the parent.
634  simb::MCParticle* parent = parentEntry->second;
635  parent->AddDaughter(particleID);
636  }
637 
638  private:
639  sim::ParticleList* particleList;
640  };
641 
642  //----------------------------------------------------------------------------
643  simb::GeneratedParticleIndex_t
644  ParticleListActionService::GetPrimaryTruthIndex(int trackId) const
645  {
646  auto const it = fPrimaryTruthMap.find(trackId);
647  return it == fPrimaryTruthMap.end() ? simb::NoGeneratedParticleIndex : it->second;
648  }
649 
650  //----------------------------------------------------------------------------
651  // Yields the ParticleList accumulated during the current event.
652  sim::ParticleList&&
653  ParticleListActionService::YieldList()
654  {
655  // check if the ParticleNavigator has entries, and if
656  // so grab the highest track id value from it to
657  // add to the fTrackIDOffset
658  int highestID = 0;
659  for (auto pn = fParticleList.begin(); pn != fParticleList.end(); pn++)
660  if ((*pn).first > highestID) highestID = (*pn).first;
661 
662  //Only change the fTrackIDOffset if there is in fact a particle to add to the event
663  if ((fParticleList.size()) != 0) {
664  fTrackIDOffset = highestID + 1;
665  mf::LogDebug("YieldList:fTrackIDOffset")
666  << "highestID = " << highestID << "\nfTrackIDOffset= " << fTrackIDOffset;
667  }
668 
669  return std::move(fParticleList);
670  } // ParticleList&& ParticleListActionService::YieldList()
671 
672  //----------------------------------------------------------------------------
673  void
674  ParticleListActionService::AddPointToCurrentParticle(TLorentzVector const& pos,
675  TLorentzVector const& mom,
676  std::string const& process)
677  {
678  // Add the first point in the trajectory.
679  fCurrentParticle.particle->AddTrajectoryPoint(pos, mom, process, fKeepTransportation);
680  }
681 
682  // Called at the end of each event. Call detectors to convert hits for the
683  // event and pass the call on to the action objects.
684  void
685  ParticleListActionService::endOfEventAction(const G4Event*)
686  {
687  // -- End of Run Report
688  if (!fNotStoredCounterUMap.empty()) { // -- Only if there is something to report
689  std::stringstream sscounter;
690  sscounter << "Not Stored Process summary:";
691  for (auto const& [process, count] : fNotStoredCounterUMap) {
692  sscounter << "\n\t" << process << " : " << count;
693  }
694  mf::LogInfo("ParticleListActionService") << sscounter.str();
695  }
696 
697  partCol_ = std::make_unique<std::vector<simb::MCParticle>>();
698  droppedCol_ = std::make_unique<std::map<int,std::set<int>>>();
699  tpassn_ =
700  std::make_unique<art::Assns<simb::MCTruth, simb::MCParticle, sim::GeneratedParticleInfo>>();
701  // Set up the utility class for the "for_each" algorithm. (We only
702  // need a separate set-up for the utility class because we need to
703  // give it the pointer to the particle list. We're using the STL
704  // "for_each" instead of the C++ "for loop" because it's supposed
705  // to be faster.
706  std::for_each(
707  fParticleList.begin(), fParticleList.end(), UpdateDaughterInformation{fParticleList});
708 
709  MF_LOG_INFO("endOfEventAction") << "MCTruth Handles Size: " << fMCLists->size();
710 
711  unsigned int nGeneratedParticles = 0;
712  unsigned int nMCTruths = 0;
713  sim::ParticleList particleList = YieldList();
714  for (size_t mcl = 0; mcl < fMCLists->size(); ++mcl) {
715  auto const& mclistHandle = (*fMCLists)[mcl];
716  MF_LOG_INFO("endOfEventAction") << "mclistHandle Size: " << mclistHandle->size();
717  for (size_t m = 0; m < mclistHandle->size(); ++m) {
718  art::Ptr<simb::MCTruth> mct(mclistHandle, m);
719  MF_LOG_INFO("endOfEventAction") << "Found " << mct->NParticles() << " particles";
720  for (simb::MCParticle* p : particleList | ranges::views::values) {
721  auto gen_index = fMCTIndexMap[p->TrackId()];
722  if (gen_index != nMCTruths) continue;
723  assert(p->NumberTrajectoryPoints() != 0ull);
724  ++nGeneratedParticles;
725  sim::GeneratedParticleInfo const truthInfo{GetPrimaryTruthIndex(p->TrackId())};
726  if (!truthInfo.hasGeneratedParticleIndex() && p->Mother() == 0) {
727  MF_LOG_WARNING("endOfEventAction") << "No GeneratedParticleIndex()!";
728  // this means it's primary but with no information; logic error!!
729  throw art::Exception(art::errors::LogicError)
730  << "Failed to match primary particle:\n"
731  << "\nwith particles from the truth record '" << mclistHandle.provenance()->inputTag()
732  << "':\n\n";
733  }
734  partCol_->push_back(std::move(*p));
735  art::Ptr<simb::MCParticle> mcp_ptr{pid_, partCol_->size() - 1, productGetter_};
736  tpassn_->addSingle(mct, mcp_ptr, truthInfo);
737  }
738  mf::LogDebug("Offset") << "nGeneratedParticles = " << nGeneratedParticles;
739  for(auto const& dropped : fdroppedTracksMap )
740  {
741  droppedCol_->insert(dropped);
742  }
743  ++nMCTruths;
744  }
745  }
746  fTrackIDOffset = 0;
747  }
748 } // namespace LArG4
util::quantities::gigaelectronvolt GeV
fKeepTransportation(p.get< bool >("KeepTransportation", false))
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
Definition: UtilFunc.cxx:42
ParticleListActionService(fhicl::ParameterSet const &)
pdgs p
Definition: selectors.fcl:22
static constexpr bool
process_name use argoneut_mc_hitfinder track
double velocity_step
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
auto cend(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:579
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
fSparsifyMargin(p.get< double >("SparsifyMargin", 0.015))
fNotStoredPhysics(p.get< std::vector< std::string >>("NotStoredPhysics",{}))
static const int NoParticleId
Definition: sim.h:21
A value measured in the specified unit.
Definition: quantities.h:566
Use Geant4&#39;s user &quot;hooks&quot; to maintain a list of particles generated by Geant4.
MF_LOG_WARNING("SimWire")<< "SimWire is an example module that works for the "<< "MicroBooNE detector. Each experiment should implement "<< "its own version of this module to simulate electronics "<< "response."
double velocity_G4
fKeepEMShowerDaughters(p.get< bool >("keepEMShowerDaughters", true))
fSparsifyTrajectories(p.get< bool >("SparsifyTrajectories", false))
double globalTime
fKeepSecondToLast(p.get< bool >("KeepSecondToLast", false))
fkeepOnlyPrimaryFullTraj(p.get< bool >("keepOnlyPrimaryFullTrajectories", false))
Contains information about a generated particle.
float mass
Definition: dedx.py:47
UpdateDaughterInformation(sim::ParticleList &p)
std::size_t count(Cont const &cont)
void operator()(sim::ParticleList::value_type &particleListEntry) const
Tools and modules for checking out the basics of the Monte Carlo.
BEGIN_PROLOG don t mess with this pandoraTrackGausCryoW true
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG process_name