LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
larg4SingleGen_module.cc
Go to the documentation of this file.
1 
10 // C++ includes.
11 #include <cctype> // std::tolower()
12 #include <cmath>
13 #include <initializer_list>
14 #include <iterator>
15 #include <map>
16 #include <memory>
17 #include <sstream>
18 #include <string>
19 
20 // Framework includes
25 #include "cetlib/filesystem.h"
26 #include "cetlib/search_path.h"
27 #include "cetlib_except/exception.h"
28 #include "fhiclcpp/types/Atom.h"
29 #include "fhiclcpp/types/Comment.h"
30 #include "fhiclcpp/types/Name.h"
34 
35 // nurandom includes
37 
38 // nusimdata includes
41 
42 #include "TDatabasePDG.h"
43 #include "TFile.h"
44 #include "TH1.h"
45 #include "TH2.h"
46 #include "TVector3.h"
47 
48 #include "CLHEP/Random/RandFlat.h"
49 #include "CLHEP/Random/RandGaussQ.h"
50 
51 namespace evgen {
52 
55 
56  public:
57  struct Config {
58  using Name = fhicl::Name;
60 
62  Name("ParticleSelectionMode"),
63  Comment("generate one particle, or one particle per PDG ID: " +
65 
67  Name("PadOutVectors"),
68  Comment("if true, all per-PDG-ID quantities must contain only one value, which is then "
69  "used for all PDG IDs")};
70 
72  Comment("PDG ID of the particles to be generated; this is the key "
73  "for the other options marked as \"per PDG ID\"")};
74 
76  Name("PDist"),
77  Comment("momentum distribution type: " + presentOptions(DistributionNames)),
79 
81  Comment("central momentum (GeV/c) to generate"),
82  [this]() { return !fromHistogram(PDist()); }};
83 
85  Comment("variation in momenta (GeV/c)"),
86  [this]() { return !fromHistogram(PDist()); }};
87 
89  Name("X0"),
90  Comment("central x position (cm) in world coordinates [per PDG ID]")};
91 
93  Name("Y0"),
94  Comment("central y position (cm) in world coordinates [per PDG ID]")};
95 
97  Name("Z0"),
98  Comment("central z position (cm) in world coordinates [per PDG ID]")};
99 
100  fhicl::Sequence<double> T0{Name("T0"), Comment("central time (s) [per PDG ID]")};
101 
103  Name("SigmaX"),
104  Comment("variation (radius or RMS) in x position (cm) [per PDG ID]")};
105 
107  Name("SigmaY"),
108  Comment("variation (radius or RMS) in y position (cm) [per PDG ID]")};
109 
111  Name("SigmaZ"),
112  Comment("variation (radius or RMS) in z position (cm) [per PDG ID]")};
113 
115  Name("SigmaT"),
116  Comment("variation (semi-interval or RMS) in time (s) [per PDG ID]")};
117 
119  Comment("distribution of starting position: " +
121 
123  Name("TDist"),
124  Comment("time distribution type: " + presentOptions(DistributionNames, true, {kHIST}))};
125 
127  Name("SingleVertex"),
128  Comment("if true, all particles are produced at the same location"),
129  false};
130 
132  Comment("angle from Z axis on world X-Z plane (degrees)")};
133 
135  Comment("angle from Z axis on world Y-Z plane (degrees)")};
136 
138  Comment("variation in angle in X-Z plane (degrees)")};
139 
141  Comment("variation in angle in Y-Z plane (degrees)")};
142 
144  Name("AngleDist"),
145  Comment("angular distribution type: " + presentOptions(DistributionNames)),
147 
149  Name("HistogramFile"),
150  Comment("ROOT file containing the required distributions for the generation"),
151  [this]() { return fromHistogram(AngleDist()) || fromHistogram(PDist()); }};
152 
154  Name("PHist"),
155  Comment("name of the histograms of momentum distributions"),
156  [this]() { return fromHistogram(PDist()); }};
157 
159  Name("ThetaXzYzHist"),
160  Comment("name of the histograms of angular (X-Z and Y-Z) distribution"),
161  [this]() { return fromHistogram(AngleDist()); }};
162 
164  Name("Seed"),
165  Comment("override the random number generator seed")};
166 
167  private:
169  bool fromHistogram(std::string const& key) const;
170 
171  }; // struct Config
172 
174 
175  explicit larg4SingleGen(Parameters const& config);
176 
177  // This is called for each event.
178  void produce(art::Event& evt) override;
179 
180  private:
182  static const std::map<int, std::string> ParticleSelectionModeNames;
184  static const std::map<int, std::string> DistributionNames;
185 
186  void SampleOne(unsigned int i, simb::MCTruth& mct);
187  void SampleMany(simb::MCTruth& mct);
188  void Sample(simb::MCTruth& mct);
189  void printVecs(std::vector<std::string> const& list);
190  bool PadVector(std::vector<double>& vec);
191  double SelectFromHist(const TH1& h);
192  void SelectFromHist(const TH2& h, double& x, double& y);
193 
196 
197  static constexpr int kSelectAllParts = 0;
198  static constexpr int kSelectOneRandPart =
199  1;
200 
204 
205  static constexpr int kUNIF = 0;
206  static constexpr int kGAUS = 1;
207  static constexpr int kHIST = 2;
208 
210  int fMode;
211  bool fPadOutVectors;
214  std::vector<int> fPDG;
218  std::vector<double> fP0;
219  std::vector<double> fSigmaP;
220  int fPDist;
221  std::vector<double> fX0;
222  std::vector<double> fY0;
223  std::vector<double> fZ0;
224  std::vector<double> fT0;
225  std::vector<double> fSigmaX;
226  std::vector<double> fSigmaY;
227  std::vector<double> fSigmaZ;
228  std::vector<double> fSigmaT;
229  int fPosDist;
230  int fTDist;
232  std::vector<double> fTheta0XZ;
233  std::vector<double> fTheta0YZ;
234  std::vector<double> fSigmaThetaXZ;
235  std::vector<double> fSigmaThetaYZ;
237  std::string fHistFileName;
238  std::vector<std::string> fPHist;
239  std::vector<std::string> fThetaXzYzHist;
240 
241  std::vector<std::unique_ptr<TH1>> hPHist;
242  std::vector<std::unique_ptr<TH2>>
244  // FYI - thetaxz and thetayz are related to standard polar angles as follows:
245  // thetaxz = atan2(math.sin(theta) * cos(phi), cos(theta))
246  // thetayz = asin(sin(theta) * sin(phi));
247 
248  CLHEP::HepRandomEngine& fEngine;
249 
251  static std::map<int, std::string> makeParticleSelectionModeNames();
252 
254  static std::map<int, std::string> makeDistributionNames();
255 
257  void beginRun(art::Run& run) override;
258 
260  void setup();
261 
284  template <typename OptionList>
285  static auto selectOption(std::string Option, OptionList const& allowedOptions)
286  -> decltype(auto);
287 
301  template <typename OptionList>
302  static std::string presentOptions(
303  OptionList const& allowedOptions,
304  bool printKey,
305  std::initializer_list<typename OptionList::value_type::first_type> exclude);
306 
307  template <typename OptionList>
308  static std::string presentOptions(OptionList const& allowedOptions, bool printKey = true)
309  {
310  return presentOptions(allowedOptions, printKey, {});
311  }
312 
314  template <typename OptionList>
315  static std::string optionName(typename OptionList::value_type::first_type optionKey,
316  OptionList const& allowedOptions,
317  std::string defName = "<unknown>");
318 
319  }; // class larg4SingleGen
320 }
321 
322 namespace evgen {
323 
325  {
326  std::map<int, std::string> names;
327  names[int(kSelectAllParts)] = "all";
328  names[int(kSelectOneRandPart)] = "singleRandom";
329  return names;
330  } // larg4SingleGen::makeParticleSelectionModeNames()
331 
332  std::map<int, std::string> larg4SingleGen::makeDistributionNames()
333  {
334  std::map<int, std::string> names;
335  names[int(kUNIF)] = "uniform";
336  names[int(kGAUS)] = "Gaussian";
337  names[int(kHIST)] = "histograms";
338  return names;
339  } // larg4SingleGen::makeDistributionNames()
340 
341  const std::map<int, std::string> larg4SingleGen::ParticleSelectionModeNames =
343  const std::map<int, std::string> larg4SingleGen::DistributionNames =
345 
346  template <typename OptionList>
347  auto larg4SingleGen::selectOption(std::string Option, OptionList const& allowedOptions)
348  -> decltype(auto)
349  {
350  using key_type = typename OptionList::value_type::first_type;
351  using tolower_type = int (*)(int);
352  auto toLower = [](auto const& S) {
353  std::string s;
354  s.reserve(S.size());
355  std::transform(S.cbegin(), S.cend(), std::back_inserter(s), (tolower_type)&std::tolower);
356  return s;
357  };
358  auto option = toLower(Option);
359  for (auto const& candidate : allowedOptions) {
360  if (toLower(candidate.second) == option) return candidate.first;
361  }
362  try {
363  std::size_t end;
364  key_type num = std::stoi(Option, &end);
365  if (allowedOptions.count(num) && (end == Option.length())) return num;
366  }
367  catch (std::invalid_argument const&) {
368  }
369  throw std::runtime_error("Option '" + Option + "' not supported.");
370  } // larg4SingleGen::selectOption()
371 
372  template <typename OptionList>
374  OptionList const& allowedOptions,
375  bool printKey /* = true */,
376  std::initializer_list<typename OptionList::value_type::first_type> exclude /* = {} */
377  )
378  {
379  std::string msg;
380 
381  unsigned int n = 0;
382  for (auto const& option : allowedOptions) {
383  auto const& key = option.first;
384  if (std::find(exclude.begin(), exclude.end(), key) != exclude.end()) continue;
385  if (n++ > 0) msg += ", ";
386  msg += '\"' + std::string(option.second) + '\"';
387  if (printKey) msg += " (" + std::to_string(key) + ")";
388  } // for
389  return msg;
390  } // larg4SingleGen::presentOptions()
391 
392  template <typename OptionList>
393  std::string larg4SingleGen::optionName(typename OptionList::value_type::first_type optionKey,
394  OptionList const& allowedOptions,
395  std::string defName /* = "<unknown>" */
396  )
397  {
398  auto iOption = allowedOptions.find(optionKey);
399  return (iOption != allowedOptions.end()) ? iOption->second : defName;
400  } // larg4SingleGen::optionName()
401 
402  //____________________________________________________________________________
403  bool larg4SingleGen::Config::fromHistogram(std::string const& key) const
404  {
406  } // larg4SingleGen::Config::fromHistogram()
407 
408  //____________________________________________________________________________
410  : EDProducer{config}
411  , fMode(selectOption(config().ParticleSelectionMode(), ParticleSelectionModeNames))
412  , fPadOutVectors(config().PadOutVectors())
413  , fPDG(config().PDG())
414  , fP0(config().P0())
415  , fSigmaP(config().SigmaP())
416  , fPDist(selectOption(config().PDist(), DistributionNames))
417  , fX0(config().X0())
418  , fY0(config().Y0())
419  , fZ0(config().Z0())
420  , fT0(config().T0())
421  , fSigmaX(config().SigmaX())
422  , fSigmaY(config().SigmaY())
423  , fSigmaZ(config().SigmaZ())
424  , fSigmaT(config().SigmaT())
425  , fPosDist(selectOption(config().PosDist(), DistributionNames))
426  , fTDist(selectOption(config().TDist(), DistributionNames))
427  , fSingleVertex(config().SingleVertex())
428  , fTheta0XZ(config().Theta0XZ())
429  , fTheta0YZ(config().Theta0YZ())
430  , fSigmaThetaXZ(config().SigmaThetaXZ())
431  , fSigmaThetaYZ(config().SigmaThetaYZ())
432  , fAngleDist(selectOption(config().AngleDist(), DistributionNames))
433  , fHistFileName(config().HistogramFile())
434  , fPHist(config().PHist())
435  , fThetaXzYzHist(config().ThetaXzYzHist())
436  , fEngine(createEngine(0))
437  {
438  setup();
440  if (config().Seed(seed)) { fEngine.setSeed(seed, 0 /* dummy? */); }
441 
442  produces<std::vector<simb::MCTruth>>();
443  // produces< sumdata::RunData, art::InRun >();
444  }
445 
446  //____________________________________________________________________________
448 
449  //____________________________________________________________________________
451  {
452  // do not put seed in reconfigure because we don't want to reset
453  // the seed midstream
454  std::vector<std::string> vlist(15);
455  vlist[0] = "PDG";
456  vlist[1] = "P0";
457  vlist[2] = "SigmaP";
458  vlist[3] = "X0";
459  vlist[4] = "Y0";
460  vlist[5] = "Z0";
461  vlist[6] = "SigmaX";
462  vlist[7] = "SigmaY";
463  vlist[8] = "SigmaZ";
464  vlist[9] = "Theta0XZ";
465  vlist[10] = "Theta0YZ";
466  vlist[11] = "SigmaThetaXZ";
467  vlist[12] = "SigmaThetaYZ";
468  vlist[13] = "T0";
469  vlist[14] = "SigmaT";
470 
471  // begin tests for multiple particle error possibilities
472  std::string list;
473  if (fPDist != kHIST) {
474  if (!this->PadVector(fP0)) { list.append(vlist[1].append(", \n")); }
475  if (!this->PadVector(fSigmaP)) { list.append(vlist[2].append(", \n")); }
476  }
477  if (!this->PadVector(fX0)) { list.append(vlist[3].append(", \n")); }
478  if (!this->PadVector(fY0)) { list.append(vlist[4].append(", \n")); }
479  if (!this->PadVector(fZ0)) { list.append(vlist[5].append(", \n")); }
480  if (!this->PadVector(fSigmaX)) { list.append(vlist[6].append(", \n")); }
481  if (!this->PadVector(fSigmaY)) { list.append(vlist[7].append(", \n")); }
482  if (!this->PadVector(fSigmaZ)) { list.append(vlist[8].append(", \n")); }
483  if (!this->PadVector(fTheta0XZ)) { list.append(vlist[9].append(", \n")); }
484  if (!this->PadVector(fTheta0YZ)) { list.append(vlist[10].append(", \n")); }
485  if (!this->PadVector(fSigmaThetaXZ)) { list.append(vlist[11].append(", \n")); }
486  if (!this->PadVector(fSigmaThetaYZ)) { list.append(vlist[12].append(" \n")); }
487  if (!this->PadVector(fT0)) { list.append(vlist[13].append(", \n")); }
488  if (!this->PadVector(fSigmaT)) { list.append(vlist[14].append(", \n")); }
489 
490  if (list.size() > 0)
491  throw cet::exception("larg4SingleGen")
492  << "The " << list << "\n vector(s) defined in the fhicl files has/have "
493  << "a different size than the PDG vector "
494  << "\n and it has (they have) more than one value, "
495  << "\n disallowing sensible padding "
496  << " and/or you have set fPadOutVectors to false. \n";
497 
498  if (fPDG.size() > 1 && fPadOutVectors) this->printVecs(vlist);
499 
500  // If needed, get histograms for momentum and angle distributions
501  TFile* histFile = nullptr;
502  if (!fHistFileName.empty()) {
503  if (fHistFileName[0] == '/') {
504  // We have an absolute path, use given name exactly.
505  if (cet::file_exists(fHistFileName)) {
506  histFile = new TFile(fHistFileName.c_str());
507  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
508  delete histFile;
509  histFile = nullptr;
511  << "Cannot open ROOT file specified in parameter HistogramFile: \"" << fHistFileName
512  << "\"";
513  }
514  }
515  else {
517  << "ROOT file specified in parameter HistogramFile: \"" << fHistFileName
518  << "\" does not exist!";
519  }
520  }
521  else {
522  // We have a relative path, search starting from current directory.
523  std::string relative_filename{"./"};
524  relative_filename += fHistFileName;
525  if (cet::file_exists(relative_filename)) {
526  histFile = new TFile(relative_filename.c_str());
527  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
528  delete histFile;
529  histFile = nullptr;
531  << "Cannot open ROOT file found using relative path and originally specified in "
532  "parameter HistogramFile: \""
533  << relative_filename << '"';
534  }
535  }
536  else {
537  cet::search_path sp{"FW_SEARCH_PATH"};
538  std::string found_filename;
539  auto found = sp.find_file(fHistFileName, found_filename);
540  if (!found) {
542  << "Cannot find ROOT file in current directory nor on FW_SEARCH_PATH specified in "
543  "parameter HistogramFile: \""
544  << fHistFileName << '"';
545  }
546  histFile = new TFile(found_filename.c_str());
547  if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
548  delete histFile;
549  histFile = nullptr;
551  << "Cannot open ROOT file found on FW_SEARCH_PATH and originally specified in "
552  "parameter HistogramFile: \""
553  << found_filename << '"';
554  }
555  }
556  }
557  }
558 
559  //
560  // deal with position distribution
561  //
562  switch (fPosDist) {
563  case kGAUS:
564  case kUNIF: break; // supported, no further action needed
565  default:
567  << "Position distribution of type '" << optionName(fPosDist, DistributionNames) << "' ("
568  << std::to_string(fPosDist) << ") is not supported.";
569  } // switch(fPosDist)
570 
571  //
572  // deal with time distribution
573  //
574  switch (fTDist) {
575  case kGAUS:
576  case kUNIF: break; // supported, no further action needed
577  default:
579  << "Time distribution of type '" << optionName(fTDist, DistributionNames) << "' ("
580  << std::to_string(fTDist) << ") is not supported.";
581  } // switch(fTDist)
582 
583  //
584  // deal with momentum distribution
585  //
586  switch (fPDist) {
587  case kHIST:
588  if (fPHist.size() != fPDG.size()) {
590  << fPHist.size() << " momentum histograms to describe " << fPDG.size()
591  << " particle types...";
592  }
593  hPHist.reserve(fPHist.size());
594  for (auto const& histName : fPHist) {
595  TH1* pHist = dynamic_cast<TH1*>(histFile->Get(histName.c_str()));
596  if (!pHist) {
598  << "Failed to read momentum histogram '" << histName << "' from '"
599  << histFile->GetPath() << "\'";
600  }
601  pHist->SetDirectory(nullptr); // make it independent of the input file
602  hPHist.emplace_back(pHist);
603  } // for
604  break;
605  default: // supported, no further action needed
606  break;
607  } // switch(fPDist)
608 
609  switch (fAngleDist) {
610  case kHIST:
611  if (fThetaXzYzHist.size() != fPDG.size()) {
613  << fThetaXzYzHist.size() << " direction histograms to describe " << fPDG.size()
614  << " particle types...";
615  }
616  hThetaXzYzHist.reserve(fThetaXzYzHist.size());
617  for (auto const& histName : fThetaXzYzHist) {
618  TH2* pHist = dynamic_cast<TH2*>(histFile->Get(histName.c_str()));
619  if (!pHist) {
621  << "Failed to read direction histogram '" << histName << "' from '"
622  << histFile->GetPath() << "\'";
623  }
624  pHist->SetDirectory(nullptr); // make it independent of the input file
625  hThetaXzYzHist.emplace_back(pHist);
626  } // for
627  default: // supported, no further action needed
628  break;
629  } // switch(fAngleDist)
630 
631  delete histFile;
632  }
633 
634  //____________________________________________________________________________
635  bool larg4SingleGen::PadVector(std::vector<double>& vec)
636  {
637  // check if the vec has the same size as fPDG
638  if (vec.size() != fPDG.size()) {
639  // if not padding out the vectors always cause an
640  // exception to be thrown if the vector in question
641  // is not the same size as the fPDG vector
642  // the exception is thrown in the reconfigure method
643  // that calls this one
644  if (!fPadOutVectors)
645  return false;
646  else if (fPadOutVectors) {
647  // if padding of vectors is desired but the vector in
648  // question has more than one entry it isn't clear
649  // what the padded values should be so cause
650  // an exception
651  if (vec.size() != 1) return false;
652 
653  // pad it out
654  vec.resize(fPDG.size(), vec[0]);
655 
656  } // end if padding out vectors
657  } // end if the vector size is not the same as fPDG
658 
659  return true;
660  }
661 
662  //____________________________________________________________________________
664  {
665 
667  std::unique_ptr<std::vector<simb::MCTruth>> truthcol(new std::vector<simb::MCTruth>);
668 
669  simb::MCTruth truth;
671  Sample(truth);
672 
673  MF_LOG_DEBUG("larg4SingleGen") << truth;
674 
675  truthcol->push_back(truth);
676 
677  evt.put(std::move(truthcol));
678 
679  return;
680  }
681 
682  //____________________________________________________________________________
683  // Draw the type, momentum and position of a single particle from the
684  // FCIHL description
685  void larg4SingleGen::SampleOne(unsigned int i, simb::MCTruth& mct)
686  {
687 
688  CLHEP::RandFlat flat(fEngine);
689  CLHEP::RandGaussQ gauss(fEngine);
690 
691  // Choose momentum
692  double p = 0.0;
693  double m = 0.0;
694  if (fPDist == kGAUS) { p = gauss.fire(fP0[i], fSigmaP[i]); }
695  else if (fPDist == kHIST) {
696  p = SelectFromHist(*(hPHist[i]));
697  }
698  else { // if (fPDist == kUNIF) {
699  p = fP0[i] + fSigmaP[i] * (2.0 * flat.fire() - 1.0);
700  }
701  // else {std::cout << "do not understand the value of PDist!";}
702 
703  static TDatabasePDG pdgt;
704  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
705  if (pdgp) m = pdgp->Mass();
706 
707  // Choose position
708  TVector3 x;
709  if (fPosDist == kGAUS) {
710  x[0] = gauss.fire(fX0[i], fSigmaX[i]);
711  ;
712  x[1] = gauss.fire(fY0[i], fSigmaY[i]);
713  x[2] = gauss.fire(fZ0[i], fSigmaZ[i]);
714  }
715  else {
716  x[0] = fX0[i] + fSigmaX[i] * (2.0 * flat.fire() - 1.0);
717  x[1] = fY0[i] + fSigmaY[i] * (2.0 * flat.fire() - 1.0);
718  x[2] = fZ0[i] + fSigmaZ[i] * (2.0 * flat.fire() - 1.0);
719  }
720 
721  double t = 0.;
722  if (fTDist == kGAUS) { t = gauss.fire(fT0[i], fSigmaT[i]); }
723  else {
724  t = fT0[i] + fSigmaT[i] * (2.0 * flat.fire() - 1.0);
725  }
726 
727  TLorentzVector pos(x[0], x[1], x[2], t);
728 
729  // Choose angles
730  double thxz = 0;
731  double thyz = 0;
732 
733  double thyzrads = 0;
734  double thyzradsplussigma = 0;
735  double thyzradsminussigma = 0;
736 
737  if (fAngleDist == kGAUS) {
738  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
739  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
740  }
741  else if (fAngleDist == kHIST) { // Select thetaxz and thetayz from histogram
742  double thetaxz = 0;
743  double thetayz = 0;
744  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
745  thxz = (180. / M_PI) * thetaxz;
746  thyz = (180. / M_PI) * thetayz;
747  }
748  else {
749 
750  // Choose angles flat in phase space, which is flat in theta_xz
751  // and flat in sin(theta_yz).
752 
753  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i] * (2.0 * flat.fire() - 1.0);
754 
755  thyzrads = std::asin(std::sin(
756  (M_PI / 180.) *
757  (fTheta0YZ
758  [i]))); //Taking asin of sin gives value between -Pi/2 and Pi/2 regardless of user input
759  thyzradsplussigma =
760  TMath::Min((thyzrads + ((M_PI / 180.) * fabs(fSigmaThetaYZ[i]))), M_PI / 2.);
761  thyzradsminussigma =
762  TMath::Max((thyzrads - ((M_PI / 180.) * fabs(fSigmaThetaYZ[i]))), -M_PI / 2.);
763 
764  //uncomment line to print angular variation info
765  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
766 
767  double sinthyzmin = std::sin(thyzradsminussigma);
768  double sinthyzmax = std::sin(thyzradsplussigma);
769  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
770  thyz = (180. / M_PI) * std::asin(sinthyz);
771  }
772 
773  double thxzrad = thxz * M_PI / 180.0;
774  double thyzrad = thyz * M_PI / 180.0;
775 
776  TLorentzVector pvec(p * std::cos(thyzrad) * std::sin(thxzrad),
777  p * std::sin(thyzrad),
778  p * std::cos(thxzrad) * std::cos(thyzrad),
779  std::sqrt(p * p + m * m));
780 
781  // set track id to -i as these are all primary particles and have id <= 0
782  int trackid = -1 * (i + 1);
783  std::string primary("primary");
784 
785  simb::MCParticle part(trackid, fPDG[i], primary);
786  part.AddTrajectoryPoint(pos, pvec);
787 
788  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
789  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
790  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
791 
792  mct.Add(part);
793  }
794 
795  //____________________________________________________________________________
796  // Draw the type, momentum and position for all particles from the
797  // FCIHL description. Start positions will all match but momenta and angles drawn from
798  // distributions defined in the fhicls
800  {
801 
802  CLHEP::RandFlat flat(fEngine);
803  CLHEP::RandGaussQ gauss(fEngine);
804 
805  // Choose position
806  TVector3 x;
807  if (fPosDist == kGAUS) {
808  x[0] = gauss.fire(fX0[0], fSigmaX[0]);
809  ;
810  x[1] = gauss.fire(fY0[0], fSigmaY[0]);
811  x[2] = gauss.fire(fZ0[0], fSigmaZ[0]);
812  }
813  else {
814  x[0] = fX0[0] + fSigmaX[0] * (2.0 * flat.fire() - 1.0);
815  x[1] = fY0[0] + fSigmaY[0] * (2.0 * flat.fire() - 1.0);
816  x[2] = fZ0[0] + fSigmaZ[0] * (2.0 * flat.fire() - 1.0);
817  }
818 
819  double t = 0.;
820  if (fTDist == kGAUS) { t = gauss.fire(fT0[0], fSigmaT[0]); }
821  else {
822  t = fT0[0] + fSigmaT[0] * (2.0 * flat.fire() - 1.0);
823  }
824 
825  TLorentzVector pos(x[0], x[1], x[2], t);
826 
827  // loop through particles and select momenta and angles
828  for (unsigned int i(0); i < fPDG.size(); ++i) {
829  // Choose momentum
830  double p = 0.0;
831  double m = 0.0;
832  if (fPDist == kGAUS) { p = gauss.fire(fP0[i], fSigmaP[i]); }
833  else if (fPDist == kHIST) {
834  p = SelectFromHist(*(hPHist[i]));
835  }
836  else {
837  p = fP0[i] + fSigmaP[i] * (2.0 * flat.fire() - 1.0);
838  }
839 
840  static TDatabasePDG pdgt;
841  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
842  if (pdgp) m = pdgp->Mass();
843 
844  // Choose angles
845  double thxz = 0;
846  double thyz = 0;
847 
848  double thyzrads = 0;
849  double thyzradsplussigma = 0;
850  double thyzradsminussigma = 0;
851 
852  if (fAngleDist == kGAUS) {
853  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
854  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
855  }
856  else if (fAngleDist == kHIST) {
857  double thetaxz = 0;
858  double thetayz = 0;
859  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
860  thxz = (180. / M_PI) * thetaxz;
861  thyz = (180. / M_PI) * thetayz;
862  }
863  else {
864 
865  // Choose angles flat in phase space, which is flat in theta_xz
866  // and flat in sin(theta_yz).
867 
868  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i] * (2.0 * flat.fire() - 1.0);
869 
870  thyzrads = std::asin(std::sin(
871  (M_PI / 180.) *
872  (fTheta0YZ
873  [i]))); //Taking asin of sin gives value between -Pi/2 and Pi/2 regardless of user input
874  thyzradsplussigma =
875  TMath::Min((thyzrads + ((M_PI / 180.) * fabs(fSigmaThetaYZ[i]))), M_PI / 2.);
876  thyzradsminussigma =
877  TMath::Max((thyzrads - ((M_PI / 180.) * fabs(fSigmaThetaYZ[i]))), -M_PI / 2.);
878 
879  //uncomment line to print angular variation info
880  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
881 
882  double sinthyzmin = std::sin(thyzradsminussigma);
883  double sinthyzmax = std::sin(thyzradsplussigma);
884  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
885  thyz = (180. / M_PI) * std::asin(sinthyz);
886  }
887 
888  double thxzrad = thxz * M_PI / 180.0;
889  double thyzrad = thyz * M_PI / 180.0;
890 
891  TLorentzVector pvec(p * std::cos(thyzrad) * std::sin(thxzrad),
892  p * std::sin(thyzrad),
893  p * std::cos(thxzrad) * std::cos(thyzrad),
894  std::sqrt(p * p + m * m));
895 
896  // set track id to -i as these are all primary particles and have id <= 0
897  int trackid = -1 * (i + 1);
898  std::string primary("primary");
899 
900  simb::MCParticle part(trackid, fPDG[i], primary);
901  part.AddTrajectoryPoint(pos, pvec);
902 
903  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
904  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
905  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
906  mct.Add(part);
907  }
908  }
909 
910  //____________________________________________________________________________
912  {
913 
914  switch (fMode) {
915  case 0: // List generation mode: every event will have one of each
916  // particle species in the fPDG array
917  if (fSingleVertex) { SampleMany(mct); }
918  else {
919  for (unsigned int i = 0; i < fPDG.size(); ++i) {
920  SampleOne(i, mct);
921  } //end loop over particles
922  }
923  break;
924  case 1: // Random selection mode: every event will exactly one particle
925  // selected randomly from the fPDG array
926  {
927  CLHEP::RandFlat flat(fEngine);
928 
929  unsigned int i = flat.fireInt(fPDG.size());
930  SampleOne(i, mct);
931  } break;
932  default:
933  mf::LogWarning("UnrecognizeOption")
934  << "larg4SingleGen does not recognize ParticleSelectionMode " << fMode;
935  break;
936  } // switch on fMode
937 
938  return;
939  }
940 
941  //____________________________________________________________________________
942  void larg4SingleGen::printVecs(std::vector<std::string> const& list)
943  {
944 
945  mf::LogInfo("larg4SingleGen")
946  << " You are using vector values for larg4SingleGen configuration.\n "
947  << " Some of the configuration vectors may have been padded out ,"
948  << " because they (weren't) as long as the pdg vector"
949  << " in your configuration. \n"
950  << " The new input particle configuration is:\n";
951 
952  std::string values;
953  for (size_t i = 0; i <= 1; ++i) { // list.size(); ++i){
954 
955  values.append(list[i]);
956  values.append(": [ ");
957 
958  for (size_t e = 0; e < fPDG.size(); ++e) {
959  std::stringstream buf;
960  buf.width(10);
961  if (i == 0) buf << fPDG[e] << ", ";
962  buf.precision(5);
963  if (i == 1) buf << fP0[e] << ", ";
964  if (i == 2) buf << fSigmaP[e] << ", ";
965  if (i == 3) buf << fX0[e] << ", ";
966  if (i == 4) buf << fY0[e] << ", ";
967  if (i == 5) buf << fZ0[e] << ", ";
968  if (i == 6) buf << fSigmaX[e] << ", ";
969  if (i == 7) buf << fSigmaY[e] << ", ";
970  if (i == 8) buf << fSigmaZ[e] << ", ";
971  if (i == 9) buf << fTheta0XZ[e] << ", ";
972  if (i == 10) buf << fTheta0YZ[e] << ", ";
973  if (i == 11) buf << fSigmaThetaXZ[e] << ", ";
974  if (i == 12) buf << fSigmaThetaYZ[e] << ", ";
975  if (i == 13) buf << fT0[e] << ", ";
976  if (i == 14) buf << fSigmaT[e] << ", ";
977  values.append(buf.str());
978  }
979 
980  values.erase(values.find_last_of(","));
981  values.append(" ] \n");
982 
983  } // end loop over vector names in list
984 
985  mf::LogInfo("larg4SingleGen") << values;
986 
987  return;
988  }
989 
990  //____________________________________________________________________________
991  double larg4SingleGen::SelectFromHist(const TH1& h) // select from a 1D histogram
992  {
993  CLHEP::RandFlat flat(fEngine);
994 
995  double throw_value = h.Integral() * flat.fire();
996  double cum_value(0);
997  for (int i(0); i < h.GetNbinsX() + 1; ++i) {
998  cum_value += h.GetBinContent(i);
999  if (throw_value < cum_value) { return flat.fire() * h.GetBinWidth(i) + h.GetBinLowEdge(i); }
1000  }
1001  return throw_value; // for some reason we've gone through all bins and failed?
1002  }
1003  //____________________________________________________________________________
1005  double& x,
1006  double& y) // select from a 2D histogram
1007  {
1008  CLHEP::RandFlat flat(fEngine);
1009 
1010  double throw_value = h.Integral() * flat.fire();
1011  double cum_value(0);
1012  for (int i(0); i < h.GetNbinsX() + 1; ++i) {
1013  for (int j(0); j < h.GetNbinsY() + 1; ++j) {
1014  cum_value += h.GetBinContent(i, j);
1015  if (throw_value < cum_value) {
1016  x = flat.fire() * h.GetXaxis()->GetBinWidth(i) + h.GetXaxis()->GetBinLowEdge(i);
1017  y = flat.fire() * h.GetYaxis()->GetBinWidth(j) + h.GetYaxis()->GetBinLowEdge(j);
1018  return;
1019  }
1020  }
1021  }
1022  return; // for some reason we've gone through all bins and failed?
1023  }
1024  //____________________________________________________________________________
1025 
1026 } //end namespace evgen
1027 
1028 namespace evgen {
1029 
1031 
1032 } //end namespace evgen
Float_t x
Definition: compare.C:6
static constexpr int kHIST
int fTDist
How to distribute t (gaus, or uniform)
void setup()
Performs checks and initialization based on the current configuration.
base_engine_t & createEngine(seed_t seed)
fhicl::Sequence< double > SigmaThetaYZ
std::vector< double > fSigmaZ
Variation in z position (cm)
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.
std::vector< std::unique_ptr< TH2 > > hThetaXzYzHist
actual TH1 for momentum distributions
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:258
std::vector< double > fSigmaP
Variation in momenta (GeV/c)
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
fhicl::Atom< std::string > PDist
static constexpr int kSelectAllParts
One particle per entry is generated.
fhicl::Sequence< double > Theta0YZ
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Float_t y
Definition: compare.C:6
fhicl::Sequence< double > SigmaT
std::vector< double > fSigmaThetaYZ
Variation in angle in YZ plane.
fhicl::Sequence< double > SigmaZ
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
static constexpr int kSelectOneRandPart
fhicl::Sequence< double > SigmaX
Particle class.
Definition: Run.h:37
double SelectFromHist(const TH1 &h)
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
static std::map< int, std::string > makeDistributionNames()
Returns a vector with the name of distribution keywords.
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
CLHEP::HepRandomEngine & fEngine
actual TH2 for angle distributions - Xz on x axis .
larg4SingleGen(Parameters const &config)
std::vector< int > fPDG
PDG code of particles to generate.
fhicl::Sequence< double > X0
void SampleOne(unsigned int i, simb::MCTruth &mct)
TString part[npart]
Definition: Style.C:32
static std::map< int, std::string > makeParticleSelectionModeNames()
Returns a vector with the name of particle selection mode keywords.
fhicl::Atom< std::string > TDist
std::vector< double > fSigmaT
Variation in t position (s)
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
static const std::map< int, std::string > DistributionNames
Names of all distribution modes.
long seed
Definition: chem4.cc:67
decltype(auto) constexpr to_string(T &&obj)
ADL-aware version of std::to_string.
std::vector< double > fZ0
Central z position (cm) in world coordinates.
single particles thrown at the detector
Definition: MCTruth.h:26
decltype(auto) values(Coll &&coll)
Range-for loop helper iterating across the values of the specified collection.
std::vector< double > fT0
Central t position (s) in world coordinates.
bool PadVector(std::vector< double > &vec)
bool fSingleVertex
if true - all particles produced at the same location
bool fromHistogram(std::string const &key) const
Returns whether the specified mode is an histogram distribution.
fhicl::Atom< std::string > HistogramFile
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
fhicl::Atom< std::string > PosDist
fhicl::OptionalAtom< rndm::NuRandomService::seed_t > Seed
void beginRun(art::Run &run) override
Act on begin of run: write "RunData" information (sumdata::RunData).
fhicl::Sequence< double > P0
std::vector< double > fTheta0XZ
Angle in XZ plane (degrees)
static std::string presentOptions(OptionList const &allowedOptions, bool printKey=true)
static auto selectOption(std::string Option, OptionList const &allowedOptions) -> decltype(auto)
Parses an option string and returns the corresponding option number.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
void SampleMany(simb::MCTruth &mct)
fhicl::Sequence< double > SigmaThetaXZ
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
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
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
#define MF_LOG_DEBUG(id)
module to produce single or multiple specified particles in the detector
std::vector< double > fSigmaY
Variation in y position (cm)
static constexpr int kUNIF
Uniform distribution.
Char_t n[5]
std::vector< double > fY0
Central y position (cm) in world coordinates.
TCEvent evt
Definition: DataStructs.cxx:8
fhicl::Sequence< std::string > ThetaXzYzHist
Event generator information.
Definition: MCTruth.h:32
std::vector< double > fSigmaX
Variation in x position (cm)
Float_t e
Definition: plot.C:35
void Sample(simb::MCTruth &mct)
std::vector< std::string > fThetaXzYzHist
name of histogram for thetaxz/thetayz distribution
fhicl::Sequence< double > Theta0XZ
fhicl::Sequence< double > Y0
std::vector< std::string > fPHist
name of histogram of momenta
Event Generation using GENIE, cosmics or single particles.
fhicl::Sequence< double > Z0
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.
std::vector< double > fP0
Central momentum (GeV/c) to generate.
std::vector< double > fSigmaThetaXZ
Variation in angle in XZ plane.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
int fPDist
How to distribute momenta (gaus or uniform)
art::detail::EngineCreator::seed_t seed_t
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