LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
MCDumpers.h
Go to the documentation of this file.
1 
11 #ifndef LARDATAALG_MCDUMPERS_MCDUMPERS_H
12 #define LARDATAALG_MCDUMPERS_MCDUMPERS_H
13 
14 // LArSoft libraries
16 #include "larcorealg/CoreUtils/DumpUtils.h" // lar::dump namespace
17 #include "larcorealg/Geometry/geo_vectors_utils_TVector.h" // geo::vect::dump
18 
19 // nutools libraries
24 
25 // ROOT libraries
26 #include "TLorentzVector.h"
27 
28 // C/C++ standard libraries
29 #include <string>
30 #include <utility> // std::forward()
31 
32 
33 namespace sim {
34 
36  namespace dump {
37 
38  //--------------------------------------------------------------------------
40 
53  template <typename Stream>
54  void DumpMCParticle(
55  Stream&& out, simb::MCParticle const& particle,
56  std::string indent, std::string firstIndent
57  );
58 
59  template <typename Stream>
60  void DumpMCParticle
61  (Stream&& out, simb::MCParticle const& particle, std::string indent = "")
62  { DumpMCParticle(std::forward<Stream>(out), particle, indent, indent); }
63 
65 
66 
67  //--------------------------------------------------------------------------
69 
86  template <typename Stream>
88  Stream&& out, simb::MCTrajectory const& trajectory,
89  unsigned int pointsPerLine, std::string indent
90  );
91 
92  template <typename Stream>
94  (Stream&& out, simb::MCTrajectory const& trajectory)
95  {
96  DumpMCParticleTrajectory(std::forward<Stream>(out), trajectory, 0U, "");
97  }
98 
99  // @}
100 
101 
102  //--------------------------------------------------------------------------
103  // @{
118  template <typename Stream>
119  void DumpMCNeutrino(
120  Stream&& out, simb::MCNeutrino const& neutrino,
121  std::string indent, std::string firstIndent
122  );
123 
124  template <typename Stream>
125  void DumpMCNeutrino
126  (Stream&& out, simb::MCNeutrino const& neutrino, std::string indent = "")
127  { DumpMCNeutrino(std::forward<Stream>(out), neutrino, indent, indent); }
128 
129  // @}
130 
131 
132  //--------------------------------------------------------------------------
133  // @{
156  template <typename Stream>
157  void DumpMCTruth(
158  Stream&& out, simb::MCTruth const& truth, unsigned int pointsPerLine,
159  std::string indent, std::string firstIndent
160  );
161 
162  template <typename Stream>
164  Stream&& out, simb::MCTruth const& truth, unsigned int pointsPerLine,
165  std::string indent = ""
166  )
167  {
169  (std::forward<Stream>(out), truth, pointsPerLine, indent, indent);
170  }
171 
172  template <typename Stream>
174  Stream&& out, simb::MCTruth const& truth,
175  std::string indent, std::string firstIndent
176  )
177  { DumpMCTruth(std::forward<Stream>(out), truth, 0, indent, firstIndent); }
178 
179  template <typename Stream>
181  Stream&& out, simb::MCTruth const& truth, std::string indent = ""
182  )
183  { DumpMCTruth(std::forward<Stream>(out), truth, indent, indent); }
184 
185  // @}
186 
187 
188  //--------------------------------------------------------------------------
189  // @{
204  template <typename Stream>
205  void DumpGTruth(
206  Stream&& out, simb::GTruth const& truth,
207  std::string indent, std::string firstIndent
208  );
209 
210  template <typename Stream>
211  void DumpGTruth
212  (Stream&& out, simb::GTruth const& truth, std::string indent = "")
213  { DumpGTruth(std::forward<Stream>(out), truth, indent, indent); }
214 
215  // @}
216 
217  //--------------------------------------------------------------------------
218 
219  } // namespace dump
220 
221 } // namespace sim
222 
223 
224 //------------------------------------------------------------------------------
225 //--- template implementation
226 //------------------------------------------------------------------------------
227 template <typename Stream>
229  Stream&& out, simb::MCParticle const& particle,
230  std::string indent, std::string firstIndent
231 ) {
232  out << firstIndent
233  << "ID=" << particle.TrackId() << ": " << ParticleName(particle.PdgCode())
234  << " mass=" << particle.Mass() << " GeV/c2 "
235  << " status=" << particle.StatusCode()
236  << " (" << ParticleStatusName(particle.StatusCode()) << ")"
237  ;
238  if (particle.Weight() != 1.0) out << " weight=" << particle.Weight();
239  if (particle.Rescatter()) {
240  out << " rescattered (" << particle.Rescatter()
241  << ") at vertex " << particle.GetGvtx();
242  }
243  out << "\n" << indent << "created via "
244  << (particle.Process().empty()? "magics": particle.Process());
245  if (particle.Mother() == 0) out << " by the gods";
246  else out << " from ID=" << particle.Mother();
247 
248  const unsigned int nDaughters = particle.NumberDaughters();
249  const unsigned int nPoints = particle.NumberTrajectoryPoints();
250  if (nPoints > 0) {
251  TLorentzVector const& start = particle.Position();
252  TLorentzVector const& start_mom = particle.Momentum();
253  out << " at " << start << " cm with momentum " << start_mom << " GeV/c";
254  }
255  if (particle.Polarization().Mag2() != 0.) {
256  out
257  << " with polarization " << lar::dump::vector3D(particle.Polarization());
258  }
259  // does this particle seem to stop? (by decay, or because of extended traj.)
260  if ((nPoints > 1) || (nDaughters > 0)) {
261  out << "\n" << indent << ((nDaughters > 0)? "ends": "stops") << " by "
262  << (particle.EndProcess().empty()? "magics": particle.EndProcess());
263  if (nPoints > 1) {
264  TLorentzVector const& stop = particle.EndPosition();
265  TLorentzVector const& stop_mom = particle.EndMomentum();
266  out << " at " << stop << " cm with momentum " << stop_mom << " GeV/c";
267  }
268  if (nDaughters > 0) {
269  out << " into ";
270  if (nDaughters == 1)
271  out << "particle ID=" << particle.FirstDaughter();
272  else {
273  out << nDaughters << " particles from ID=" << particle.FirstDaughter()
274  << " to ID=" << particle.LastDaughter();
275  }
276  } // if daughters
277  } // if it stops
278  if (nPoints > 1) {
279  simb::MCTrajectory const& traj = particle.Trajectory();
280  out << "\n" << indent << "comes with a trajectory " << traj.TotalLength()
281  << " cm long in " << nPoints << " points";
282  } // if has points
283 
284 } // sim::dump::DumpMCParticle()
285 
286 
287 //------------------------------------------------------------------------------
288 template <typename Stream>
290  Stream&& out, simb::MCTrajectory const& trajectory,
291  unsigned int pointsPerLine, std::string indent
292 ) {
293  unsigned int page = 0;
294  for (auto const& pair: trajectory) {
295  if ((pointsPerLine > 0) && (page-- == 0)) {
296  out << "\n" << indent << " ";
297  page = pointsPerLine - 1;
298  }
299  else out << " -- ";
300 
301  TLorentzVector const& pos = pair.first;
302  out << pos;
303  } // for trajectory points
304 
305 } // sim::dump::DumpMCParticleTrajectory()
306 
307 
308 //------------------------------------------------------------------------------
309 template <typename Stream>
311  Stream&& out, simb::MCNeutrino const& nu,
312  std::string indent, std::string firstIndent
313 ) {
314 
315  out << firstIndent
316  << "from " << TruthCCNCname(nu.CCNC())
318  << ", mode: " << nu.Mode() << " (" << TruthReactionMode(nu.Mode()) << ")"
319  << '\n' << indent
320  << "target: " << nu.Target() << " (" << ParticleName(nu.Target()) << ")"
321  ;
322  if (nu.HitNuc() != 0) {
323  out << ", hit nucleon: " << nu.HitNuc()
324  << " (" << ParticleName(nu.HitNuc()) << ")";
325  }
326  if (nu.HitQuark() != 0) {
327  out << ", hit quark: " << nu.HitQuark()
328  << " (" << ParticleName(nu.HitQuark()) << ")";
329  }
330  out
331  << '\n' << indent
332  << "x=" << nu.X() << " y=" << nu.Y() << " w=" << nu.W()
333  << " Q^2=" << nu.QSqr() << " GeV^2; theta=" << nu.Theta()
334  << " rad pT=" << nu.Pt() << " GeV/c"
335  ;
336  out << '\n' << indent << "neutrino: ";
337  DumpMCParticle(std::forward<Stream>(out), nu.Nu(), indent + " ", "");
338  out << '\n' << indent << "outgoing lepton: ";
339  DumpMCParticle(std::forward<Stream>(out), nu.Lepton(), indent + " ", "");
340 
341 } // sim::dump::DumpMCNeutrino()
342 
343 
344 //------------------------------------------------------------------------------
345 template <typename Stream>
347  Stream&& out, simb::MCTruth const& truth, unsigned int pointsPerLine,
348  std::string indent, std::string firstIndent
349 ) {
350  unsigned int const nParticles = truth.NParticles();
351  out << firstIndent
352  << nParticles << " particles from "
353  << TruthOriginName(truth.Origin());
354  if (truth.NeutrinoSet()) {
355  out << '\n' << indent << "neutrino information: ";
357  (std::forward<Stream>(out), truth.GetNeutrino(), indent + " ", "");
358  }
359  for (unsigned int i = 0; i < nParticles; ++i) {
360  out << '\n' << indent << "[#" << i << "] ";
361  simb::MCParticle const& particle = truth.GetParticle(i);
362  DumpMCParticle(std::forward<Stream>(out), particle, indent + " ", "");
363 
364  const unsigned int nPoints = particle.NumberTrajectoryPoints();
365  if ((nPoints > 0) && (pointsPerLine > 0)) {
366  out << ":";
368  std::forward<Stream>(out), particle.Trajectory(),
369  pointsPerLine, indent + " "
370  );
371  } // if has points
372  } // for all particles
373 
374 } // sim::dump::DumpMCTruth()
375 
376 
377 //------------------------------------------------------------------------------
378 template <typename Stream>
380  Stream&& out, simb::GTruth const& truth,
381  std::string indent, std::string firstIndent
382 ) {
383 
384  unsigned int const nCharged
385  = truth.fNumPiPlus + truth.fNumPiMinus + truth.fNumProton;
386  unsigned int const nNeutral = truth.fNumPi0 + truth.fNumNeutron;
387  unsigned int const nPions
388  = truth.fNumPiPlus + truth.fNumPiMinus + truth.fNumPi0;
389  unsigned int const nNucleons = truth.fNumProton + truth.fNumNeutron;
390  unsigned int const nTotalParticles = nCharged + nNeutral;
391 
392  out << firstIndent
393  << "interaction code: " << truth.fGint
394  << ", neutrino scattering code: " << truth.fGscatter
395  << " at " << truth.fVertex
396  << "\n" << indent
397  << "probe: " << ParticleName(truth.fProbePDG)
398  << " with cp=" << truth.fProbeP4
399  << " hit nucleon with cp=" << truth.fHitNucP4 << " GeV"
400  << " (" << (truth.fIsSeaQuark? "": "not a ") << "sea quark)"
401  << " in target: " << ParticleName(truth.ftgtPDG)
402  << " (Z: " << truth.ftgtZ << ", A: " << truth.ftgtA << ")"
403  << "\n" << indent
404  << "event interaction weight (genie internal): " << truth.fweight
405  << ", interaction probability: " << truth.fprobability
406  << ", cross section: " << truth.fXsec
407  << ", differential cross section: " << truth.fDiffXsec
408  << "\n" << indent
409  << "particles after reaction, before FSI: "
410  << truth.fNumPiPlus << " pi+"
411  << ", " << truth.fNumPiMinus << " pi-"
412  << ", " << truth.fNumPi0 << " pi0"
413  << ", " << truth.fNumProton << " p/pbar"
414  << ", " << truth.fNumNeutron << " n/nbar"
415  << "\n" << indent
416  << " total " << nTotalParticles << " particles after reaction before FSI"
417  ": " << nCharged << "/" << nNeutral << " charged/neutral"
418  ", " << nPions << " pions, " << nNucleons << " nucleons"
419  << "\n" << indent << "process "
420  << (truth.fIsCharm? "with": "without") << " charmed hadron";
421  if (truth.fResNum == -1) out << ", no resonance";
422  else out << ", resonance: #" << truth.fResNum;
423  out
424  << "\n" << indent
425  << "internal (on shell) genie kinematics: Q^2: " << truth.fgQ2 << " GeV^2"
426  << " q^2: " << truth.fgq2 << " GeV^2"
427  << ", w: " << truth.fgW << " GeV^2"
428  << ", t: " << truth.fgT << " GeV^2"
429  << ", x: " << truth.fgX
430  << ", y: " << truth.fgY
431  << "\n" << indent
432  << "FShadSyst: " << truth.fFShadSystP4
433  ;
434 
435 } // sim::DumpGTruth::DumpTruth()
436 
437 
438 //------------------------------------------------------------------------------
439 
440 
441 #endif // LARDATAALG_MCDUMPERS_MCDUMPERS_H
double fgW
Definition: GTruth.h:48
int fGint
interaction code
Definition: GTruth.h:25
unsigned int NumberTrajectoryPoints() const
Definition: MCParticle.h:222
const TVector3 & Polarization() const
Definition: MCParticle.h:218
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:223
int PdgCode() const
Definition: MCParticle.h:216
int CCNC() const
Definition: MCNeutrino.h:152
Specializations of geo_vectors_utils.h for ROOT old vector types.
double fgq2
Definition: GTruth.h:47
double QSqr() const
Definition: MCNeutrino.h:161
double Theta() const
angle between incoming and outgoing leptons, in radians
Definition: MCNeutrino.cxx:63
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:74
double fgX
Definition: GTruth.h:50
int ftgtA
Definition: GTruth.h:58
auto vector3D(Vector3D const &v)
Returns a manipulator which will print the specified vector.
Definition: DumpUtils.h:301
void DumpGTruth(Stream &&out, simb::GTruth const &truth, std::string indent, std::string firstIndent)
Dumps the content of the GENIE truth in the output stream.
Definition: MCDumpers.h:379
Trajectory class.
std::string TruthInteractionTypeName(int type)
const TLorentzVector & EndPosition() const
Definition: MCParticle.h:229
double Weight() const
Definition: MCParticle.h:258
int HitQuark() const
Definition: MCNeutrino.h:157
int FirstDaughter() const
Definition: MCParticle.h:254
const simb::MCTrajectory & Trajectory() const
Definition: MCParticle.h:257
int Mother() const
Definition: MCParticle.h:217
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:150
simb::Origin_t Origin() const
Definition: MCTruth.h:71
double Mass() const
Definition: MCParticle.h:243
double Pt() const
transverse momentum of interaction, in GeV/c
Definition: MCNeutrino.cxx:74
int fNumNeutron
number of neutrons after reaction, before FSI
Definition: GTruth.h:40
void DumpMCParticleTrajectory(Stream &&out, simb::MCTrajectory const &trajectory, unsigned int pointsPerLine, std::string indent)
Dumps the specified particle trajectory into the output stream.
Definition: MCDumpers.h:289
int HitNuc() const
Definition: MCNeutrino.h:156
int ftgtZ
Definition: GTruth.h:57
double fXsec
cross section of interaction
Definition: GTruth.h:32
int fNumPiPlus
number of pi pluses after reaction, before FSI
Definition: GTruth.h:36
int StatusCode() const
Definition: MCParticle.h:215
int fNumPiMinus
number of pi minuses after reaction, before FSI
Definition: GTruth.h:37
std::string TruthCCNCname(int ccnc)
int NParticles() const
Definition: MCTruth.h:72
std::string Process() const
Definition: MCParticle.h:219
TLorentzVector GetGvtx() const
Definition: MCParticle.h:249
Particle class.
int NumberDaughters() const
Definition: MCParticle.h:221
int TrackId() const
Definition: MCParticle.h:214
TLorentzVector fProbeP4
Definition: GTruth.h:63
int fResNum
resonance number
Definition: GTruth.h:42
int fNumProton
number of protons after reaction, before FSI
Definition: GTruth.h:39
int InteractionType() const
Definition: MCNeutrino.h:154
double fprobability
interaction probability
Definition: GTruth.h:31
std::string TruthReactionMode(int mode)
Returns the "mode" of the reaction (a lesser version of interaction type).
const simb::MCParticle & Lepton() const
Definition: MCNeutrino.h:151
Utilities to dump objects into a stream.
double W() const
Definition: MCNeutrino.h:158
int fProbePDG
Definition: GTruth.h:62
std::string EndProcess() const
Definition: MCParticle.h:220
double Y() const
Definition: MCNeutrino.h:160
int fGscatter
neutrino scattering code
Definition: GTruth.h:26
int fNumPi0
number of pi0 after reaction, before FSI
Definition: GTruth.h:38
std::string indent(std::size_t const i)
std::string ParticleName(int pigid)
Returns a string with the name of particle the specified with PDG ID.
void DumpMCTruth(Stream &&out, simb::MCTruth const &truth, unsigned int pointsPerLine, std::string indent, std::string firstIndent)
Dumps the content of the specified MC truth in the output stream.
Definition: MCDumpers.h:346
bool fIsCharm
did the interaction produce a charmed hadron?
Definition: GTruth.h:41
double fweight
event interaction weight (genie internal)
Definition: GTruth.h:30
double X() const
Definition: MCNeutrino.h:159
Monte Carlo Simulation.
TLorentzVector fHitNucP4
Definition: GTruth.h:56
const simb::MCParticle & GetParticle(int i) const
Definition: MCTruth.h:73
int Target() const
Definition: MCNeutrino.h:155
std::string TruthOriginName(simb::Origin_t origin)
Returns a string representing the specified process origin.
double TotalLength() const
int ftgtPDG
PDG of Target Nucleus, nucleon only if free.
Definition: GTruth.h:59
double fgQ2
< these are for the internal (on shell) genie kinematics
Definition: GTruth.h:46
Utility functions to print MC truth information.
int LastDaughter() const
Definition: MCParticle.h:255
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:224
TLorentzVector fFShadSystP4
Definition: GTruth.h:52
int Rescatter() const
Definition: MCParticle.h:256
double fgT
Definition: GTruth.h:49
void DumpMCNeutrino(Stream &&out, simb::MCNeutrino const &neutrino, std::string indent, std::string firstIndent)
Dumps the content of the specified neutrino in the output stream.
Definition: MCDumpers.h:310
bool NeutrinoSet() const
Definition: MCTruth.h:75
Event generator information.
Definition: MCTruth.h:30
Event generator information.
Definition: MCNeutrino.h:18
int Mode() const
Definition: MCNeutrino.h:153
bool fIsSeaQuark
Definition: GTruth.h:55
std::string ParticleStatusName(int code)
Describes the status of a particle (simb::MCParticle::StatusCode()).
TLorentzVector fVertex
Definition: GTruth.h:64
const TLorentzVector & EndMomentum() const
Definition: MCParticle.h:244
double fgY
Definition: GTruth.h:51
double fDiffXsec
differential cross section of interaction
Definition: GTruth.h:33
void DumpMCParticle(Stream &&out, simb::MCParticle const &particle, std::string indent, std::string firstIndent)
Dumps the content of the specified particle in the output stream.
Definition: MCDumpers.h:228