All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SingleGen_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file SingleGen_module.cc
3 /// \brief Generator for cosmic-rays
4 ///
5 /// Module designed to produce a set list of particles for a MC event
6 ///
7 /// \author brebel@fnal.gov
8 ////////////////////////////////////////////////////////////////////////
9 
10 // C++ includes.
11 #include <sstream>
12 #include <string>
13 #include <cmath>
14 #include <memory>
15 #include <iterator>
16 #include <map>
17 #include <initializer_list>
18 #include <cctype> // std::tolower()
19 
20 
21 // Framework includes
22 #include "art/Framework/Core/EDProducer.h"
23 #include "art/Framework/Principal/Event.h"
24 #include "art/Framework/Principal/Run.h"
25 #include "fhiclcpp/types/Name.h"
26 #include "fhiclcpp/types/Comment.h"
27 #include "fhiclcpp/types/Atom.h"
28 #include "fhiclcpp/types/OptionalAtom.h"
29 #include "fhiclcpp/types/Sequence.h"
30 #include "art/Framework/Services/Registry/ServiceHandle.h"
31 #include "art/Framework/Core/ModuleMacros.h"
32 #include "messagefacility/MessageLogger/MessageLogger.h"
33 #include "cetlib_except/exception.h"
34 #include "cetlib/filesystem.h"
35 #include "cetlib/search_path.h"
36 
37 // nurandom includes
38 #include "nurandom/RandomUtils/NuRandomService.h"
39 
40 // nusimdata includes
41 #include "nusimdata/SimulationBase/MCTruth.h"
42 #include "nusimdata/SimulationBase/MCParticle.h"
43 
44 // LArSoft libraries
49 
50 
51 #include "TVector3.h"
52 #include "TDatabasePDG.h"
53 #include "TFile.h"
54 #include "TH1.h"
55 #include "TH2.h"
56 
57 #include "CLHEP/Random/RandFlat.h"
58 #include "CLHEP/Random/RandGaussQ.h"
59 
60 namespace evgen {
61 
62  /// module to produce single or multiple specified particles in the detector
63  class SingleGen : public art::EDProducer {
64 
65  public:
66 
67  struct Config {
68  using Name = fhicl::Name;
69  using Comment = fhicl::Comment;
70 
71  fhicl::Atom<std::string> ParticleSelectionMode{
72  Name("ParticleSelectionMode"),
73  Comment("generate one particle, or one particle per PDG ID: " + presentOptions(ParticleSelectionModeNames))
74  };
75 
76  fhicl::Atom<bool> PadOutVectors{
77  Name("PadOutVectors"),
78  Comment("if true, all per-PDG-ID quantities must contain only one value, which is then used for all PDG IDs")
79  };
80 
81  fhicl::Sequence<int> PDG{
82  Name("PDG"),
83  Comment("PDG ID of the particles to be generated; this is the key for the other options marked as \"per PDG ID\"")
84  };
85 
86  fhicl::Atom<std::string> PDist{
87  Name("PDist"),
88  Comment("momentum distribution type: " + presentOptions(DistributionNames)),
90  };
91 
92  fhicl::Sequence<double> P0{
93  Name("P0"),
94  Comment("central momentum (GeV/c) to generate"),
95  [this](){ return !fromHistogram(PDist()); }
96  };
97 
98  fhicl::Sequence<double> SigmaP{
99  Name("SigmaP"),
100  Comment("variation in momenta (GeV/c)"),
101  [this](){ return !fromHistogram(PDist()); }
102  };
103 
104  fhicl::Sequence<double> X0{
105  Name("X0"),
106  Comment("central x position (cm) in world coordinates [per PDG ID]")
107  };
108 
109  fhicl::Sequence<double> Y0{
110  Name("Y0"),
111  Comment("central y position (cm) in world coordinates [per PDG ID]")
112  };
113 
114  fhicl::Sequence<double> Z0{
115  Name("Z0"),
116  Comment("central z position (cm) in world coordinates [per PDG ID]")
117  };
118 
119  fhicl::Sequence<double> T0{
120  Name("T0"),
121  Comment("central time (s) [per PDG ID]")
122  };
123 
124  fhicl::Sequence<double> SigmaX{
125  Name("SigmaX"),
126  Comment("variation (radius or RMS) in x position (cm) [per PDG ID]")
127  };
128 
129  fhicl::Sequence<double> SigmaY{
130  Name("SigmaY"),
131  Comment("variation (radius or RMS) in y position (cm) [per PDG ID]")
132  };
133 
134  fhicl::Sequence<double> SigmaZ{
135  Name("SigmaZ"),
136  Comment("variation (radius or RMS) in z position (cm) [per PDG ID]")
137  };
138 
139  fhicl::Sequence<double> SigmaT{
140  Name("SigmaT"),
141  Comment("variation (semi-interval or RMS) in time (s) [per PDG ID]")
142  };
143 
144  fhicl::Atom<std::string> PosDist{
145  Name("PosDist"),
146  Comment("distribution of starting position: " + presentOptions(DistributionNames, true, { kHIST }))
147  };
148 
149  fhicl::Atom<std::string> TDist{
150  Name("TDist"),
151  Comment("time distribution type: " + presentOptions(DistributionNames, true, { kHIST }))
152  };
153 
154  fhicl::Atom<bool> SingleVertex{
155  Name("SingleVertex"),
156  Comment("if true, all particles are produced at the same location"),
157  false
158  };
159 
160  fhicl::Sequence<double> Theta0XZ{
161  Name("Theta0XZ"),
162  Comment("angle from Z axis on world X-Z plane (degrees)")
163  };
164 
165  fhicl::Sequence<double> Theta0YZ{
166  Name("Theta0YZ"),
167  Comment("angle from Z axis on world Y-Z plane (degrees)")
168  };
169 
170  fhicl::Sequence<double> SigmaThetaXZ{
171  Name("SigmaThetaXZ"),
172  Comment("variation in angle in X-Z plane (degrees)")
173  };
174 
175  fhicl::Sequence<double> SigmaThetaYZ{
176  Name("SigmaThetaYZ"),
177  Comment("variation in angle in Y-Z plane (degrees)")
178  };
179 
180  fhicl::Atom<std::string> AngleDist{
181  Name("AngleDist"),
182  Comment("angular distribution type: " + presentOptions(DistributionNames)),
184  };
185 
186  fhicl::Atom<std::string> HistogramFile{
187  Name("HistogramFile"),
188  Comment("ROOT file containing the required distributions for the generation"),
189  [this](){ return fromHistogram(AngleDist()) || fromHistogram(PDist()); }
190  };
191 
192  fhicl::Sequence<std::string> PHist{
193  Name("PHist"),
194  Comment("name of the histograms of momentum distributions"),
195  [this](){ return fromHistogram(PDist()); }
196  };
197 
198  /*
199  fhicl::Sequence<std::string> ThetaPhiHist{
200  Name("ThetaPhiHist"),
201  Comment("name of the histograms of angular (theta/phi) distribution"),
202  [this](){ return fromHistogram(AngleDist()); }
203  };
204  */
205  fhicl::Sequence<std::string> ThetaXzYzHist{
206  Name("ThetaXzYzHist"),
207  Comment("name of the histograms of angular (X-Z and Y-Z) distribution"),
208  [this](){ return fromHistogram(AngleDist()); }
209  };
210 
211  fhicl::OptionalAtom<rndm::NuRandomService::seed_t> Seed{
212  Name("Seed"),
213  Comment("override the random number generator seed")
214  };
215 
216 
217  private:
218 
219  /// Returns whether the specified mode is an histogram distribution.
220  bool fromHistogram(std::string const& key) const;
221 
222  }; // struct Config
223 
224 
225  using Parameters = art::EDProducer::Table<Config>;
226 
227 
228  explicit SingleGen(Parameters const& config);
229 
230  // This is called for each event.
231  void produce(art::Event& evt) override;
232 
233  private:
234 
235  /// Names of all particle selection modes.
236  static const std::map<int, std::string> ParticleSelectionModeNames;
237  /// Names of all distribution modes.
238  static const std::map<int, std::string> DistributionNames;
239 
240  void SampleOne(unsigned int i,
241  simb::MCTruth &mct);
242  void SampleMany(simb::MCTruth &mct);
243  void Sample(simb::MCTruth &mct);
244  void printVecs(std::vector<std::string> const& list);
245  bool PadVector(std::vector<double> &vec);
246  double SelectFromHist(const TH1& h);
247  void SelectFromHist(const TH2& h, double &x, double &y);
248 
249  /// @{
250  /// @name Constants for particle type extraction mode (`ParticleSelectionMode` parameter).
251 
252  static constexpr int kSelectAllParts = 0; ///< One particle per entry is generated
253  static constexpr int kSelectOneRandPart = 1; ///< One particle is generated, extracted from the provided options.
254  /// @}
255 
256  /// @{
257  /// @name Constants for kinematic distribution options.
258 
259  static constexpr int kUNIF = 0; ///< Uniform distribution.
260  static constexpr int kGAUS = 1; ///< Gaussian distribution.
261  static constexpr int kHIST = 2; ///< Distribution from histograms.
262  /// @}
263 
264  int fMode; ///< Particle Selection Mode
265  ///< 0--generate a list of all particles,
266  ///< 1--generate a single particle selected randomly from the list
267  bool fPadOutVectors; ///< Select to pad out configuration vectors if they are not of
268  ///< of the same length as PDG
269  ///< false: don't pad out - all values need to specified
270  ///< true: pad out - default values assumed and printed out
271  std::vector<int> fPDG; ///< PDG code of particles to generate
272  std::vector<double> fP0; ///< Central momentum (GeV/c) to generate
273  std::vector<double> fSigmaP; ///< Variation in momenta (GeV/c)
274  int fPDist; ///< How to distribute momenta (gaus or uniform)
275  std::vector<double> fX0; ///< Central x position (cm) in world coordinates
276  std::vector<double> fY0; ///< Central y position (cm) in world coordinates
277  std::vector<double> fZ0; ///< Central z position (cm) in world coordinates
278  std::vector<double> fT0; ///< Central t position (s) in world coordinates
279  std::vector<double> fSigmaX; ///< Variation in x position (cm)
280  std::vector<double> fSigmaY; ///< Variation in y position (cm)
281  std::vector<double> fSigmaZ; ///< Variation in z position (cm)
282  std::vector<double> fSigmaT; ///< Variation in t position (s)
283  int fPosDist; ///< How to distribute xyz (gaus, or uniform)
284  int fTDist; ///< How to distribute t (gaus, or uniform)
285  bool fSingleVertex; ///< if true - all particles produced at the same location
286  std::vector<double> fTheta0XZ; ///< Angle in XZ plane (degrees)
287  std::vector<double> fTheta0YZ; ///< Angle in YZ plane (degrees)
288  std::vector<double> fSigmaThetaXZ; ///< Variation in angle in XZ plane
289  std::vector<double> fSigmaThetaYZ; ///< Variation in angle in YZ plane
290  int fAngleDist; ///< How to distribute angles (gaus, uniform)
291  std::string fHistFileName; ///< Filename containing histogram of momenta
292  std::vector<std::string> fPHist; ///< name of histogram of momenta
293  std::vector<std::string> fThetaXzYzHist; ///< name of histogram for thetaxz/thetayz distribution
294 
295  std::vector<std::unique_ptr<TH1>> hPHist ; /// actual TH1 for momentum distributions
296  std::vector<std::unique_ptr<TH2>> hThetaXzYzHist ; /// actual TH2 for angle distributions - Xz on x axis .
297  // FYI - thetaxz and thetayz are related to standard polar angles as follows:
298  // thetaxz = atan2(math.sin(theta) * cos(phi), cos(theta))
299  // thetayz = asin(sin(theta) * sin(phi));
300 
301  CLHEP::HepRandomEngine& fEngine; ///< art-managed random-number engine
302 
303 
304  /// Returns a vector with the name of particle selection mode keywords.
305  static std::map<int, std::string> makeParticleSelectionModeNames();
306 
307  /// Returns a vector with the name of distribution keywords.
308  static std::map<int, std::string> makeDistributionNames();
309 
310 
311  /// Act on begin of run: write "RunData" information (`sumdata::RunData`).
312  void beginRun(art::Run& run) override;
313 
314 
315  /// Performs checks and initialization based on the current configuration.
316  void setup();
317 
318  /**
319  * @brief Parses an option string and returns the corresponding option number.
320  * @tparam OptionList type of list of options (e.g. `std::map<int, std::string>`)
321  * @param Option the string of the option to be parsed
322  * @param allowedOptions list of valid options, as key/name pairs
323  * @return the key of the `Option` string from `allowedOptions`
324  * @throws std::runtime_error if `Option` is not in the option list
325  *
326  * The option string `Option` represent a single one among the supported
327  * options as defined in `allowedOptions`. The option string can be either
328  * one of the option names (the matching is not case-sensitive) or the
329  * number of the option itself.
330  *
331  * `OptionList` requirements
332  * --------------------------
333  *
334  * `OptionList` must behave like a sequence with forward iterators.
335  * Each element must behave as a pair, whose first element is the option key
336  * and the second element is the option name, equivalent to a string in that
337  * it must be forward-iterable and its elements can be converted by
338  * `std::tolower()`. The key type has no requirements beside being copiable.
339  */
340  template <typename OptionList>
341  static auto selectOption
342  (std::string Option, OptionList const& allowedOptions) -> decltype(auto);
343 
344  /**
345  * @brief Returns a string describing all options in the list
346  * @tparam OptionList type of list of options (e.g. `std::map<int, std::string>`)
347  * @param allowedOptions the list of allowed options
348  * @param printKey whether to print the key of the option beside its name
349  * @param excludeKeys list of keys to be ignored (none by default)
350  * @return a string with all options in a line
351  *
352  * The result string is a list of option names, separated by commas, like in
353  * `"'apple', 'orange', 'banana'"`. If `printKey` is `true`, the key of each
354  * option is also written in parentheses, like in
355  * `"'apple' (1), 'orange' (7), 'banana' (2)"`.
356  */
357  template <typename OptionList>
358  static std::string presentOptions(
359  OptionList const& allowedOptions, bool printKey,
360  std::initializer_list<typename OptionList::value_type::first_type> exclude
361  );
362 
363  template <typename OptionList>
364  static std::string presentOptions
365  (OptionList const& allowedOptions, bool printKey = true)
366  { return presentOptions(allowedOptions, printKey, {}); }
367 
368 
369  /// Returns the name of the specified option key, or `defName` if not known.
370  template <typename OptionList>
371  static std::string optionName(
372  typename OptionList::value_type::first_type optionKey,
373  OptionList const& allowedOptions,
374  std::string defName = "<unknown>"
375  );
376 
377  }; // class SingleGen
378 }
379 
380 namespace evgen{
381 
382  std::map<int, std::string> SingleGen::makeParticleSelectionModeNames() {
383  std::map<int, std::string> names;
384  names[int(kSelectAllParts )] = "all";
385  names[int(kSelectOneRandPart)] = "singleRandom";
386  return names;
387  } // SingleGen::makeParticleSelectionModeNames()
388 
389  std::map<int, std::string> SingleGen::makeDistributionNames() {
390  std::map<int, std::string> names;
391  names[int(kUNIF)] = "uniform";
392  names[int(kGAUS)] = "Gaussian";
393  names[int(kHIST)] = "histograms";
394  return names;
395  } // SingleGen::makeDistributionNames()
396 
397  const std::map<int, std::string> SingleGen::ParticleSelectionModeNames
399  const std::map<int, std::string> SingleGen::DistributionNames
401 
402 
403  template <typename OptionList>
405  (std::string Option, OptionList const& allowedOptions) -> decltype(auto)
406  {
407  using key_type = typename OptionList::value_type::first_type;
408  using tolower_type = int(*)(int);
409  auto toLower = [](auto const& S)
410  {
411  std::string s;
412  s.reserve(S.size());
413  std::transform(S.cbegin(), S.cend(), std::back_inserter(s),
414  (tolower_type) &std::tolower);
415  return s;
416  };
417  auto option = toLower(Option);
418  for (auto const& candidate: allowedOptions) {
419  if (toLower(candidate.second) == option) return candidate.first;
420  }
421  try {
422  std::size_t end;
423  key_type num = std::stoi(Option, &end);
424  if (allowedOptions.count(num) && (end == Option.length())) return num;
425  }
426  catch (std::invalid_argument const&) {}
427  throw std::runtime_error("Option '" + Option + "' not supported.");
428  } // SingleGen::selectOption()
429 
430 
431  template <typename OptionList>
433  OptionList const& allowedOptions, bool printKey /* = true */,
434  std::initializer_list<typename OptionList::value_type::first_type> exclude /* = {} */
435  ) {
436  std::string msg;
437 
438  unsigned int n = 0;
439  for (auto const& option: allowedOptions) {
440  auto const& key = option.first;
441  if (std::find(exclude.begin(), exclude.end(), key) != exclude.end())
442  continue;
443  if (n++ > 0) msg += ", ";
444  msg += '\"' + std::string(option.second) + '\"';
445  if (printKey)
446  msg += " (" + std::to_string(key) + ")";
447  } // for
448  return msg;
449  } // SingleGen::presentOptions()
450 
451 
452  template <typename OptionList>
454  typename OptionList::value_type::first_type optionKey,
455  OptionList const& allowedOptions,
456  std::string defName /* = "<unknown>" */
457  ) {
458  auto iOption = allowedOptions.find(optionKey);
459  return (iOption != allowedOptions.end())? iOption->second: defName;
460  } // SingleGen::optionName()
461 
462 
463  //____________________________________________________________________________
464  bool SingleGen::Config::fromHistogram(std::string const& key) const {
466  } // SingleGen::Config::fromHistogram()
467 
468  //____________________________________________________________________________
470  : EDProducer{config}
471  , fMode (selectOption(config().ParticleSelectionMode(), ParticleSelectionModeNames))
472  , fPadOutVectors(config().PadOutVectors())
473  , fPDG (config().PDG())
474  , fP0 (config().P0())
475  , fSigmaP (config().SigmaP())
476  , fPDist (selectOption(config().PDist(), DistributionNames))
477  , fX0 (config().X0())
478  , fY0 (config().Y0())
479  , fZ0 (config().Z0())
480  , fT0 (config().T0())
481  , fSigmaX (config().SigmaX())
482  , fSigmaY (config().SigmaY())
483  , fSigmaZ (config().SigmaZ())
484  , fSigmaT (config().SigmaT())
485  , fPosDist (selectOption(config().PosDist(), DistributionNames))
486  , fTDist (selectOption(config().TDist(), DistributionNames))
487  , fSingleVertex (config().SingleVertex())
488  , fTheta0XZ (config().Theta0XZ())
489  , fTheta0YZ (config().Theta0YZ())
490  , fSigmaThetaXZ (config().SigmaThetaXZ())
491  , fSigmaThetaYZ (config().SigmaThetaYZ())
492  , fAngleDist (selectOption(config().AngleDist(), DistributionNames))
493  , fHistFileName (config().HistogramFile())
494  , fPHist (config().PHist())
495  , fThetaXzYzHist(config().ThetaXzYzHist())
496  , fEngine(art::ServiceHandle<rndm::NuRandomService>()->createEngine(*this))
497  {
498  setup();
499  rndm::NuRandomService::seed_t seed;
500  if (config().Seed(seed)) {
501  fEngine.setSeed(seed, 0 /* dummy? */);
502  }
503 
504  produces< std::vector<simb::MCTruth> >();
505  produces< sumdata::RunData, art::InRun >();
506 
507  }
508 
509 
510  //____________________________________________________________________________
511  void SingleGen::beginRun(art::Run& run)
512  {
513  geo::GeometryCore const& geom = *(lar::providerFrom<geo::Geometry>());
514  run.put
515  (std::make_unique<sumdata::RunData>(geom.DetectorName()), art::fullRun());
516  }
517 
518 
519  //____________________________________________________________________________
521  {
522  // do not put seed in reconfigure because we don't want to reset
523  // the seed midstream
524  std::vector<std::string> vlist(15);
525  vlist[0] = "PDG";
526  vlist[1] = "P0";
527  vlist[2] = "SigmaP";
528  vlist[3] = "X0";
529  vlist[4] = "Y0";
530  vlist[5] = "Z0";
531  vlist[6] = "SigmaX";
532  vlist[7] = "SigmaY";
533  vlist[8] = "SigmaZ";
534  vlist[9] = "Theta0XZ";
535  vlist[10] = "Theta0YZ";
536  vlist[11] = "SigmaThetaXZ";
537  vlist[12] = "SigmaThetaYZ";
538  vlist[13] = "T0";
539  vlist[14] = "SigmaT";
540 
541 // vlist[15] = "PHist";
542 // vlist[16] = "ThetaHist";
543 // vlist[17] = "PhiHist";
544 
545  // begin tests for multiple particle error possibilities
546  std::string list;
547  if (fPDist != kHIST) {
548  if( !this->PadVector(fP0 ) ){ list.append(vlist[1].append(", \n")); }
549  if( !this->PadVector(fSigmaP ) ){ list.append(vlist[2].append(", \n")); }
550  }
551  if( !this->PadVector(fX0 ) ){ list.append(vlist[3].append(", \n")); }
552  if( !this->PadVector(fY0 ) ){ list.append(vlist[4].append(", \n")); }
553  if( !this->PadVector(fZ0 ) ){ list.append(vlist[5].append(", \n")); }
554  if( !this->PadVector(fSigmaX ) ){ list.append(vlist[6].append(", \n")); }
555  if( !this->PadVector(fSigmaY ) ){ list.append(vlist[7].append(", \n")); }
556  if( !this->PadVector(fSigmaZ ) ){ list.append(vlist[8].append(", \n")); }
557  if( !this->PadVector(fTheta0XZ ) ){ list.append(vlist[9].append(", \n")); }
558  if( !this->PadVector(fTheta0YZ ) ){ list.append(vlist[10].append(", \n")); }
559  if( !this->PadVector(fSigmaThetaXZ) ){ list.append(vlist[11].append(", \n")); }
560  if( !this->PadVector(fSigmaThetaYZ) ){ list.append(vlist[12].append(" \n")); }
561  if( !this->PadVector(fT0 ) ){ list.append(vlist[13].append(", \n")); }
562  if( !this->PadVector(fSigmaT ) ){ list.append(vlist[14].append(", \n")); }
563 
564 
565 
566  if(list.size() > 0)
567  throw cet::exception("SingleGen") << "The "<< list
568  << "\n vector(s) defined in the fhicl files has/have "
569  << "a different size than the PDG vector "
570  << "\n and it has (they have) more than one value, "
571  << "\n disallowing sensible padding "
572  << " and/or you have set fPadOutVectors to false. \n";
573 
574  if(fPDG.size() > 1 && fPadOutVectors) this->printVecs(vlist);
575 
576  // If needed, get histograms for momentum and angle distributions
577  TFile* histFile = nullptr;
578  if (!fHistFileName.empty()) {
579  if (fHistFileName[0] == '/') {
580  // We have an absolute path, use given name exactly.
581  if (cet::file_exists(fHistFileName)) {
582  histFile = new TFile(fHistFileName.c_str());
583  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
584  delete histFile;
585  histFile = nullptr;
586  throw art::Exception(art::errors::NotFound) << "Cannot open ROOT file specified in parameter HistogramFile: \"" << fHistFileName << "\"";
587  }
588  }
589  else {
590  throw art::Exception(art::errors::NotFound) << "ROOT file specified in parameter HistogramFile: \"" << fHistFileName << "\" does not exist!";
591  }
592  }
593  else {
594  // We have a relative path, search starting from current directory.
595  std::string relative_filename{"./"};
596  relative_filename += fHistFileName;
597  if (cet::file_exists(relative_filename)) {
598  histFile = new TFile(relative_filename.c_str());
599  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
600  delete histFile;
601  histFile = nullptr;
602  throw art::Exception(art::errors::NotFound) << "Cannot open ROOT file found using relative path and originally specified in parameter HistogramFile: \"" << relative_filename << '"';
603  }
604  }
605  else {
606  cet::search_path sp{"FW_SEARCH_PATH"};
607  std::string found_filename;
608  auto found = sp.find_file(fHistFileName, found_filename);
609  if (!found) {
610  throw art::Exception(art::errors::NotFound) << "Cannot find ROOT file in current directory nor on FW_SEARCH_PATH specified in parameter HistogramFile: \"" << fHistFileName << '"';
611  }
612  histFile = new TFile(found_filename.c_str());
613  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
614  delete histFile;
615  histFile = nullptr;
616  throw art::Exception(art::errors::NotFound) << "Cannot open ROOT file found on FW_SEARCH_PATH and originally specified in parameter HistogramFile: \"" << found_filename << '"';
617  }
618  }
619  }
620  }
621 
622  //
623  // deal with position distribution
624  //
625  switch (fPosDist) {
626  case kGAUS: case kUNIF: break; // supported, no further action needed
627  default:
628  throw art::Exception(art::errors::Configuration)
629  << "Position distribution of type '"
631  << "' (" << std::to_string(fPosDist) << ") is not supported.";
632  } // switch(fPosDist)
633 
634  //
635  // deal with time distribution
636  //
637  switch (fTDist) {
638  case kGAUS: case kUNIF: break; // supported, no further action needed
639  default:
640  throw art::Exception(art::errors::Configuration)
641  << "Time distribution of type '"
643  << "' (" << std::to_string(fTDist) << ") is not supported.";
644  } // switch(fTDist)
645 
646  //
647  // deal with momentum distribution
648  //
649  switch (fPDist) {
650  case kHIST:
651  if (fPHist.size() != fPDG.size()) {
652  throw art::Exception(art::errors::Configuration)
653  << fPHist.size() << " momentum histograms to describe " << fPDG.size() << " particle types...";
654  }
655  hPHist.reserve(fPHist.size());
656  for (auto const& histName: fPHist) {
657  TH1* pHist = dynamic_cast<TH1*>(histFile->Get(histName.c_str()));
658  if (!pHist) {
659  throw art::Exception(art::errors::NotFound)
660  << "Failed to read momentum histogram '" << histName << "' from '" << histFile->GetPath() << "\'";
661  }
662  pHist->SetDirectory(nullptr); // make it independent of the input file
663  hPHist.emplace_back(pHist);
664  } // for
665  break;
666  default: // supported, no further action needed
667  break;
668  } // switch(fPDist)
669 
670  switch (fAngleDist) {
671  case kHIST:
672  if (fThetaXzYzHist.size() != fPDG.size()) {
673  throw art::Exception(art::errors::Configuration)
674  << fThetaXzYzHist.size() << " direction histograms to describe " << fPDG.size() << " particle types...";
675  }
676  hThetaXzYzHist.reserve(fThetaXzYzHist.size());
677  for (auto const& histName: fThetaXzYzHist) {
678  TH2* pHist = dynamic_cast<TH2*>(histFile->Get(histName.c_str()));
679  if (!pHist) {
680  throw art::Exception(art::errors::NotFound)
681  << "Failed to read direction histogram '" << histName << "' from '" << histFile->GetPath() << "\'";
682  }
683  pHist->SetDirectory(nullptr); // make it independent of the input file
684  hThetaXzYzHist.emplace_back(pHist);
685  } // for
686  default: // supported, no further action needed
687  break;
688  } // switch(fAngleDist)
689 
690  delete histFile;
691 
692  }
693 
694  //____________________________________________________________________________
695  bool SingleGen::PadVector(std::vector<double> &vec)
696  {
697  // check if the vec has the same size as fPDG
698  if( vec.size() != fPDG.size() ){
699  // if not padding out the vectors always cause an
700  // exception to be thrown if the vector in question
701  // is not the same size as the fPDG vector
702  // the exception is thrown in the reconfigure method
703  // that calls this one
704  if (!fPadOutVectors) return false;
705  else if( fPadOutVectors){
706  // if padding of vectors is desired but the vector in
707  // question has more than one entry it isn't clear
708  // what the padded values should be so cause
709  // an exception
710  if(vec.size() != 1) return false;
711 
712  // pad it out
713  vec.resize(fPDG.size(), vec[0]);
714 
715  }// end if padding out vectors
716  }// end if the vector size is not the same as fPDG
717 
718  return true;
719  }
720 
721  //____________________________________________________________________________
722  void SingleGen::produce(art::Event& evt)
723  {
724 
725  ///unique_ptr allows ownership to be transferred to the art::Event after the put statement
726  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
727 
728  simb::MCTruth truth;
729  truth.SetOrigin(simb::kSingleParticle);
730  Sample(truth);
731 
732  MF_LOG_DEBUG("SingleGen") << truth;
733 
734  truthcol->push_back(truth);
735 
736  evt.put(std::move(truthcol));
737 
738  return;
739  }
740 
741  //____________________________________________________________________________
742  // Draw the type, momentum and position of a single particle from the
743  // FCIHL description
744  void SingleGen::SampleOne(unsigned int i, simb::MCTruth &mct){
745 
746  CLHEP::RandFlat flat(fEngine);
747  CLHEP::RandGaussQ gauss(fEngine);
748 
749  // Choose momentum
750  double p = 0.0;
751  double m = 0.0;
752  if (fPDist == kGAUS) {
753  p = gauss.fire(fP0[i], fSigmaP[i]);
754  }
755  else if (fPDist == kHIST){
756  p = SelectFromHist(*(hPHist[i]));
757  }
758  else{// if (fPDist == kUNIF) {
759  p = fP0[i] + fSigmaP[i]*(2.0*flat.fire()-1.0);
760  }
761 // else {std::cout << "do not understand the value of PDist!";}
762 
763  static TDatabasePDG pdgt;
764  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
765  if (pdgp) m = pdgp->Mass();
766 
767  // Choose position
768  TVector3 x;
769  if (fPosDist == kGAUS) {
770  x[0] = gauss.fire(fX0[i], fSigmaX[i]);;
771  x[1] = gauss.fire(fY0[i], fSigmaY[i]);
772  x[2] = gauss.fire(fZ0[i], fSigmaZ[i]);
773  }
774  else {
775  x[0] = fX0[i] + fSigmaX[i]*(2.0*flat.fire()-1.0);
776  x[1] = fY0[i] + fSigmaY[i]*(2.0*flat.fire()-1.0);
777  x[2] = fZ0[i] + fSigmaZ[i]*(2.0*flat.fire()-1.0);
778  }
779 
780  double t = 0.;
781  if(fTDist==kGAUS){
782  t = gauss.fire(fT0[i], fSigmaT[i]);
783  }
784  else{
785  t = fT0[i] + fSigmaT[i]*(2.0*flat.fire()-1.0);
786  }
787 
788  TLorentzVector pos(x[0], x[1], x[2], t);
789 
790  // Choose angles
791  double thxz = 0;
792  double thyz = 0;
793 
794  double thyzrads = 0;
795  double thyzradsplussigma = 0;
796  double thyzradsminussigma = 0;
797 
798  if (fAngleDist == kGAUS) {
799  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
800  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
801  }
802  else if (fAngleDist == kHIST){ // Select thetaxz and thetayz from histogram
803  double thetaxz = 0;
804  double thetayz = 0;
805  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
806  thxz = (180./M_PI)*thetaxz;
807  thyz = (180./M_PI)*thetayz;
808  }
809  else {
810 
811  // Choose angles flat in phase space, which is flat in theta_xz
812  // and flat in sin(theta_yz).
813 
814  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i]*(2.0*flat.fire()-1.0);
815 
816  thyzrads = std::asin(std::sin((M_PI/180.)*(fTheta0YZ[i]))); //Taking asin of sin gives value between -Pi/2 and Pi/2 regardless of user input
817  thyzradsplussigma = TMath::Min((thyzrads + ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), M_PI/2.);
818  thyzradsminussigma = TMath::Max((thyzrads - ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), -M_PI/2.);
819 
820  //uncomment line to print angular variation info
821  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
822 
823  double sinthyzmin = std::sin(thyzradsminussigma);
824  double sinthyzmax = std::sin(thyzradsplussigma);
825  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
826  thyz = (180. / M_PI) * std::asin(sinthyz);
827  }
828 
829  double thxzrad=thxz*M_PI/180.0;
830  double thyzrad=thyz*M_PI/180.0;
831 
832  TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
833  p*std::sin(thyzrad),
834  p*std::cos(thxzrad)*std::cos(thyzrad),
835  std::sqrt(p*p+m*m));
836 
837  // set track id to -i as these are all primary particles and have id <= 0
838  int trackid = -1*(i+1);
839  std::string primary("primary");
840 
841  simb::MCParticle part(trackid, fPDG[i], primary);
842  part.AddTrajectoryPoint(pos, pvec);
843 
844  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
845  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
846  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
847 
848  mct.Add(part);
849  }
850 
851  //____________________________________________________________________________
852  // Draw the type, momentum and position for all particles from the
853  // FCIHL description. Start positions will all match but momenta and angles drawn from
854  // distributions defined in the fhicls
855  void SingleGen::SampleMany(simb::MCTruth &mct){
856 
857  CLHEP::RandFlat flat(fEngine);
858  CLHEP::RandGaussQ gauss(fEngine);
859 
860  // Choose position
861  TVector3 x;
862  if (fPosDist == kGAUS) {
863  x[0] = gauss.fire(fX0[0], fSigmaX[0]);;
864  x[1] = gauss.fire(fY0[0], fSigmaY[0]);
865  x[2] = gauss.fire(fZ0[0], fSigmaZ[0]);
866  }
867  else {
868  x[0] = fX0[0] + fSigmaX[0]*(2.0*flat.fire()-1.0);
869  x[1] = fY0[0] + fSigmaY[0]*(2.0*flat.fire()-1.0);
870  x[2] = fZ0[0] + fSigmaZ[0]*(2.0*flat.fire()-1.0);
871  }
872 
873  double t = 0.;
874  if(fTDist==kGAUS){
875  t = gauss.fire(fT0[0], fSigmaT[0]);
876  }
877  else{
878  t = fT0[0] + fSigmaT[0]*(2.0*flat.fire()-1.0);
879  }
880 
881  TLorentzVector pos(x[0], x[1], x[2], t);
882 
883  // loop through particles and select momenta and angles
884  for (unsigned int i(0); i<fPDG.size(); ++i){
885  // Choose momentum
886  double p = 0.0;
887  double m = 0.0;
888  if (fPDist == kGAUS) {
889  p = gauss.fire(fP0[i], fSigmaP[i]);
890  }
891  else if (fPDist == kHIST){
892  p = SelectFromHist(*(hPHist[i]));
893  }
894  else {
895  p = fP0[i] + fSigmaP[i]*(2.0*flat.fire()-1.0);
896  }
897 
898  static TDatabasePDG pdgt;
899  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
900  if (pdgp) m = pdgp->Mass();
901 
902 
903  // Choose angles
904  double thxz = 0;
905  double thyz = 0;
906 
907  double thyzrads = 0;
908  double thyzradsplussigma = 0;
909  double thyzradsminussigma = 0;
910 
911  if (fAngleDist == kGAUS) {
912  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
913  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
914  }
915  else if (fAngleDist == kHIST){
916  double thetaxz = 0;
917  double thetayz = 0;
918  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
919  thxz = (180./M_PI)*thetaxz;
920  thyz = (180./M_PI)*thetayz;
921  }
922  else {
923 
924  // Choose angles flat in phase space, which is flat in theta_xz
925  // and flat in sin(theta_yz).
926 
927  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i]*(2.0*flat.fire()-1.0);
928 
929  thyzrads = std::asin(std::sin((M_PI/180.)*(fTheta0YZ[i]))); //Taking asin of sin gives value between -Pi/2 and Pi/2 regardless of user input
930  thyzradsplussigma = TMath::Min((thyzrads + ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), M_PI/2.);
931  thyzradsminussigma = TMath::Max((thyzrads - ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), -M_PI/2.);
932 
933  //uncomment line to print angular variation info
934  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
935 
936  double sinthyzmin = std::sin(thyzradsminussigma);
937  double sinthyzmax = std::sin(thyzradsplussigma);
938  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
939  thyz = (180. / M_PI) * std::asin(sinthyz);
940  }
941 
942  double thxzrad=thxz*M_PI/180.0;
943  double thyzrad=thyz*M_PI/180.0;
944 
945  TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
946  p*std::sin(thyzrad),
947  p*std::cos(thxzrad)*std::cos(thyzrad),
948  std::sqrt(p*p+m*m));
949 
950  // set track id to -i as these are all primary particles and have id <= 0
951  int trackid = -1*(i+1);
952  std::string primary("primary");
953 
954  simb::MCParticle part(trackid, fPDG[i], primary);
955  part.AddTrajectoryPoint(pos, pvec);
956 
957  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
958  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
959  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
960  mct.Add(part);
961  }
962  }
963 
964 
965  //____________________________________________________________________________
966  void SingleGen::Sample(simb::MCTruth &mct)
967  {
968 
969  switch (fMode) {
970  case 0: // List generation mode: every event will have one of each
971  // particle species in the fPDG array
972  if (fSingleVertex){
973  SampleMany(mct);
974  }
975  else{
976  for (unsigned int i=0; i<fPDG.size(); ++i) {
977  SampleOne(i,mct);
978  }//end loop over particles
979  }
980  break;
981  case 1: // Random selection mode: every event will exactly one particle
982  // selected randomly from the fPDG array
983  {
984  CLHEP::RandFlat flat(fEngine);
985 
986  unsigned int i=flat.fireInt(fPDG.size());
987  SampleOne(i,mct);
988  }
989  break;
990  default:
991  mf::LogWarning("UnrecognizeOption") << "SingleGen does not recognize ParticleSelectionMode "
992  << fMode;
993  break;
994  } // switch on fMode
995 
996  return;
997  }
998 
999  //____________________________________________________________________________
1000  void SingleGen::printVecs(std::vector<std::string> const& list)
1001  {
1002 
1003  mf::LogInfo("SingleGen") << " You are using vector values for SingleGen configuration.\n "
1004  << " Some of the configuration vectors may have been padded out ,"
1005  << " because they (weren't) as long as the pdg vector"
1006  << " in your configuration. \n"
1007  << " The new input particle configuration is:\n" ;
1008 
1009  std::string values;
1010  for(size_t i = 0; i <=1; ++i){// list.size(); ++i){
1011 
1012  values.append(list[i]);
1013  values.append(": [ ");
1014 
1015  for(size_t e = 0; e < fPDG.size(); ++e){
1016  std::stringstream buf;
1017  buf.width(10);
1018  if(i == 0 ) buf << fPDG[e] << ", ";
1019  buf.precision(5);
1020  if(i == 1 ) buf << fP0[e] << ", ";
1021  if(i == 2 ) buf << fSigmaP[e] << ", ";
1022  if(i == 3 ) buf << fX0[e] << ", ";
1023  if(i == 4 ) buf << fY0[e] << ", ";
1024  if(i == 5 ) buf << fZ0[e] << ", ";
1025  if(i == 6 ) buf << fSigmaX[e] << ", ";
1026  if(i == 7 ) buf << fSigmaY[e] << ", ";
1027  if(i == 8 ) buf << fSigmaZ[e] << ", ";
1028  if(i == 9 ) buf << fTheta0XZ[e] << ", ";
1029  if(i == 10) buf << fTheta0YZ[e] << ", ";
1030  if(i == 11) buf << fSigmaThetaXZ[e] << ", ";
1031  if(i == 12) buf << fSigmaThetaYZ[e] << ", ";
1032  if(i == 13) buf << fT0[e] << ", ";
1033  if(i == 14) buf << fSigmaT[e] << ", ";
1034  values.append(buf.str());
1035  }
1036 
1037  values.erase(values.find_last_of(","));
1038  values.append(" ] \n");
1039 
1040  }// end loop over vector names in list
1041 
1042  mf::LogInfo("SingleGen") << values;
1043 
1044  return;
1045  }
1046 
1047 
1048  //____________________________________________________________________________
1049  double SingleGen::SelectFromHist(const TH1& h) // select from a 1D histogram
1050  {
1051  CLHEP::RandFlat flat(fEngine);
1052 
1053  double throw_value = h.Integral() * flat.fire();
1054  double cum_value(0);
1055  for (int i(0); i < h.GetNbinsX()+1; ++i){
1056  cum_value += h.GetBinContent(i);
1057  if (throw_value < cum_value){
1058  return flat.fire()*h.GetBinWidth(i) + h.GetBinLowEdge(i);
1059  }
1060  }
1061  return throw_value; // for some reason we've gone through all bins and failed?
1062  }
1063  //____________________________________________________________________________
1064  void SingleGen::SelectFromHist(const TH2& h, double &x, double &y) // select from a 2D histogram
1065  {
1066  CLHEP::RandFlat flat(fEngine);
1067 
1068  double throw_value = h.Integral() * flat.fire();
1069  double cum_value(0);
1070  for (int i(0); i < h.GetNbinsX()+1; ++i){
1071  for (int j(0); j < h.GetNbinsY()+1; ++j){
1072  cum_value += h.GetBinContent(i, j);
1073  if (throw_value < cum_value){
1074  x = flat.fire()*h.GetXaxis()->GetBinWidth(i) + h.GetXaxis()->GetBinLowEdge(i);
1075  y = flat.fire()*h.GetYaxis()->GetBinWidth(j) + h.GetYaxis()->GetBinLowEdge(j);
1076  return;
1077  }
1078  }
1079  }
1080  return; // for some reason we've gone through all bins and failed?
1081  }
1082  //____________________________________________________________________________
1083 
1084 
1085 }//end namespace evgen
1086 
1087 namespace evgen{
1088 
1089  DEFINE_ART_MODULE(SingleGen)
1090 
1091 }//end namespace evgen
int fPDist
How to distribute momenta (gaus or uniform)
static const std::map< int, std::string > DistributionNames
Names of all distribution modes.
process_name physics producers generator physics producers generator SigmaZ
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator Theta0XZ
Definition: gen_protons.fcl:45
process_name physics producers generator PDist
int fTDist
How to distribute t (gaus, or uniform)
module to produce single or multiple specified particles in the detector
process_name stream1 can override from command line with o or output services DetectorPropertiesService services DetectorPropertiesService services DetectorPropertiesService services DetectorPropertiesService physics analyzers pmtresponse NeutronTrackingCut services LArG4Parameters gaussian physics producers generator PDG
then source grid fermiapp products dune setup_dune_fermiapp sh exit else echo No setup file found exit fi setup
see a below echo S(symbol in a section other than those above)
Utilities related to art service access.
bool PadVector(std::vector< double > &vec)
static std::string presentOptions(OptionList const &allowedOptions, bool printKey, std::initializer_list< typename OptionList::value_type::first_type > exclude)
Returns a string describing all options in the list.
process_name physics producers generator hPHist_pi physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator Theta0YZ
static constexpr Sample_t transform(Sample_t sample)
fhicl::Sequence< double > SigmaT
createEngine fPadOutVectors
static std::map< int, std::string > makeParticleSelectionModeNames()
Returns a vector with the name of particle selection mode keywords.
process_name opflash particleana ie x
createEngine fT0
static std::map< int, std::string > makeDistributionNames()
Returns a vector with the name of distribution keywords.
fhicl::Sequence< double > SigmaP
process_name opdaq physics producers generator PosDist
Definition: gen_protons.fcl:45
std::vector< std::unique_ptr< TH2 > > hThetaXzYzHist
actual TH1 for momentum distributions
std::vector< double > fSigmaT
Variation in t position (s)
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator SigmaP
Definition: gen_protons.fcl:45
createEngine fTDist
fhicl::Sequence< double > SigmaThetaXZ
process_name opdaq physics producers generator physics producers generator physics producers generator Z0
Definition: gen_protons.fcl:45
pdgs p
Definition: selectors.fcl:22
std::vector< std::unique_ptr< TH1 > > hPHist
fhicl::OptionalAtom< rndm::NuRandomService::seed_t > Seed
process_name physics producers generator hPHist_pi physics producers generator physics producers generator SigmaX
process_name opdaq physics producers generator physics producers generator Y0
Definition: gen_protons.fcl:45
std::vector< std::string > fThetaXzYzHist
name of histogram for thetaxz/thetayz distribution
fhicl::Sequence< double > Y0
void SampleOne(unsigned int i, simb::MCTruth &mct)
bool fSingleVertex
if true - all particles produced at the same location
fhicl::Sequence< double > SigmaThetaYZ
fhicl::Sequence< double > SigmaZ
fhicl::Sequence< std::string > ThetaXzYzHist
tuple m
now if test mode generate materials, CRT shell, world, gdml header else just generate CRT shell for u...
int fPosDist
How to distribute xyz (gaus, or uniform)
void printVecs(std::vector< std::string > const &list)
fhicl::Sequence< double > P0
fhicl::Sequence< double > SigmaY
process_name physics producers generator hPHist_pi physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator AngleDist
createEngine fPDG
std::vector< double > fSigmaZ
Variation in z position (cm)
fhicl::Sequence< double > X0
static constexpr int kSelectOneRandPart
fhicl::Sequence< double > Theta0YZ
while getopts h
Access the description of detector geometry.
fhicl::Sequence< double > T0
fhicl::Atom< std::string > HistogramFile
standard_singlep gaussian distribution ie ie ie gaussian TDist
Definition: multigen.fcl:18
std::vector< double > fT0
Central t position (s) in world coordinates.
fhicl::Atom< std::string > AngleDist
void produce(art::Event &evt) override
standard_singlep gaussian distribution X0
Definition: multigen.fcl:8
process_name opflash particleana ie ie y
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
std::vector< double > fP0
Central momentum (GeV/c) to generate.
std::vector< double > fTheta0YZ
Angle in YZ plane (degrees)
createEngine fSigmaT
createEngine fMode
fhicl::Atom< bool > SingleVertex
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
fEngine(art::ServiceHandle< rndm::NuRandomService >() ->createEngine(*this, pset,"Seed"))
fhicl::Sequence< double > SigmaX
unsigned int seed
BEGIN_PROLOG vertical distance to the surface Name
fhicl::Atom< std::string > ParticleSelectionMode
static constexpr int kSelectAllParts
One particle per entry is generated.
double SelectFromHist(const TH1 &h)
void Sample(simb::MCTruth &mct)
fhicl::Sequence< int > PDG
Description of geometry of one entire detector.
static auto selectOption(std::string Option, OptionList const &allowedOptions) -> decltype(auto)
Parses an option string and returns the corresponding option number.
fhicl::Sequence< double > Theta0XZ
fhicl::Sequence< double > Z0
process_name physics producers generator physics producers generator physics producers generator SigmaT
std::vector< double > fY0
Central y position (cm) in world coordinates.
std::vector< int > fPDG
PDG code of particles to generate.
static const std::vector< std::string > names
std::vector< double > fSigmaThetaYZ
Variation in angle in YZ plane.
std::vector< double > fSigmaP
Variation in momenta (GeV/c)
fhicl::Sequence< std::string > PHist
std::vector< double > fZ0
Central z position (cm) in world coordinates.
std::vector< std::string > fPHist
name of histogram of momenta
std::vector< double > fSigmaX
Variation in x position (cm)
std::string to_string(WindowPattern const &pattern)
fhicl::Atom< std::string > PosDist
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
static constexpr int kGAUS
Gaussian distribution.
bool fromHistogram(std::string const &key) const
Returns whether the specified mode is an histogram distribution.
std::vector< double > fX0
Central x position (cm) in world coordinates.
do i e
static std::string optionName(typename OptionList::value_type::first_type optionKey, OptionList const &allowedOptions, std::string defName="<unknown>")
Returns the name of the specified option key, or defName if not known.
echo Invalid option
Definition: TrainMVA.sh:17
static const std::map< int, std::string > ParticleSelectionModeNames
Names of all particle selection modes.
std::string fHistFileName
Filename containing histogram of momenta.
int fAngleDist
How to distribute angles (gaus, uniform)
SingleGen(Parameters const &config)
static constexpr int kUNIF
Uniform distribution.
std::vector< double > fTheta0XZ
Angle in XZ plane (degrees)
TCEvent evt
Definition: DataStructs.cxx:8
void beginRun(art::Run &run) override
Act on begin of run: write &quot;RunData&quot; information (sumdata::RunData).
physics producers generator PadOutVectors
process_name physics producers generator hPHist_pi physics producers generator P0
static constexpr int kHIST
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator SigmaThetaYZ
Definition: gen_protons.fcl:45
void setup()
Performs checks and initialization based on the current configuration.
fhicl::Atom< std::string > PDist
fhicl::Atom< bool > PadOutVectors
process_name physics producers generator SigmaY
art::EDProducer::Table< Config > Parameters
list
Definition: file_to_url.sh:28
std::vector< double > fSigmaThetaXZ
Variation in angle in XZ plane.
void SampleMany(simb::MCTruth &mct)
art framework interface to geometry description
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator SigmaThetaXZ
Definition: gen_protons.fcl:45
fhicl::Atom< std::string > TDist
CLHEP::HepRandomEngine & fEngine
actual TH2 for angle distributions - Xz on x axis .
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator T0
Definition: gen_protons.fcl:45
std::vector< double > fSigmaY
Variation in y position (cm)