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