LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
ParticleListAction.h
Go to the documentation of this file.
1 
10 //
15 
16 #ifndef LArG4_ParticleListAction_h
17 #define LArG4_ParticleListAction_h
18 
19 #include "larsim/LArG4/ParticleFilters.h" // larg4::PositionInVolumeFilter
20 #include "lardataobj/Simulation/sim.h" // sim::GeneratorIndex_t, ...
21 #include "nutools/ParticleNavigation/ParticleList.h" // larg4::PositionInVolumeFilter
22 
24 #include "nusimdata/SimulationBase/simb.h" // simb::GeneratedParticleIndex_t
26 
27 #include "Geant4/globals.hh"
28 #include <map>
29 #include <limits> // std::numeric_limits<>
30 #include <cstddef> // std::size_t
31 
32 // Forward declarations.
33 class G4Event;
34 class G4Track;
35 class G4Step;
36 
37 namespace sim {
38  class ParticleList;
39 }
40 
41 namespace larg4 {
42 
44  {
45  public:
47 
48  struct ParticleInfo_t {
49 
50  simb::MCParticle* particle = nullptr;
51  bool keep = false;
54 
56  void clear()
57  {
58  particle = nullptr;
59  keep = false;
60  truthIndex = simb::NoGeneratedParticleIndex;
61  }
62 
64  bool hasParticle() const { return particle; }
65 
67  bool isPrimary() const { return simb::isGeneratedParticleIndex(truthIndex); }
68 
70  bool keepParticle() const { return hasParticle() && keep; }
71 
73  GeneratedParticleIndex_t truthInfoIndex() const { return truthIndex; }
74 
75  }; // ParticleInfo_t
76 
77  // Standard constructors and destructors;
78  ParticleListAction(double energyCut, bool storeTrajectories=false, bool keepEMShowerDaughters=false);
79  virtual ~ParticleListAction();
80 
81  // UserActions method that we'll override, to obtain access to
82  // Geant4's particle tracks and trajectories.
83  virtual void BeginOfEventAction(const G4Event*);
84  virtual void EndOfEventAction (const G4Event*);
85  virtual void PreTrackingAction (const G4Track*);
86  virtual void PostTrackingAction(const G4Track*);
87  virtual void SteppingAction (const G4Step* );
88 
90  void ParticleFilter(std::unique_ptr<PositionInVolumeFilter>&& filter)
91  { fFilter = std::move(filter); }
92 
93 
94  // TrackID of the current particle, EveID if the particle is from an EM shower
95  static int GetCurrentTrackID() { return fCurrentTrackID; }
96  static int GetCurrentPdgCode() { return fCurrentPdgCode; }
97 
98  void ResetTrackIDOffset() { fTrackIDOffset = 0; }
99 
100  // Returns the ParticleList accumulated during the current event.
101  const sim::ParticleList* GetList() const;
102 
105  std::map<int, GeneratedParticleIndex_t> const& GetPrimaryTruthMap() const
106  { return fPrimaryTruthMap; }
107 
109  GeneratedParticleIndex_t GetPrimaryTruthIndex(int trackId) const;
110 
111  // Yields the ParticleList accumulated during the current event.
112  sim::ParticleList&& YieldList();
113 
115  static bool isDropped(simb::MCParticle const* p);
116 
117  private:
118 
119  // this method will loop over the fParentIDMap to get the
120  // parentage of the provided trackid
121  int GetParentage(int trackid) const;
122 
123  G4double fenergyCut;
124  ParticleInfo_t fCurrentParticle;
126  sim::ParticleList* fparticleList;
128  G4bool fstoreTrajectories;
130  std::map<int, int> fParentIDMap;
131  static int fCurrentTrackID;
132  static int fCurrentPdgCode;
134  static int fTrackIDOffset;
135  bool fKeepEMShowerDaughters;
137 
138  std::unique_ptr<PositionInVolumeFilter> fFilter;
139 
141  std::map<int, GeneratedParticleIndex_t> fPrimaryTruthMap;
142 
144  void AddPointToCurrentParticle(TLorentzVector const& pos,
145  TLorentzVector const& mom,
146  std::string const& process);
147 
148  };
149 
150 } // namespace LArG4
151 
152 #endif // LArG4_ParticleListAction_h
std::unique_ptr< PositionInVolumeFilter > fFilter
filter for particles to be kept
std::map< int, GeneratedParticleIndex_t > fPrimaryTruthMap
Map: particle track ID -> index of primary information in MC truth.
bool keepParticle() const
Rerturns whether there is a particle known to be kept.
Geant4 interface.
Particle class.
GeneratedParticleIndex_t truthInfoIndex() const
Returns the index of the particle in the generator truth record.
std::map< int, int > fParentIDMap
key is current track ID, value is parent ID
Framework includes.
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)
see below
simb::GeneratedParticleIndex_t GeneratedParticleIndex_t
bool isGeneratedParticleIndex(GeneratedParticleIndex_t index)
Returns whether the specified one is an acceptable generator index.
Definition: simb.h:37
Monte Carlo Simulation.
bool isPrimary() const
Returns whether there is a particle.
void ParticleFilter(std::unique_ptr< PositionInVolumeFilter > &&filter)
Grabs a particle filter.
bool hasParticle() const
Returns whether there is a particle.
Defines classes to filter particles based on their trajectory.
Particle list in DetSim contains Monte Carlo particle information.
Tools and modules for checking out the basics of the Monte Carlo.
Common type definitions for data products (and a bit beyond).
std::map< int, GeneratedParticleIndex_t > const & GetPrimaryTruthMap() const
std::size_t GeneratedParticleIndex_t
Type of particle index in the generator truth record (simb::MCTruth).
Definition: simb.h:30