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