All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AnalysisExample_module.cc
Go to the documentation of this file.
1 /**
2  * @file AnalysisExample_module.cc
3  * @brief A basic "skeleton" to read in art::Event records from a file,
4  * access their information, and do something with them.
5  * @ingroup AnalysisExample
6  * @author William Seligman (seligman@nevis.columbia.edu)
7  *
8  * See
9  * <https://cdcvs.fnal.gov/redmine/projects/larsoft/wiki/Using_art_in_LArSoft>
10  * for a description of the ART classes used here.
11  *
12  * Almost everything you see in the code below may have to be changed
13  * by you to suit your task. The example task is to make histograms
14  * and n-tuples related to @f$ dE/dx @f$ of particle tracks in the detector.
15  *
16  * As you try to understand why things are done a certain way in this
17  * example ("What's all this stuff about 'auto const&'?"), it will help
18  * to read ADDITIONAL_NOTES.md in the same directory as this file.
19  *
20  * Also note that, despite our efforts, the documentation and the practices in
21  * this code may fall out of date. In doubt, ask!
22  *
23  * The last revision of this code was done in August 2017 with LArSoft
24  * v06_44_00.
25  *
26  * This is in-source documentation in Doxygen format. Doxygen is a
27  * utility that creates web pages from specially-formatted comments in
28  * program code. If your package is ever added to an official LArSoft
29  * release, the doxygen-style comments will be added to the LArSoft
30  * doxygen pages at <http://nusoft.fnal.gov/larsoft/doxsvn/html/>.
31  *
32  * You can see the doxygenated version of the following lines at
33  * <http://nusoft.fnal.gov/larsoft/doxsvn/html/classlar_1_1example_1_1AnalysisExample.html>
34  *
35  * Doxygen tips:
36  *
37  * In general, comments that start with "///" or in C++ comment blocks
38  * (like this one) will appear in the web pages created by Doxygen.
39  *
40  * A comment starting with "///" will be associated wth the next code
41  * line that follows it, while one starting with "///<" will be
42  * associated with the one above it.
43  *
44  * For more on doxygen, see the documentation at
45  * <http://www.stack.nl/~dimitri/doxygen/manual/index.html>
46  *
47  */
48 
49 // Always include headers defining everything you use. Starting from
50 // LArSoft and going up the software layers (nusimdata, art, etc.)
51 // ending with C++ is standard.
52 
53 // LArSoft includes
63 #include "nusimdata/SimulationBase/MCParticle.h"
64 #include "nusimdata/SimulationBase/MCTruth.h"
65 
66 // Framework includes
67 #include "art/Framework/Core/EDAnalyzer.h"
68 #include "art/Framework/Core/ModuleMacros.h"
69 #include "art/Framework/Principal/Event.h"
70 #include "art/Framework/Principal/Handle.h"
71 #include "art/Framework/Services/Registry/ServiceHandle.h"
72 #include "art_root_io/TFileService.h"
73 #include "canvas/Persistency/Common/FindManyP.h"
74 #include "canvas/Utilities/Exception.h"
75 
76 // Utility libraries
77 #include "cetlib/pow.h" // cet::sum_of_squares()
78 #include "fhiclcpp/types/Atom.h"
79 #include "fhiclcpp/types/Table.h"
80 #include "messagefacility/MessageLogger/MessageLogger.h"
81 
82 // ROOT includes. Note: To look up the properties of the ROOT classes,
83 // use the ROOT web site; e.g.,
84 // <https://root.cern.ch/doc/master/annotated.html>
85 #include "TH1.h"
86 #include "TLorentzVector.h"
87 #include "TTree.h"
88 #include "TVector3.h"
89 
90 // C++ includes
91 #include <cmath>
92 #include <map>
93 
94 namespace {
95 
96  // This is a local namespace (as opposed to the one below, which has
97  // the nested name lar::example::).
98  //
99  // The stuff you declare in an local namespace will not be
100  // visible beyond this file (more technically, this "translation
101  // unit", which is everything that gets included in the .cc file from
102  // now on). In this way, any functions you define in this namespace
103  // won't affect the environment of other modules.
104 
105  // We will define functions at the end, but we declare them here so
106  // that the module can freely use them.
107 
108  /// Utility function to get the diagonal of the detector
109  /// @ingroup AnalysisExample
110  double DetectorDiagonal(geo::GeometryCore const& geom);
111 
112  /// Comparison routine for using std::lower/upper_bound to search
113  /// TDCIDE vectors.
114  /// @ingroup AnalysisExample
115  bool TDCIDETimeCompare(const sim::TDCIDE&, const sim::TDCIDE&);
116 
117 } // local namespace
118 
119 // It is good programming practice to put your code into a namespace,
120 // to make sure that no method or function you define will conflict
121 // with anything in any other module with the same name. Here we
122 // follow the LArSoft standard of declaring a main namespace of "lar"
123 // with a nested namespace of "example" because we're in the
124 // larexamples product. If an outside package wanted to call this
125 // module, it would have to use something like
126 // lar::example::AnalysisExample.
127 
128 namespace lar {
129  namespace example {
130 
131  // BEGIN AnalysisExample group
132  // -----------------------------------------------
133  /// @ingroup AnalysisExample
134  /// @{
135  //-----------------------------------------------------------------------
136  //-----------------------------------------------------------------------
137  // class definition
138  /**
139  * @brief Example analyzer module.
140  * @see @ref AnalysisExample "analysis module example overview"
141  *
142  * This class extracts information from the generated and reconstructed
143  * particles.
144  *
145  * It produces histograms for the simulated particles in the input file:
146  * - PDG ID (flavor) of all particles
147  * - momentum of the primary particles selected to have a specific PDG ID
148  * - length of the selected particle trajectory
149  *
150  * It also produces two ROOT trees.
151  *
152  * The first ROOT tree contains information on the simulated
153  * particles, including "dEdx", a binned histogram of collected
154  * charge as function of track range.
155  *
156  * The second ROOT tree contains information on the reconstructed hits.
157  *
158  * Configuration parameters
159  * =========================
160  *
161  * - *SimulationLabel* (string, default: "largeant"): tag of the input data
162  * product with the detector simulation information (typically an instance
163  * of the LArG4 module)
164  *
165  * - *HitLabel* (string, mandatory): tag of the input data product with
166  * reconstructed hits
167  *
168  * - *ClusterLabel* (string, mandatory): tag of the input data product with
169  * reconstructed clusters
170  *
171  * - *PDGcode* (integer, mandatory): particle type (PDG ID) to be selected;
172  * only primary particles of this type will be considered
173  *
174  * - *BinSize* (real, mandatory): @f$ dx @f$ [cm] used for the @f$ dE/dx @f$
175  * calculation
176  *
177  */
178  class AnalysisExample : public art::EDAnalyzer {
179  public:
180  // This structure describes the configuration parameters of the
181  // module. Any missing or unknown parameters will generate a
182  // configuration error.
183  //
184  // With an additional trick (see below) it allows configuration
185  // documentation to be displayed as a command-line option (see
186  // below).
187  //
188  // Note that, in this example, the Name() string (that is, the
189  // name you call a parameter in the FHiCL configuration file) and
190  // the data member name in the Config struct (that is, the name
191  // you access that parameter in your C++ code) match. This is not
192  // required, but it makes it easier to remember them.
193  //
194  // More details at:
195  // https://cdcvs.fnal.gov/redmine/projects/fhicl-cpp/wiki/Configuration_validation_and_fhiclcpp_types
196  //
197  struct Config {
198 
199  // Save some typing:
200  using Name = fhicl::Name;
201  using Comment = fhicl::Comment;
202 
203  // One Atom for each parameter
204  fhicl::Atom<art::InputTag> SimulationLabel{
205  Name("SimulationLabel"),
206  Comment("tag of the input data product with the detector simulation "
207  "information")};
208 
209  fhicl::Atom<art::InputTag> HitLabel{
210  Name("HitLabel"),
211  Comment("tag of the input data product with reconstructed hits")};
212 
213  fhicl::Atom<art::InputTag> ClusterLabel{
214  Name("ClusterLabel"),
215  Comment("tag of the input data product with reconstructed clusters")};
216 
217  fhicl::Atom<int> PDGcode{
218  Name("PDGcode"),
219  Comment("particle type (PDG ID) of the primary particle to be selected")};
220 
221  fhicl::Atom<double> BinSize{Name("BinSize"),
222  Comment("dx [cm] used for the dE/dx calculation")};
223 
224  }; // Config
225 
226  // If we define "Parameters" in this way, art will know to use the
227  // "Comment" items above for the module description. See what
228  // this command does:
229  //
230  // lar --print-description AnalysisExample
231  //
232  // The details of the art side of this trick are at:
233  //
234  // https://cdcvs.fnal.gov/redmine/projects/art/wiki/Configuration_validation_and_description
235  //
236  using Parameters = art::EDAnalyzer::Table<Config>;
237 
238  // -------------------------------------------------------------------
239  // -------------------------------------------------------------------
240  // Standard constructor for an ART module with configuration validation;
241  // we don't need a special destructor here.
242 
243  /// Constructor: configures the module (see the Config structure above)
244  explicit AnalysisExample(Parameters const& config);
245 
246  // The following methods have a standard meaning and a standard signature
247  // inherited from the framework (art::EDAnalyzer class).
248 
249  // - The "virtual" keyword is a reminder that the function we are
250  // dealing with is, in fact, virtual. You don't need to
251  // understand it now, but it's very important when you write a
252  // new algorithm.
253 
254  // - The "override" keyword, new in C++ 2011, is an important
255  // safety measure that ensures that the method we are going to
256  // write will override a matching method in the base class. For
257  // example, if we mistyped it as
258 
259  // virtual void beginJob() const;
260 
261  // (adding "const"), the compiler will be very happy with it,
262  // but art will not touch it, because art needs a "void
263  // beginJob()" (non-const) and it will find one in the base
264  // class (void art::EDAnalyzer::beginJob()) and will silently
265  // use that one instead. If you accidentally use:
266 
267  // virtual void beginJob() const override;
268 
269  // the compiler will immediately complain with us that this
270  // method is overriding nothing, hinting to some mistake (the
271  // spurious "const" in this case).
272 
273  // This method is called once, at the start of the job. In this
274  // example, it will define the histograms and n-tuples we'll
275  // write.
276  virtual void beginJob() override;
277 
278  // This method is called once, at the start of each run. It's a
279  // good place to read databases or files that may have
280  // run-dependent information.
281  virtual void beginRun(const art::Run& run) override;
282 
283  // The analysis routine, called once per event.
284  virtual void analyze(const art::Event& event) override;
285 
286  private:
287  // The stuff below is the part you'll most likely have to change to
288  // go from this custom example to your own task.
289 
290  // The parameters we'll read from the .fcl file.
291  art::InputTag fSimulationProducerLabel; ///< The name of the producer that tracked
292  ///< simulated particles through the detector
293  art::InputTag fHitProducerLabel; ///< The name of the producer that created hits
294  art::InputTag fClusterProducerLabel; ///< The name of the producer that
295  ///< created clusters
296  int fSelectedPDG; ///< PDG code of particle we'll focus on
297  double fBinSize; ///< For dE/dx work: the value of dx.
298 
299  // Pointers to the histograms we'll create.
300  TH1D* fPDGCodeHist; ///< PDG code of all particles
301  TH1D* fMomentumHist; ///< momentum [GeV] of all selected particles
302  TH1D* fTrackLengthHist; ///< true length [cm] of all selected particles
303 
304  // The n-tuples we'll create.
305  TTree* fSimulationNtuple; ///< tuple with simulated data
306  TTree* fReconstructionNtuple; ///< tuple with reconstructed data
307 
308  // The comment lines with the @ symbols define groups in doxygen.
309  /// @name The variables that will go into both n-tuples.
310  /// @{
311  int fEvent; ///< number of the event being processed
312  int fRun; ///< number of the run being processed
313  int fSubRun; ///< number of the sub-run being processed
314  /// @}
315 
316  /// @name The variables that will go into the simulation n-tuple.
317  /// @{
318  int fSimPDG; ///< PDG ID of the particle being processed
319  int fSimTrackID; ///< GEANT ID of the particle being processed
320 
321  // Arrays for 4-vectors: (x,y,z,t) and (Px,Py,Pz,E).
322  // Note: old-style C++ arrays are considered obsolete. However,
323  // to create simple n-tuples, we still need to use them.
324  double fStartXYZT[4]; ///< (x,y,z,t) of the true start of the particle
325  double fEndXYZT[4]; ///< (x,y,z,t) of the true end of the particle
326  double fStartPE[4]; ///< (Px,Py,Pz,E) at the true start of the particle
327  double fEndPE[4]; ///< (Px,Py,Pz,E) at the true end of the particle
328 
329  /// Number of dE/dx bins in a given track.
331 
332  /// The vector that will be used to accumulate dE/dx values as a function
333  /// of range.
334  std::vector<double> fSimdEdxBins;
335  /// @}
336 
337  /// @name Variables used in the reconstruction n-tuple
338  /// @{
339  int fRecoPDG; ///< PDG ID of the particle being processed
340  int fRecoTrackID; ///< GEANT ID of the particle being processed
341 
342  /// Number of dE/dx bins in a given track.
344 
345  /// The vector that will be used to accumulate dE/dx values as a function
346  /// of range.
347  std::vector<double> fRecodEdxBins;
348 
349  /// @}
350 
351  // Other variables that will be shared between different methods.
352  geo::GeometryCore const* fGeometryService; ///< pointer to Geometry provider
353  double fElectronsToGeV; ///< conversion factor
354  int fTriggerOffset; ///< (units of ticks) time of expected neutrino event
355 
356  }; // class AnalysisExample
357 
358  /// @}
359  // END AnalysisExample group
360  // -------------------------------------------------
361 
362  //-----------------------------------------------------------------------
363  //-----------------------------------------------------------------------
364  // class implementation
365 
366  //-----------------------------------------------------------------------
367  // Constructor
368  //
369  // Note that config is a Table<Config>, and to access the Config
370  // value we need to use an operator: "config()". In the same way,
371  // each element in Config is an Atom<Type>, so to access the type we
372  // again use the call operator, e.g. "SimulationLabel()".
373  //
375  : EDAnalyzer(config)
376  , fSimulationProducerLabel(config().SimulationLabel())
377  , fHitProducerLabel(config().HitLabel())
378  , fClusterProducerLabel(config().ClusterLabel())
379  , fSelectedPDG(config().PDGcode())
380  , fBinSize(config().BinSize())
381  {
382  // Get a pointer to the geometry service provider.
383  fGeometryService = lar::providerFrom<geo::Geometry>();
384 
385  auto const clock_data =
386  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
387  fTriggerOffset = trigger_offset(clock_data);
388 
389  // Since art 2.8, you can and should tell beforehand, here in the
390  // constructor, all the data the module is going to read ("consumes") or
391  // might read
392  // ("may_consume"). Diligence here will in the future help the framework
393  // execute modules in parallel, making sure the order is correct.
394  consumes<std::vector<simb::MCParticle>>(fSimulationProducerLabel);
395  consumes<std::vector<sim::SimChannel>>(fSimulationProducerLabel);
396  consumes<art::Assns<simb::MCTruth, simb::MCParticle>>(fSimulationProducerLabel);
397  consumes<std::vector<recob::Hit>>(fHitProducerLabel);
398  consumes<std::vector<recob::Cluster>>(fClusterProducerLabel);
399  consumes<art::Assns<recob::Cluster, recob::Hit>>(fHitProducerLabel);
400  }
401 
402  //-----------------------------------------------------------------------
403  void
405  {
406  // Get the detector length, to determine the maximum bin edge of one
407  // of the histograms.
408  const double detectorLength = DetectorDiagonal(*fGeometryService);
409 
410  // Access ART's TFileService, which will handle creating and writing
411  // histograms and n-tuples for us.
412  art::ServiceHandle<art::TFileService const> tfs;
413 
414  // For TFileService, the arguments to 'make<whatever>' are the
415  // same as those passed to the 'whatever' constructor, provided
416  // 'whatever' is a ROOT class that TFileService recognizes.
417 
418  // Define the histograms. Putting semi-colons around the title
419  // causes it to be displayed as the x-axis label if the histogram
420  // is drawn (the format is "title;label on abscissae;label on ordinates").
421  fPDGCodeHist = tfs->make<TH1D>("pdgcodes", ";PDG Code;", 5000, -2500, 2500);
422  fMomentumHist = tfs->make<TH1D>("mom", ";particle Momentum (GeV);", 100, 0., 10.);
424  tfs->make<TH1D>("length", ";particle track length (cm);", 200, 0, detectorLength);
425 
426  // Define our n-tuples, which are limited forms of ROOT
427  // TTrees. Start with the TTree itself.
429  tfs->make<TTree>("AnalysisExampleSimulation", "AnalysisExampleSimulation");
431  tfs->make<TTree>("AnalysisExampleReconstruction", "AnalysisExampleReconstruction");
432 
433  // Define the branches (columns) of our simulation n-tuple. To
434  // write a variable, we give the address of the variable to
435  // TTree::Branch.
436  fSimulationNtuple->Branch("Event", &fEvent, "Event/I");
437  fSimulationNtuple->Branch("SubRun", &fSubRun, "SubRun/I");
438  fSimulationNtuple->Branch("Run", &fRun, "Run/I");
439  fSimulationNtuple->Branch("TrackID", &fSimTrackID, "TrackID/I");
440  fSimulationNtuple->Branch("PDG", &fSimPDG, "PDG/I");
441  // When we write arrays, we give the address of the array to
442  // TTree::Branch; in C++ this is simply the array name.
443  fSimulationNtuple->Branch("StartXYZT", fStartXYZT, "StartXYZT[4]/D");
444  fSimulationNtuple->Branch("EndXYZT", fEndXYZT, "EndXYZT[4]/D");
445  fSimulationNtuple->Branch("StartPE", fStartPE, "StartPE[4]/D");
446  fSimulationNtuple->Branch("EndPE", fEndPE, "EndPE[4]/D");
447  // For a variable-length array: include the number of bins.
448  fSimulationNtuple->Branch("NdEdx", &fSimNdEdxBins, "NdEdx/I");
449  // ROOT branches can contain std::vector objects.
450  fSimulationNtuple->Branch("dEdx", &fSimdEdxBins);
451 
452  // A similar definition for the reconstruction n-tuple. Note that we
453  // use some of the same variables in both n-tuples.
454  fReconstructionNtuple->Branch("Event", &fEvent, "Event/I");
455  fReconstructionNtuple->Branch("SubRun", &fSubRun, "SubRun/I");
456  fReconstructionNtuple->Branch("Run", &fRun, "Run/I");
457  fReconstructionNtuple->Branch("TrackID", &fRecoTrackID, "TrackID/I");
458  fReconstructionNtuple->Branch("PDG", &fRecoPDG, "PDG/I");
459  fReconstructionNtuple->Branch("NdEdx", &fRecoNdEdxBins, "NdEdx/I");
460  fReconstructionNtuple->Branch("dEdx", &fRecodEdxBins);
461  }
462 
463  //-----------------------------------------------------------------------
464  // art expects this function to have a art::Run argument; C++
465  // expects us to use all the arguments we are given, or it will
466  // generate an "unused variable" warning. But we don't actually need
467  // nor use the art::Run object in this example. The trick to prevent
468  // that warning is to omit (or comment out) the name of the
469  // parameter.
470 
471  void
472  AnalysisExample::beginRun(const art::Run& /*run*/)
473  {
474  // How to convert from number of electrons to GeV. The ultimate
475  // source of this conversion factor is
476  // ${LARCOREOBJ_INC}/larcoreobj/SimpleTypesAndConstants/PhysicalConstants.h.
477  // But sim::LArG4Parameters might in principle ask a database for it.
478  art::ServiceHandle<sim::LArG4Parameters const> larParameters;
479  fElectronsToGeV = 1. / larParameters->GeVToElectrons();
480  }
481 
482  //-----------------------------------------------------------------------
483  void
484  AnalysisExample::analyze(const art::Event& event)
485  {
486  // Start by fetching some basic event information for our n-tuple.
487  fEvent = event.id().event();
488  fRun = event.run();
489  fSubRun = event.subRun();
490 
491  // This is the standard method of reading multiple objects
492  // associated with the same event; see
493  // <https://cdcvs.fnal.gov/redmine/projects/larsoft/wiki/Using_art_in_LArSoft>
494  // for more information.
495  //
496  // Define a "handle" to point to a vector of the objects.
497  art::Handle<std::vector<simb::MCParticle>> particleHandle;
498 
499  // Then tell the event to fill the vector with all the objects of
500  // that type produced by a particular producer.
501  //
502  // Note that if I don't test whether this is successful, and there
503  // aren't any simb::MCParticle objects, then the first time we
504  // access particleHandle, art will display a "ProductNotFound"
505  // exception message and, depending on art settings, it may skip
506  // all processing for the rest of this event (including any
507  // subsequent analysis steps) or stop the execution.
508  if (!event.getByLabel(fSimulationProducerLabel, particleHandle)) {
509  // If we have no MCParticles at all in an event, then we're in
510  // big trouble. Throw an exception to force this module to
511  // stop. Try to include enough information for the user to
512  // figure out what's going on. Note that when we throw a
513  // cet::exception, the run and event number will automatically
514  // be displayed.
515  //
516  // __LINE__ and __FILE__ are values computed by the compiler.
517 
518  throw cet::exception("AnalysisExample")
519  << " No simb::MCParticle objects in this event - "
520  << " Line " << __LINE__ << " in file " << __FILE__ << std::endl;
521  }
522 
523  // Get all the simulated channels for the event. These channels
524  // include the energy deposited for each simulated track.
525  //
526  // Here we use a different method to access objects:
527  // art::ValidHandle. Using this method, if there aren't any
528  // objects of the given type (sim::SimChannel in this case) in the
529  // input file for this art::Event, art will throw a
530  // ProductNotFound exception.
531  //
532  // The "auto" type means that the C++ compiler will determine the
533  // appropriate type for us, based on the return type of
534  // art::Event::getValidHandle<T>(). The "auto" keyword is a great
535  // timesaver, especially with frameworks like LArSoft which often
536  // have complicated data product structures.
537 
538  auto simChannelHandle =
539  event.getValidHandle<std::vector<sim::SimChannel>>(fSimulationProducerLabel);
540 
541  //
542  // Let's compute the variables for the simulation n-tuple first.
543  //
544 
545  // The MCParticle objects are not necessarily in any particular
546  // order. Since we may have to search the list of particles, let's
547  // put them into a map, a sorted container that will make
548  // searching fast and easy. To save both space and time, the map
549  // will not contain a copy of the MCParticle, but a pointer to it.
550  std::map<int, const simb::MCParticle*> particleMap;
551 
552  // This is a "range-based for loop" in the 2011 version of C++; do
553  // a web search on "c++ range based for loop" for more
554  // information. Here's how it breaks down:
555 
556  // - A range-based for loop operates on a container.
557  // "particleHandle" is not a container; it's a pointer to a
558  // container. If we want C++ to "see" a container, we have to
559  // dereference the pointer, like this: *particleHandle.
560 
561  // - The loop variable that is set to each object in the container
562  // is named "particle". As for the loop variable's type:
563 
564  // - To save a little bit of typing, and to make the code easier
565  // to maintain, we're going to let the C++ compiler deduce the
566  // type of what's in the container (simb::MCParticle objects
567  // in this case), so we use "auto".
568 
569  // - We do _not_ want to change the contents of the container,
570  // so we use the "const" keyword to make sure.
571 
572  // - We don't need to copy each object from the container into
573  // the variable "particle". It's sufficient to refer to the
574  // object by its address. So we use the reference operator "&"
575  // to tell the compiler to just copy the address, not the
576  // entire object.
577 
578  // It sounds complicated, but the end result is that we loop over
579  // the list of particles in the art::Event in the most efficient
580  // way possible.
581 
582  for (auto const& particle : (*particleHandle)) {
583  // For the methods you can call to get particle information,
584  // see ${NUSIMDATA_INC}/nusimdata/SimulationBase/MCParticle.h.
585  fSimTrackID = particle.TrackId();
586 
587  // Add the address of the MCParticle to the map, with the
588  // track ID as the key.
589  particleMap[fSimTrackID] = &particle;
590 
591  // Histogram the PDG code of every particle in the event.
592  fSimPDG = particle.PdgCode();
593  fPDGCodeHist->Fill(fSimPDG);
594 
595  // For this example, we want to fill the n-tuples and histograms
596  // only with information from the primary particles in the
597  // event, whose PDG codes match a value supplied in the .fcl file.
598  if (particle.Process() != "primary" || fSimPDG != fSelectedPDG) continue;
599 
600  // A particle has a trajectory, consisting of a set of
601  // 4-positions and 4-mommenta.
602  const size_t numberTrajectoryPoints = particle.NumberTrajectoryPoints();
603 
604  // For trajectories, as for vectors and arrays, the first
605  // point is #0, not #1.
606  const int last = numberTrajectoryPoints - 1;
607  const TLorentzVector& positionStart = particle.Position(0);
608  const TLorentzVector& positionEnd = particle.Position(last);
609  const TLorentzVector& momentumStart = particle.Momentum(0);
610  const TLorentzVector& momentumEnd = particle.Momentum(last);
611 
612  // Make a histogram of the starting momentum.
613  fMomentumHist->Fill(momentumStart.P());
614 
615  // Fill arrays with the 4-values. (Don't be fooled by
616  // the name of the method; it just puts the numbers from
617  // the 4-vector into the array.)
618  positionStart.GetXYZT(fStartXYZT);
619  positionEnd.GetXYZT(fEndXYZT);
620  momentumStart.GetXYZT(fStartPE);
621  momentumEnd.GetXYZT(fEndPE);
622 
623  // Use a polar-coordinate view of the 4-vectors to
624  // get the track length.
625  const double trackLength = (positionEnd - positionStart).Rho();
626 
627  // Let's print some debug information in the job output to see
628  // that everything is fine. MF_LOG_DEBUG() is a messagefacility
629  // macro that prints stuff when the message level is set to
630  // standard_debug in the .fcl file.
631  MF_LOG_DEBUG("AnalysisExample") << "Track length: " << trackLength << " cm";
632 
633  // Fill a histogram of the track length.
634  fTrackLengthHist->Fill(trackLength);
635 
636  MF_LOG_DEBUG("AnalysisExample")
637  << "track ID=" << fSimTrackID << " (PDG ID: " << fSimPDG << ") " << trackLength
638  << " cm long, momentum " << momentumStart.P() << " GeV/c, has " << numberTrajectoryPoints
639  << " trajectory points";
640 
641  // Determine the number of dE/dx bins for the n-tuple.
642  fSimNdEdxBins = int(trackLength / fBinSize) + 1;
643  // Initialize the vector of dE/dx bins to be empty.
644  fSimdEdxBins.clear();
645 
646  // To look at the energy deposited by this particle's track,
647  // we loop over the SimChannel objects in the event.
648  for (auto const& channel : (*simChannelHandle)) {
649  // Get the numeric ID associated with this channel. (The
650  // channel number is a 32-bit unsigned int, which normally
651  // requires a special data type. Let's use "auto" so we
652  // don't have to remember "raw::ChannelID_t".
653  auto const channelNumber = channel.Channel();
654 
655  // A little care: There is more than one plane that reacts
656  // to charge in the TPC. We only want to include the
657  // energy from the collection plane. Note:
658  // geo::kCollection is defined in
659  // ${LARCOREOBJ_INC}/larcoreobj/SimpleTypesAndConstants/geo_types.h
660  if (fGeometryService->SignalType(channelNumber) != geo::kCollection) continue;
661 
662  // Each channel has a map inside it that connects a time
663  // slice to energy deposits in the detector. We'll use
664  // "auto", but it's worth noting that the full type of
665  // this map is
666  // std::map<unsigned short, std::vector<sim::IDE>>
667  auto const& timeSlices = channel.TDCIDEMap();
668 
669  // For every time slice in this channel:
670  for (auto const& timeSlice : timeSlices) {
671  // Each entry in a map is a pair<first,second>. For
672  // the timeSlices map, the 'first' is a time slice
673  // number. The 'second' is a vector of IDE objects.
674  auto const& energyDeposits = timeSlice.second;
675 
676  // Loop over the energy deposits. An "energy deposit"
677  // object is something that knows how much
678  // charge/energy was deposited in a small volume, by
679  // which particle, and where. The type of
680  // 'energyDeposit' will be sim::IDE, which is defined
681  // in
682  // ${LARDATAOBJ_INC}/lardataobj/Simulation/SimChannel.h.
683  for (auto const& energyDeposit : energyDeposits) {
684  // Check if the track that deposited the
685  // energy matches the track of the particle.
686  if (energyDeposit.trackID != fSimTrackID) continue;
687  // Get the (x,y,z) of the energy deposit.
688  TVector3 location(energyDeposit.x, energyDeposit.y, energyDeposit.z);
689 
690  // Distance from the start of the track.
691  const double distance = (location - positionStart.Vect()).Mag();
692 
693  // Into which bin of the dE/dx array do we add the energy?
694  const unsigned int bin = (unsigned int)(distance / fBinSize);
695 
696  // Is the dE/dx array big enough to include this bin?
697  if (fSimdEdxBins.size() < bin + 1) {
698  // Increase the array size to accomodate
699  // the new bin, padding it with zeros.
700  fSimdEdxBins.resize(bin + 1, 0.);
701  }
702 
703  // Add the energy deposit to that bin. (If you look at the
704  // definition of sim::IDE, you'll see that there's another
705  // way to get the energy. Are the two methods equivalent?
706  // Compare the results and see!)
707  fSimdEdxBins[bin] += energyDeposit.numElectrons * fElectronsToGeV;
708 
709  } // For each energy deposit
710  } // For each time slice
711  } // For each SimChannel
712 
713  // At this point we've filled in the values of all the
714  // variables and arrays we want to write to the n-tuple
715  // for this particle. The following command actually
716  // writes those values.
717  fSimulationNtuple->Fill();
718 
719  } // loop over all particles in the event.
720 
721  //
722  // Reconstruction n-tuple
723  //
724 
725  // All of the above is based on objects entirely produced by the
726  // simulation. Let's try to do something based on reconstructed
727  // objects. A Hit (see ${LARDATAOBJ_INC}/lardataobj/RecoBase/Hit.h)
728  // is a 2D object in a plane.
729 
730  // This code duplicates much of the code in the previous big
731  // simulation loop, and will produce the similar results. (It
732  // won't be identical, because of shaping and other factors; not
733  // every bit of charge in a channel ends up contributing to a
734  // hit.) The point is to show different methods of accessing
735  // information, not to produce profound physics results -- that
736  // part is up to you!
737 
738  // For the rest of this method, I'm going to assume you've read
739  // the comments in previous section; I won't repeat all the C++
740  // coding tricks and whatnot.
741 
742  // Start by reading the Hits. We don't use art::ValidHandle here
743  // because there might be no hits in the input; e.g., we ran the
744  // simulation but did not run the reconstruction steps yet. If
745  // there are no hits we may as well skip the rest of this module.
746 
747  art::Handle<std::vector<recob::Hit>> hitHandle;
748  if (!event.getByLabel(fHitProducerLabel, hitHandle)) return;
749 
750  // Our goal is to accumulate the dE/dx of any particles associated
751  // with the hits that match our criteria: primary particles with
752  // the PDG code from the .fcl file. I don't know how many such
753  // particles there will be in a given event. We'll use a map, with
754  // track ID as the key, to hold the vectors of dE/dx information.
755  std::map<int, std::vector<double>> dEdxMap;
756 
757  auto const clock_data =
758  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(event);
759 
760  // For every Hit:
761  for (auto const& hit : (*hitHandle)) {
762  // The channel associated with this hit.
763  auto hitChannelNumber = hit.Channel();
764 
765  // We have a hit. For this example let's just focus on the
766  // hits in the collection plane.
767  if (fGeometryService->SignalType(hitChannelNumber) != geo::kCollection) continue;
768 
769  MF_LOG_DEBUG("AnalysisExample") << "Hit in collection plane" << std::endl;
770 
771  // In a few lines we're going to look for possible energy
772  // deposits that correspond to that hit. Determine a
773  // reasonable range of times that might correspond to those
774  // energy deposits.
775 
776  // In reconstruction, the channel waveforms are truncated. So
777  // we have to adjust the Hit TDC ticks to match those of the
778  // SimChannels, which were created during simulation.
779 
780  // Save a bit of typing, while still allowing for potential
781  // changes in the definitions of types in
782  // $LARDATAOBJ_DIR/source/lardataobj/Simulation/SimChannel.h
783 
784  typedef sim::SimChannel::StoredTDC_t TDC_t;
785  TDC_t start_tdc = clock_data.TPCTick2TDC(hit.StartTick());
786  TDC_t end_tdc = clock_data.TPCTick2TDC(hit.EndTick());
787  TDC_t hitStart_tdc = clock_data.TPCTick2TDC(hit.PeakTime() - 3. * hit.SigmaPeakTime());
788  TDC_t hitEnd_tdc = clock_data.TPCTick2TDC(hit.PeakTime() + 3. * hit.SigmaPeakTime());
789 
790  start_tdc = std::max(start_tdc, hitStart_tdc);
791  end_tdc = std::min(end_tdc, hitEnd_tdc);
792 
793  // In the simulation section, we started with particles to find
794  // channels with a matching track ID. Now we search in reverse:
795  // search the SimChannels for matching channel number, then look
796  // at the tracks inside the channel.
797 
798  for (auto const& channel : (*simChannelHandle)) {
799  auto simChannelNumber = channel.Channel();
800  if (simChannelNumber != hitChannelNumber) continue;
801 
802  MF_LOG_DEBUG("AnalysisExample")
803  << "SimChannel number = " << simChannelNumber << std::endl;
804 
805  // The time slices in this channel.
806  auto const& timeSlices = channel.TDCIDEMap();
807 
808  // We want to look over the range of time slices in this
809  // channel that correspond to the range of hit times.
810 
811  // To do this, we're going to use some fast STL search
812  // methods; STL algorithms are usually faster than the
813  // ones you might write yourself. The price for this speed
814  // is a bit of code complexity: in particular, we need to
815  // a custom comparison function, TDCIDETimeCompare, to
816  // define a "less-than" function for the searches.
817 
818  // For an introduction to STL algorithms, see
819  // <http://www.learncpp.com/cpp-tutorial/16-4-stl-algorithms-overview/>.
820  // For a list of available STL algorithms, see
821  // <http://en.cppreference.com/w/cpp/algorithm>
822 
823  // We have to create "dummy" time slices for the search.
824  sim::TDCIDE startTime;
825  sim::TDCIDE endTime;
826  startTime.first = start_tdc;
827  endTime.first = end_tdc;
828 
829  // Here are the fast searches:
830 
831  // Find a pointer to the first channel with time >= start_tdc.
832  auto const startPointer =
833  std::lower_bound(timeSlices.begin(), timeSlices.end(), startTime, TDCIDETimeCompare);
834 
835  // From that time slice, find the last channel with time < end_tdc.
836  auto const endPointer =
837  std::upper_bound(startPointer, timeSlices.end(), endTime, TDCIDETimeCompare);
838 
839  // Did we find anything? If not, skip.
840  if (startPointer == timeSlices.end() || startPointer == endPointer) continue;
841  MF_LOG_DEBUG("AnalysisExample")
842  << "Time slice start = " << (*startPointer).first << std::endl;
843 
844  // Loop over the channel times we found that match the hit
845  // times.
846  for (auto slicePointer = startPointer; slicePointer != endPointer; slicePointer++) {
847  auto const timeSlice = *slicePointer;
848  auto time = timeSlice.first;
849 
850  // How to debug a problem: Lots of print statements. There are
851  // debuggers such as gdb, but they can be tricky to use with
852  // shared libraries and don't work if you're using software
853  // that was compiled somewhere else (e.g., you're accessing
854  // LArSoft libraries via CVMFS).
855 
856  // The MF_LOG_DEBUG statements below are left over from when I
857  // was trying to solve a problem about hit timing. I could
858  // have deleted them, but decided to leave them to demonsrate
859  // what a typical debugging process looks like.
860 
861  MF_LOG_DEBUG("AnalysisExample")
862  << "Hit index = " << hit.LocalIndex() << " channel number = " << hitChannelNumber
863  << " start TDC tick = " << hit.StartTick() << " end TDC tick = " << hit.EndTick()
864  << " peak TDC tick = " << hit.PeakTime()
865  << " sigma peak time = " << hit.SigmaPeakTime()
866  << " adjusted start TDC tick = " << clock_data.TPCTick2TDC(hit.StartTick())
867  << " adjusted end TDC tick = " << clock_data.TPCTick2TDC(hit.EndTick())
868  << " adjusted peak TDC tick = " << clock_data.TPCTick2TDC(hit.PeakTime())
869  << " adjusted start_tdc = " << start_tdc << " adjusted end_tdc = " << end_tdc
870  << " time = " << time << std::endl;
871 
872  // Loop over the energy deposits.
873  auto const& energyDeposits = timeSlice.second;
874  for (auto const& energyDeposit : energyDeposits) {
875  // Remember that map of MCParticles we created
876  // near the top of this method? Now we can use
877  // it. Search the map for the track ID associated
878  // with this energy deposit. Since a map is
879  // automatically sorted, we can use a fast binary
880  // search method, 'find()'.
881 
882  // By the way, the type of "search" is an iterator
883  // (to be specific, it's an
884  // std::map<int,simb::MCParticle*>::const_iterator,
885  // which makes you thankful for the "auto"
886  // keyword). If you're going to work with C++
887  // containers, you'll have to learn about
888  // iterators eventually; do a web search on "STL
889  // iterator" to get started.
890  auto search = particleMap.find(energyDeposit.trackID);
891 
892  // Did we find this track ID in the particle map?
893  // It's possible for the answer to be "no"; some
894  // particles are too low in kinetic energy to be
895  // written by the simulation (see
896  // ${LARSIM_DIR}/job/simulationservices.fcl,
897  // parameter ParticleKineticEnergyCut).
898  if (search == particleMap.end()) continue;
899 
900  // "search" points to a pair in the map: <track ID, MCParticle*>
901  int trackID = (*search).first;
902  const simb::MCParticle& particle = *((*search).second);
903 
904  // Is this a primary particle, with a PDG code that
905  // matches the user input?
906  if (particle.Process() != "primary" || particle.PdgCode() != fSelectedPDG) continue;
907 
908  // Determine the dE/dx of this particle.
909  const TLorentzVector& positionStart = particle.Position(0);
910  TVector3 location(energyDeposit.x, energyDeposit.y, energyDeposit.z);
911  double distance = (location - positionStart.Vect()).Mag();
912  unsigned int bin = int(distance / fBinSize);
913  double energy = energyDeposit.numElectrons * fElectronsToGeV;
914 
915  // A feature of maps: if we refer to
916  // dEdxMap[trackID], and there's no such entry in
917  // the map yet, it will be automatically created
918  // with a zero-size vector. Test to see if the
919  // vector for this track ID is big enough.
920  //
921  // dEdxMap is a map, which is a slow container
922  // compared to a vector. If we are going to access
923  // the same element over and over, it is a good
924  // idea to find that element once, and then refer
925  // to that item directly. Since we don't care
926  // about the specific type of dEdxMap[trackID] (a
927  // vector, by the way), we can use "auto" to save
928  // some time. This must be a reference, since we
929  // want to change the original value in the map,
930  // and can't be constant.
931  auto& track_dEdX = dEdxMap[trackID];
932  if (track_dEdX.size() < bin + 1) {
933  // Increase the vector size, padding it with
934  // zeroes.
935  track_dEdX.resize(bin + 1, 0);
936  }
937 
938  // Add the energy to the dE/dx bins for this track.
939  track_dEdX[bin] += energy;
940 
941  } // loop over energy deposits
942  } // loop over time slices
943  } // for each SimChannel
944  } // for each Hit
945 
946  // We have a map of dE/dx vectors. Write each one of them to the
947  // reconstruction n-tuple.
948  for (const auto& dEdxEntry : dEdxMap) {
949  // Here, the map entries are <first,second>=<track ID, dE/dx vector>
950  fRecoTrackID = dEdxEntry.first;
951 
952  // This is an example of how we'd pick out the PDG code if
953  // there are multiple particle types or tracks in a single
954  // event allowed in the n-tuple.
955  fRecoPDG = particleMap[fRecoTrackID]->PdgCode();
956 
957  // Get the number of bins for this track.
958  const std::vector<double>& dEdx = dEdxEntry.second;
959  fRecoNdEdxBins = dEdx.size();
960 
961  // Copy this track's dE/dx information.
963 
964  // At this point, we've filled in all the reconstruction
965  // n-tuple's variables. Write it.
966  fReconstructionNtuple->Fill();
967  }
968 
969  // Think about the two big loops above, One starts from the
970  // particles then looks at the channels; the other starts with the
971  // hits and backtracks to the particles. What links the
972  // information in simb::MCParticle and sim::SimChannel is the
973  // track ID number assigned by the LArG4 simulation; what links
974  // sim::SimChannel and recob::Hit is the channel ID.
975 
976  // In general, that is not how objects in the LArSoft
977  // reconstruction chain are linked together. Most of them are
978  // linked using associations and the art::Assns class:
979  // <https://cdcvs.fnal.gov/redmine/projects/larsoft/wiki/Using_art_in_LArSoft#artAssns>
980 
981  // The web page may a bit difficult to understand (at least, it is
982  // for me), so let's try to simplify things:
983 
984  // - If you're doing an analysis task, you don't have to worry
985  // about creating art::Assns objects.
986 
987  // - You don't have to read the art:Assns objects on your
988  // own. There are helper classes (FindXXX) which will do that for
989  // you.
990 
991  // - There's only one helper you need: art::FindManyP. It will do
992  // what you want with a minimum of fuss.
993 
994  // - Associations are symmetric. If you see an
995  // art::Assns<ClassA,ClassB>, the helper classes will find all of
996  // ClassA associated with ClassB or all of ClassB associated with
997  // ClassA.
998 
999  // - To know what associations exist, you have to be a 'code
1000  // detective'. The important clue is to look for a 'produces' line
1001  // in the code that writes objects to an art::Event. For example,
1002  // in ${LARSIM_DIR}/source/larsim/LArG4/LArG4_module.cc, you'll see this
1003  // line:
1004 
1005  // produces< art::Assns<simb::MCTruth, simb::MCParticle> >();
1006 
1007  // That means a simulated event will have an association between
1008  // simb::MCTruth (the primary particle produced by the event
1009  // generator) and the simb::MCParticle objects (the secondary
1010  // particles produced in the detector simulation).
1011 
1012  // Let's try it. The following statement will find the
1013  // simb::MCTruth objects associated with the simb::MCParticle
1014  // objects in the event (referenced by particleHandle):
1015 
1016  const art::FindManyP<simb::MCTruth> findManyTruth(
1017  particleHandle, event, fSimulationProducerLabel);
1018 
1019  // Note that we still have to supply the module label of the step
1020  // that created the association. Also note that we did not have to
1021  // explicitly read in the simb::MCTruth objects from the
1022  // art::Event object 'event'; FindManyP did that for us.
1023 
1024  // Also note that at this point art::FindManyP has already found
1025  // all the simb::MCTruth associated with each of the particles in
1026  // particleHandle. This is a slow process, so in general you want
1027  // to do it only once. If we had a loop over the particles, we
1028  // would still do this outside that loop.
1029 
1030  // Now we can query the 'findManyTruth' object to access the
1031  // information. First, check that there wasn't some kind of error:
1032 
1033  if (!findManyTruth.isValid()) {
1034  mf::LogError("AnalysisExample")
1035  << "findManyTruth simb::MCTruth for simb::MCParticle failed!";
1036  }
1037 
1038  // I'm going to be lazy, and just look at the simb::MCTruth object
1039  // associated with the first simb::MCParticle we read in. (The
1040  // main reason I'm being lazy is that if I used the
1041  // single-particle generator in prodsingle.fcl, every particle in
1042  // the event is going to be associated with just the one primary
1043  // particle from the event generator.)
1044 
1045  size_t particle_index = 0; // look at first particle in
1046  // particleHandle's vector.
1047 
1048  // I'm using "auto" to save on typing. The result of
1049  // FindManyP::at() is a vector of pointers, in this case
1050  // simb::MCTruth*. In this case it will be a vector with just one
1051  // entry; I could have used art::FindOneP instead. (This will be a
1052  // vector of art::Ptr, which is a type of smart pointer; see
1053  // <https://cdcvs.fnal.gov/redmine/projects/larsoft/wiki/Using_art_in_LArSoft#artPtrltTgt-and-artPtrVectorltTgt>
1054  // To avoid unnecessary copying, and since art::FindManyP returns
1055  // a constant reference, use "auto const&".
1056 
1057  auto const& truth = findManyTruth.at(particle_index);
1058 
1059  // Make sure there's no problem.
1060  if (truth.empty()) {
1061  mf::LogError("AnalysisExample")
1062  << "Particle ID=" << particleHandle->at(particle_index).TrackId() << " has no primary!";
1063  }
1064 
1065  // Use the message facility to write output. I don't have to write
1066  // the event, run, or subrun numbers; the message facility takes
1067  // care of that automatically. I'm "going at warp speed" with the
1068  // vectors, pointers, and methods; see
1069  // ${NUSIMDATA_INC}/nusimdata/SimulationBase/MCTruth.h and
1070  // ${NUSIMDATA_INC}/nusimdata/SimulationBase/MCParticle.h for the
1071  // nested calls I'm using.
1072  mf::LogInfo("AnalysisExample")
1073  << "Particle ID=" << particleHandle->at(particle_index).TrackId()
1074  << " primary PDG code=" << truth[0]->GetParticle(0).PdgCode();
1075 
1076  // Let's try a slightly more realistic example. Suppose I want to
1077  // read in the clusters, and learn what hits are associated with
1078  // them. Then I could backtrack from those hits to determine the
1079  // dE/dx of the particles in the clusters. (Don't worry; I won't
1080  // put you through the n-tuple creation procedure for a third
1081  // time.)
1082 
1083  // First, read in the clusters (if there are any).
1084  art::Handle<std::vector<recob::Cluster>> clusterHandle;
1085  if (!event.getByLabel(fClusterProducerLabel, clusterHandle)) return;
1086 
1087  // Now use the associations to find which hits are associated with
1088  // which clusters. Note that this is not as trivial a query as the
1089  // one with MCTruth, since it's possible for a hit to be assigned
1090  // to more than one cluster.
1091 
1092  // We have to include fClusterProducerLabel, since that's the step
1093  // that created the art::Assns<recob::Hit,recob::Cluster> object;
1094  // look at the modules in ${LARRECO_DIR}/source/larreco/ClusterFinder/
1095  // and search for the 'produces' lines. (I did not know this before I
1096  // wrote these comments. I had to be a code detective and use UNIX
1097  // tools like 'grep' and 'find' to locate those routines.)
1098  const art::FindManyP<recob::Hit> findManyHits(clusterHandle, event, fClusterProducerLabel);
1099 
1100  if (!findManyHits.isValid()) {
1101  mf::LogError("AnalysisExample") << "findManyHits recob::Hit for recob::Cluster failed;"
1102  << " cluster label='" << fClusterProducerLabel << "'";
1103  }
1104 
1105  // Now we'll loop over the clusters to see the hits associated
1106  // with each one. Note that we're not using a range-based for
1107  // loop. That's because FindManyP::at() expects a number as an
1108  // argument, so we might as well go through the cluster objects
1109  // via numerical index instead.
1110  for (size_t cluster_index = 0; cluster_index != clusterHandle->size(); cluster_index++) {
1111  // In this case, FindManyP::at() will return a vector of
1112  // pointers to recob::Hit that corresponds to the
1113  // "cluster_index"-th cluster.
1114  auto const& hits = findManyHits.at(cluster_index);
1115 
1116  // We have a vector of pointers to the hits associated
1117  // with the cluster, but for this example I'm not going to
1118  // do anything fancy with them. I'll just print out how
1119  // many there are.
1120  mf::LogInfo("AnalysisExample") << "Cluster ID=" << clusterHandle->at(cluster_index).ID()
1121  << " has " << hits.size() << " hits";
1122  }
1123 
1124  } // AnalysisExample::analyze()
1125 
1126  // This macro has to be defined for this module to be invoked from a
1127  // .fcl file; see AnalysisExample.fcl for more information.
1128  DEFINE_ART_MODULE(AnalysisExample)
1129 
1130  } // namespace example
1131 } // namespace lar
1132 
1133 // Back to our local namespace.
1134 namespace {
1135 
1136  // Define a local function to calculate the detector diagonal.
1137  double
1138  DetectorDiagonal(geo::GeometryCore const& geom)
1139  {
1140  const double length = geom.DetLength();
1141  const double width = 2. * geom.DetHalfWidth();
1142  const double height = 2. * geom.DetHalfHeight();
1143 
1144  return std::sqrt(cet::sum_of_squares(length, width, height));
1145  } // DetectorDiagonal()
1146 
1147  // Define a comparison function to use in std::upper_bound and
1148  // std::lower_bound searches above.
1149  bool
1150  TDCIDETimeCompare(const sim::TDCIDE& lhs, const sim::TDCIDE& rhs)
1151  {
1152  return lhs.first < rhs.first;
1153  }
1154 } // local namespace
Store parameters for running LArG4.
int fSimNdEdxBins
Number of dE/dx bins in a given track.
TH1D * fPDGCodeHist
PDG code of all particles.
Utilities related to art service access.
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
fhicl::Atom< art::InputTag > SimulationLabel
int fRecoTrackID
GEANT ID of the particle being processed.
std::pair< unsigned short, std::vector< sim::IDE > > TDCIDE
List of energy deposits at the same time (on this channel)
Definition: SimChannel.h:127
Declaration of signal hit object.
art::InputTag fHitProducerLabel
The name of the producer that created hits.
double fEndXYZT[4]
(x,y,z,t) of the true end of the particle
geo::GeometryCore const * fGeometryService
pointer to Geometry provider
double fStartPE[4]
(Px,Py,Pz,E) at the true start of the particle
AnalysisExample(Parameters const &config)
Constructor: configures the module (see the Config structure above)
int fSimTrackID
GEANT ID of the particle being processed.
TTree * fReconstructionNtuple
tuple with reconstructed data
virtual void analyze(const art::Event &event) override
process_name hit
Definition: cheaterreco.fcl:51
SigType_t SignalType(geo::PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
int fRecoNdEdxBins
Number of dE/dx bins in a given track.
constexpr details::BinObj< T > bin(T value)
Returns a wrapper to print the specified data in binary format.
Access the description of detector geometry.
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
virtual void beginRun(const art::Run &run) override
art::EDAnalyzer::Table< Config > Parameters
double fStartXYZT[4]
(x,y,z,t) of the true start of the particle
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
BEGIN_PROLOG vertical distance to the surface Name
Description of geometry of one entire detector.
Declaration of cluster object.
int fRecoPDG
PDG ID of the particle being processed.
Definition of data types for geometry description.
float dEdx(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const TCSlice &slc, TP3D &tp3d)
Definition: PFPUtils.cxx:2687
int fSimPDG
PDG ID of the particle being processed.
int fRun
number of the run being processed
TH1D * fMomentumHist
momentum [GeV] of all selected particles
int fEvent
number of the event being processed
object containing MC truth information necessary for making RawDigits and doing back tracking ...
int trigger_offset(DetectorClocksData const &data)
double fEndPE[4]
(Px,Py,Pz,E) at the true end of the particle
TTree * fSimulationNtuple
tuple with simulated data
art::ServiceHandle< art::TFileService > tfs
double fElectronsToGeV
conversion factor
int fTriggerOffset
(units of ticks) time of expected neutrino event
double fBinSize
For dE/dx work: the value of dx.
art framework interface to geometry description
TH1D * fTrackLengthHist
true length [cm] of all selected particles
int fSelectedPDG
PDG code of particle we&#39;ll focus on.
TDCIDE::first_type StoredTDC_t
Type for TDC tick used in the internal representation.
Definition: SimChannel.h:149
Signal from collection planes.
Definition: geo_types.h:146