LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
SingleGen_module.cc
Go to the documentation of this file.
1 #ifndef EVGEN_SINGLEGEN
10 #define EVGEN_SINGLEGEN
11 
12 // C++ includes.
13 #include <iostream>
14 #include <sstream>
15 #include <string>
16 #include <cmath>
17 #include <memory>
18 #include <iterator>
19 #include <vector>
20 #include <map>
21 #include <initializer_list>
22 #include <cctype> // std::tolower()
23 
24 
25 // Framework includes
28 #include "fhiclcpp/ParameterSet.h"
29 #include "fhiclcpp/types/Name.h"
30 #include "fhiclcpp/types/Comment.h"
31 #include "fhiclcpp/types/Atom.h"
41 #include "cetlib_except/exception.h"
42 #include "cetlib/filesystem.h"
43 #include "cetlib/search_path.h"
44 
45 // art extensions
47 
48 // nutools includes
52 
53 // lar includes
56 
57 #include "TVector3.h"
58 #include "TDatabasePDG.h"
59 #include "TMath.h"
60 #include "TFile.h"
61 #include "TH1.h"
62 #include "TH2.h"
63 
64 #include "CLHEP/Random/RandFlat.h"
65 #include "CLHEP/Random/RandGaussQ.h"
66 
67 namespace simb { class MCTruth; }
68 
69 namespace evgen {
70 
72  class SingleGen : public art::EDProducer {
73 
74  public:
75 
76  struct Config {
77  using Name = fhicl::Name;
79 
80  fhicl::Atom<std::string> ParticleSelectionMode{
81  Name("ParticleSelectionMode"),
82  Comment("generate one particle, or one particle per PDG ID: " + presentOptions(ParticleSelectionModeNames))
83  };
84 
85  fhicl::Atom<bool> PadOutVectors{
86  Name("PadOutVectors"),
87  Comment("if true, all per-PDG-ID quantities must contain only one value, which is then used for all PDG IDs")
88  };
89 
91  Name("PDG"),
92  Comment("PDG ID of the particles to be generated; this is the key for the other options marked as \"per PDG ID\"")
93  };
94 
96  Name("PDist"),
97  Comment("momentum distribution type: " + presentOptions(DistributionNames)),
98  optionName(kHIST, DistributionNames)
99  };
100 
102  Name("P0"),
103  Comment("central momentum (GeV/c) to generate"),
104  [this](){ return !fromHistogram(PDist()); }
105  };
106 
108  Name("SigmaP"),
109  Comment("variation in momenta (GeV/c)"),
110  [this](){ return !fromHistogram(PDist()); }
111  };
112 
114  Name("X0"),
115  Comment("central x position (cm) in world coordinates [per PDG ID]")
116  };
117 
119  Name("Y0"),
120  Comment("central y position (cm) in world coordinates [per PDG ID]")
121  };
122 
124  Name("Z0"),
125  Comment("central z position (cm) in world coordinates [per PDG ID]")
126  };
127 
129  Name("T0"),
130  Comment("central time (s) [per PDG ID]")
131  };
132 
134  Name("SigmaX"),
135  Comment("variation (radius or RMS) in x position (cm) [per PDG ID]")
136  };
137 
139  Name("SigmaY"),
140  Comment("variation (radius or RMS) in y position (cm) [per PDG ID]")
141  };
142 
144  Name("SigmaZ"),
145  Comment("variation (radius or RMS) in z position (cm) [per PDG ID]")
146  };
147 
149  Name("SigmaT"),
150  Comment("variation (semi-interval or RMS) in time (s) [per PDG ID]")
151  };
152 
154  Name("PosDist"),
155  Comment("distribution of starting position: " + presentOptions(DistributionNames, true, { kHIST }))
156  };
157 
159  Name("TDist"),
160  Comment("time distribution type: " + presentOptions(DistributionNames, true, { kHIST }))
161  };
162 
163  fhicl::Atom<bool> SingleVertex{
164  Name("SingleVertex"),
165  Comment("if true, all particles are produced at the same location"),
166  false
167  };
168 
170  Name("Theta0XZ"),
171  Comment("angle from Z axis on world X-Z plane (degrees)")
172  };
173 
175  Name("Theta0YZ"),
176  Comment("angle from Z axis on world Y-Z plane (degrees)")
177  };
178 
180  Name("SigmaThetaXZ"),
181  Comment("variation in angle in X-Z plane (degrees)")
182  };
183 
185  Name("SigmaThetaYZ"),
186  Comment("variation in angle in Y-Z plane (degrees)")
187  };
188 
190  Name("AngleDist"),
191  Comment("angular distribution type: " + presentOptions(DistributionNames)),
192  optionName(kHIST, DistributionNames)
193  };
194 
196  Name("HistogramFile"),
197  Comment("ROOT file containing the required distributions for the generation"),
198  [this](){ return fromHistogram(AngleDist()) || fromHistogram(PDist()); }
199  };
200 
202  Name("PHist"),
203  Comment("name of the histograms of momentum distributions"),
204  [this](){ return fromHistogram(PDist()); }
205  };
206 
207  /*
208  fhicl::Sequence<std::string> ThetaPhiHist{
209  Name("ThetaPhiHist"),
210  Comment("name of the histograms of angular (theta/phi) distribution"),
211  [this](){ return fromHistogram(AngleDist()); }
212  };
213  */
215  Name("ThetaXzYzHist"),
216  Comment("name of the histograms of angular (X-Z and Y-Z) distribution"),
217  [this](){ return fromHistogram(AngleDist()); }
218  };
219 
221  Name("Seed"),
222  Comment("override the random number generator seed")
223  };
224 
225 
226  private:
227 
229  bool fromHistogram(std::string const& key) const;
230 
231  }; // struct Config
232 
233 
235 
236 
237  explicit SingleGen(Parameters const& config);
238 
239  // This is called for each event.
240  void produce(art::Event& evt);
241  void beginRun(art::Run& run);
242 
243  private:
244 
246  static const std::map<int, std::string> ParticleSelectionModeNames;
248  static const std::map<int, std::string> DistributionNames;
249 
250  void SampleOne(unsigned int i,
251  simb::MCTruth &mct);
252  void SampleMany(simb::MCTruth &mct);
253  void Sample(simb::MCTruth &mct);
254  void printVecs(std::vector<std::string> const& list);
255  bool PadVector(std::vector<double> &vec);
256  double SelectFromHist(const TH1& h);
257  void SelectFromHist(const TH2& h, double &x, double &y);
258 
261 
262  static constexpr int kSelectAllParts = 0;
263  static constexpr int kSelectOneRandPart = 1;
264 
268 
269  static constexpr int kUNIF = 0;
270  static constexpr int kGAUS = 1;
271  static constexpr int kHIST = 2;
272 
274  int fMode;
275  bool fPadOutVectors;
278  std::vector<int> fPDG;
282  std::vector<double> fP0;
283  std::vector<double> fSigmaP;
284  int fPDist;
285  std::vector<double> fX0;
286  std::vector<double> fY0;
287  std::vector<double> fZ0;
288  std::vector<double> fT0;
289  std::vector<double> fSigmaX;
290  std::vector<double> fSigmaY;
291  std::vector<double> fSigmaZ;
292  std::vector<double> fSigmaT;
293  int fPosDist;
294  int fTDist;
296  std::vector<double> fTheta0XZ;
297  std::vector<double> fTheta0YZ;
298  std::vector<double> fSigmaThetaXZ;
299  std::vector<double> fSigmaThetaYZ;
301  std::string fHistFileName;
302  std::vector<std::string> fPHist;
303 // std::vector<std::string> fThetaPhiHist; ///< name of histogram for theta/phi distribution
304  std::vector<std::string> fThetaXzYzHist;
305 
306  std::vector<std::unique_ptr<TH1>> hPHist ;
307 // std::vector<TH2*> hThetaPhiHist ; /// actual TH1 for theta distributions - Theta on x axis
308  std::vector<std::unique_ptr<TH2>> hThetaXzYzHist ;
309  // FYI - thetaxz and thetayz are related to standard polar angles as follows:
310  // thetaxz = atan2(math.sin(theta) * cos(phi), cos(theta))
311  // thetayz = asin(sin(theta) * sin(phi));
312 
313 
315  static std::map<int, std::string> makeParticleSelectionModeNames();
316 
318  static std::map<int, std::string> makeDistributionNames();
319 
320 
322  void setup();
323 
346  template <typename OptionList>
347  static auto selectOption
348  (std::string Option, OptionList const& allowedOptions) -> decltype(auto);
349 
363  template <typename OptionList>
364  static std::string presentOptions(
365  OptionList const& allowedOptions, bool printKey,
366  std::initializer_list<typename OptionList::value_type::first_type> exclude
367  );
368 
369  template <typename OptionList>
370  static std::string presentOptions
371  (OptionList const& allowedOptions, bool printKey = true)
372  { return presentOptions(allowedOptions, printKey, {}); }
373 
374 
376  template <typename OptionList>
377  static std::string optionName(
378  typename OptionList::value_type::first_type optionKey,
379  OptionList const& allowedOptions,
380  std::string defName = "<unknown>"
381  );
382 
383  }; // class SingleGen
384 }
385 
386 namespace evgen{
387 
388  std::map<int, std::string> SingleGen::makeParticleSelectionModeNames() {
389  std::map<int, std::string> names;
390  names[int(kSelectAllParts )] = "all";
391  names[int(kSelectOneRandPart)] = "singleRandom";
392  return names;
393  } // SingleGen::makeParticleSelectionModeNames()
394 
395  std::map<int, std::string> SingleGen::makeDistributionNames() {
396  std::map<int, std::string> names;
397  names[int(kUNIF)] = "uniform";
398  names[int(kGAUS)] = "Gaussian";
399  names[int(kHIST)] = "histograms";
400  return names;
401  } // SingleGen::makeDistributionNames()
402 
403  const std::map<int, std::string> SingleGen::ParticleSelectionModeNames
404  = SingleGen::makeParticleSelectionModeNames();
405  const std::map<int, std::string> SingleGen::DistributionNames
406  = SingleGen::makeDistributionNames();
407 
408 
409  template <typename OptionList>
410  auto SingleGen::selectOption
411  (std::string Option, OptionList const& allowedOptions) -> decltype(auto)
412  {
413  using key_type = typename OptionList::value_type::first_type;
414  using tolower_type = int(*)(int);
415  auto toLower = [](auto const& S)
416  {
417  std::string s;
418  s.reserve(S.size());
419  std::transform(S.cbegin(), S.cend(), std::back_inserter(s),
420  (tolower_type) &std::tolower);
421  return s;
422  };
423  auto option = toLower(Option);
424  for (auto const& candidate: allowedOptions) {
425  if (toLower(candidate.second) == option) return candidate.first;
426  }
427  try {
428  std::size_t end;
429  key_type num = std::stoi(Option, &end);
430  if (allowedOptions.count(num) && (end == Option.length())) return num;
431  }
432  catch (std::invalid_argument const&) {}
433  throw std::runtime_error("Option '" + Option + "' not supported.");
434  } // SingleGen::selectOption()
435 
436 
437  template <typename OptionList>
438  std::string SingleGen::presentOptions(
439  OptionList const& allowedOptions, bool printKey /* = true */,
440  std::initializer_list<typename OptionList::value_type::first_type> exclude /* = {} */
441  ) {
442  std::string msg;
443 
444  unsigned int n = 0;
445  for (auto const& option: allowedOptions) {
446  auto const& key = option.first;
447  if (std::find(exclude.begin(), exclude.end(), key) != exclude.end())
448  continue;
449  if (n++ > 0) msg += ", ";
450  msg += '\"' + std::string(option.second) + '\"';
451  if (printKey)
452  msg += " (" + std::to_string(key) + ")";
453  } // for
454  return msg;
455  } // SingleGen::presentOptions()
456 
457 
458  template <typename OptionList>
459  std::string SingleGen::optionName(
460  typename OptionList::value_type::first_type optionKey,
461  OptionList const& allowedOptions,
462  std::string defName /* = "<unknown>" */
463  ) {
464  auto iOption = allowedOptions.find(optionKey);
465  return (iOption != allowedOptions.end())? iOption->second: defName;
466  } // SingleGen::optionName()
467 
468 
469  //____________________________________________________________________________
470  bool SingleGen::Config::fromHistogram(std::string const& key) const {
471  return selectOption(PDist(), DistributionNames) == kHIST;
472  } // SingleGen::Config::fromHistogram()
473 
474  //____________________________________________________________________________
475  SingleGen::SingleGen(Parameters const& config)
476  : fMode (selectOption(config().ParticleSelectionMode(), ParticleSelectionModeNames))
477  , fPadOutVectors(config().PadOutVectors())
478  , fPDG (config().PDG())
479  , fP0 (config().P0())
480  , fSigmaP (config().SigmaP())
481  , fPDist (selectOption(config().PDist(), DistributionNames))
482  , fX0 (config().X0())
483  , fY0 (config().Y0())
484  , fZ0 (config().Z0())
485  , fT0 (config().T0())
486  , fSigmaX (config().SigmaX())
487  , fSigmaY (config().SigmaY())
488  , fSigmaZ (config().SigmaZ())
489  , fSigmaT (config().SigmaT())
490  , fPosDist (selectOption(config().PosDist(), DistributionNames))
491  , fTDist (selectOption(config().TDist(), DistributionNames))
492  , fTheta0XZ (config().Theta0XZ())
493  , fTheta0YZ (config().Theta0YZ())
494  , fSigmaThetaXZ (config().SigmaThetaXZ())
495  , fSigmaThetaYZ (config().SigmaThetaYZ())
496  , fAngleDist (selectOption(config().AngleDist(), DistributionNames))
497  , fHistFileName (config().HistogramFile())
498  , fPHist (config().PHist())
499 // , fThetaPhiHist (config().ThetaPhiHist())
500  , fThetaXzYzHist(config().ThetaXzYzHist())
501  {
502  setup();
503 
504  // create a default random engine; obtain the random seed from NuRandomService,
505  // unless overridden in configuration with key "Seed"
508  if (config().Seed(seed)) {
510  (seed, 0 /* dummy? */);
511  }
512 
513  produces< std::vector<simb::MCTruth> >();
514  produces< sumdata::RunData, art::InRun >();
515 
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:
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:
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()) {
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) {
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()) {
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) {
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  //____________________________________________________________________________
723  {
724 
725  // grab the geometry object to see what geometry we are using
727  std::unique_ptr<sumdata::RunData> runcol(new sumdata::RunData(geo->DetectorName()));
728 
729  run.put(std::move(runcol));
730 
731  return;
732  }
733 
734  //____________________________________________________________________________
736  {
737 
739  std::unique_ptr< std::vector<simb::MCTruth> > truthcol(new std::vector<simb::MCTruth>);
740 
741  simb::MCTruth truth;
743  Sample(truth);
744 
745  LOG_DEBUG("SingleGen") << truth;
746 
747  truthcol->push_back(truth);
748 
749  evt.put(std::move(truthcol));
750 
751  return;
752  }
753 
754  //____________________________________________________________________________
755  // Draw the type, momentum and position of a single particle from the
756  // FCIHL description
757  void SingleGen::SampleOne(unsigned int i, simb::MCTruth &mct){
758 
759  // get the random number generator service and make some CLHEP generators
761  CLHEP::HepRandomEngine &engine = rng->getEngine();
762  CLHEP::RandFlat flat(engine);
763  CLHEP::RandGaussQ gauss(engine);
764 
765  // Choose momentum
766  double p = 0.0;
767  double m = 0.0;
768  if (fPDist == kGAUS) {
769  p = gauss.fire(fP0[i], fSigmaP[i]);
770  }
771  else if (fPDist == kHIST){
772  p = SelectFromHist(*(hPHist[i]));
773  }
774  else{// if (fPDist == kUNIF) {
775  p = fP0[i] + fSigmaP[i]*(2.0*flat.fire()-1.0);
776  }
777 // else {std::cout << "do not understand the value of PDist!";}
778 
779  static TDatabasePDG pdgt;
780  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
781  if (pdgp) m = pdgp->Mass();
782 
783  // Choose position
784  TVector3 x;
785  if (fPosDist == kGAUS) {
786  x[0] = gauss.fire(fX0[i], fSigmaX[i]);;
787  x[1] = gauss.fire(fY0[i], fSigmaY[i]);
788  x[2] = gauss.fire(fZ0[i], fSigmaZ[i]);
789  }
790  else {
791  x[0] = fX0[i] + fSigmaX[i]*(2.0*flat.fire()-1.0);
792  x[1] = fY0[i] + fSigmaY[i]*(2.0*flat.fire()-1.0);
793  x[2] = fZ0[i] + fSigmaZ[i]*(2.0*flat.fire()-1.0);
794  }
795 
796  double t = 0.;
797  if(fTDist==kGAUS){
798  t = gauss.fire(fT0[i], fSigmaT[i]);
799  }
800  else{
801  t = fT0[i] + fSigmaT[i]*(2.0*flat.fire()-1.0);
802  }
803 
804  TLorentzVector pos(x[0], x[1], x[2], t);
805 
806  // Choose angles
807  double thxz = 0;
808  double thyz = 0;
809 
810  double thyzrads = 0;
811  double thyzradsplussigma = 0;
812  double thyzradsminussigma = 0;
813 
814  if (fAngleDist == kGAUS) {
815  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
816  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
817  }
818  else if (fAngleDist == kHIST){ // Select thetaxz and thetayz from histogram
819  double thetaxz = 0;
820  double thetayz = 0;
821  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
822  thxz = (180./M_PI)*thetaxz;
823  thyz = (180./M_PI)*thetayz;
824  }
825  else {
826 
827  // Choose angles flat in phase space, which is flat in theta_xz
828  // and flat in sin(theta_yz).
829 
830  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i]*(2.0*flat.fire()-1.0);
831 
832  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
833  thyzradsplussigma = TMath::Min((thyzrads + ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), M_PI/2.);
834  thyzradsminussigma = TMath::Max((thyzrads - ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), -M_PI/2.);
835 
836  //uncomment line to print angular variation info
837  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
838 
839  double sinthyzmin = std::sin(thyzradsminussigma);
840  double sinthyzmax = std::sin(thyzradsplussigma);
841  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
842  thyz = (180. / M_PI) * std::asin(sinthyz);
843  }
844 
845  double thxzrad=thxz*M_PI/180.0;
846  double thyzrad=thyz*M_PI/180.0;
847 
848  TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
849  p*std::sin(thyzrad),
850  p*std::cos(thxzrad)*std::cos(thyzrad),
851  std::sqrt(p*p+m*m));
852 
853  // set track id to -i as these are all primary particles and have id <= 0
854  int trackid = -1*(i+1);
855  std::string primary("primary");
856 
857  simb::MCParticle part(trackid, fPDG[i], primary);
858  part.AddTrajectoryPoint(pos, pvec);
859 
860  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
861  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
862  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
863 
864  mct.Add(part);
865  }
866 
867  //____________________________________________________________________________
868  // Draw the type, momentum and position for all particles from the
869  // FCIHL description. Start positions will all match but momenta and angles drawn from
870  // distributions defined in the fhicls
872 
873  // get the random number generator service and make some CLHEP generators
875  CLHEP::HepRandomEngine &engine = rng->getEngine();
876  CLHEP::RandFlat flat(engine);
877  CLHEP::RandGaussQ gauss(engine);
878 
879  // Choose position
880  TVector3 x;
881  if (fPosDist == kGAUS) {
882  x[0] = gauss.fire(fX0[0], fSigmaX[0]);;
883  x[1] = gauss.fire(fY0[0], fSigmaY[0]);
884  x[2] = gauss.fire(fZ0[0], fSigmaZ[0]);
885  }
886  else {
887  x[0] = fX0[0] + fSigmaX[0]*(2.0*flat.fire()-1.0);
888  x[1] = fY0[0] + fSigmaY[0]*(2.0*flat.fire()-1.0);
889  x[2] = fZ0[0] + fSigmaZ[0]*(2.0*flat.fire()-1.0);
890  }
891 
892  double t = 0.;
893  if(fTDist==kGAUS){
894  t = gauss.fire(fT0[0], fSigmaT[0]);
895  }
896  else{
897  t = fT0[0] + fSigmaT[0]*(2.0*flat.fire()-1.0);
898  }
899 
900  TLorentzVector pos(x[0], x[1], x[2], t);
901 
902  // loop through particles and select momenta and angles
903  for (unsigned int i(0); i<fPDG.size(); ++i){
904  // Choose momentum
905  double p = 0.0;
906  double m = 0.0;
907  if (fPDist == kGAUS) {
908  p = gauss.fire(fP0[i], fSigmaP[i]);
909  }
910  else if (fPDist == kHIST){
911  p = SelectFromHist(*(hPHist[i]));
912  }
913  else {
914  p = fP0[i] + fSigmaP[i]*(2.0*flat.fire()-1.0);
915  }
916 
917  static TDatabasePDG pdgt;
918  TParticlePDG* pdgp = pdgt.GetParticle(fPDG[i]);
919  if (pdgp) m = pdgp->Mass();
920 
921 
922  // Choose angles
923  double thxz = 0;
924  double thyz = 0;
925 
926  double thyzrads = 0;
927  double thyzradsplussigma = 0;
928  double thyzradsminussigma = 0;
929 
930  if (fAngleDist == kGAUS) {
931  thxz = gauss.fire(fTheta0XZ[i], fSigmaThetaXZ[i]);
932  thyz = gauss.fire(fTheta0YZ[i], fSigmaThetaYZ[i]);
933  }
934  else if (fAngleDist == kHIST){
935  double thetaxz = 0;
936  double thetayz = 0;
937  SelectFromHist(*(hThetaXzYzHist[i]), thetaxz, thetayz);
938  thxz = (180./M_PI)*thetaxz;
939  thyz = (180./M_PI)*thetayz;
940  }
941  else {
942 
943  // Choose angles flat in phase space, which is flat in theta_xz
944  // and flat in sin(theta_yz).
945 
946  thxz = fTheta0XZ[i] + fSigmaThetaXZ[i]*(2.0*flat.fire()-1.0);
947 
948  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
949  thyzradsplussigma = TMath::Min((thyzrads + ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), M_PI/2.);
950  thyzradsminussigma = TMath::Max((thyzrads - ((M_PI/180.)*fabs(fSigmaThetaYZ[i]))), -M_PI/2.);
951 
952  //uncomment line to print angular variation info
953  //std::cout << "Central angle: " << (180./M_PI)*thyzrads << " Max angle: " << (180./M_PI)*thyzradsplussigma << " Min angle: " << (180./M_PI)*thyzradsminussigma << std::endl;
954 
955  double sinthyzmin = std::sin(thyzradsminussigma);
956  double sinthyzmax = std::sin(thyzradsplussigma);
957  double sinthyz = sinthyzmin + flat.fire() * (sinthyzmax - sinthyzmin);
958  thyz = (180. / M_PI) * std::asin(sinthyz);
959  }
960 
961  double thxzrad=thxz*M_PI/180.0;
962  double thyzrad=thyz*M_PI/180.0;
963 
964  TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
965  p*std::sin(thyzrad),
966  p*std::cos(thxzrad)*std::cos(thyzrad),
967  std::sqrt(p*p+m*m));
968 
969  // set track id to -i as these are all primary particles and have id <= 0
970  int trackid = -1*(i+1);
971  std::string primary("primary");
972 
973  simb::MCParticle part(trackid, fPDG[i], primary);
974  part.AddTrajectoryPoint(pos, pvec);
975 
976  //std::cout << "Px: " << pvec.Px() << " Py: " << pvec.Py() << " Pz: " << pvec.Pz() << std::endl;
977  //std::cout << "x: " << pos.X() << " y: " << pos.Y() << " z: " << pos.Z() << " time: " << pos.T() << std::endl;
978  //std::cout << "YZ Angle: " << (thyzrad * (180./M_PI)) << " XZ Angle: " << (thxzrad * (180./M_PI)) << std::endl;
979  mct.Add(part);
980  }
981  }
982 
983 
984  //____________________________________________________________________________
986  {
987 
988  switch (fMode) {
989  case 0: // List generation mode: every event will have one of each
990  // particle species in the fPDG array
991  if (fSingleVertex){
992  SampleMany(mct);
993  }
994  else{
995  for (unsigned int i=0; i<fPDG.size(); ++i) {
996  SampleOne(i,mct);
997  }//end loop over particles
998  }
999  break;
1000  case 1: // Random selection mode: every event will exactly one particle
1001  // selected randomly from the fPDG array
1002  {
1004  CLHEP::HepRandomEngine &engine = rng->getEngine();
1005  CLHEP::RandFlat flat(engine);
1006 
1007  unsigned int i=flat.fireInt(fPDG.size());
1008  SampleOne(i,mct);
1009  }
1010  break;
1011  default:
1012  mf::LogWarning("UnrecognizeOption") << "SingleGen does not recognize ParticleSelectionMode "
1013  << fMode;
1014  break;
1015  } // switch on fMode
1016 
1017  return;
1018  }
1019 
1020  //____________________________________________________________________________
1021  void SingleGen::printVecs(std::vector<std::string> const& list)
1022  {
1023 
1024  mf::LogInfo("SingleGen") << " You are using vector values for SingleGen configuration.\n "
1025  << " Some of the configuration vectors may have been padded out ,"
1026  << " because they (weren't) as long as the pdg vector"
1027  << " in your configuration. \n"
1028  << " The new input particle configuration is:\n" ;
1029 
1030  std::string values;
1031  for(size_t i = 0; i <=1; ++i){// list.size(); ++i){
1032 
1033  values.append(list[i]);
1034  values.append(": [ ");
1035 
1036  for(size_t e = 0; e < fPDG.size(); ++e){
1037  std::stringstream buf;
1038  buf.width(10);
1039  if(i == 0 ) buf << fPDG[e] << ", ";
1040  buf.precision(5);
1041  if(i == 1 ) buf << fP0[e] << ", ";
1042  if(i == 2 ) buf << fSigmaP[e] << ", ";
1043  if(i == 3 ) buf << fX0[e] << ", ";
1044  if(i == 4 ) buf << fY0[e] << ", ";
1045  if(i == 5 ) buf << fZ0[e] << ", ";
1046  if(i == 6 ) buf << fSigmaX[e] << ", ";
1047  if(i == 7 ) buf << fSigmaY[e] << ", ";
1048  if(i == 8 ) buf << fSigmaZ[e] << ", ";
1049  if(i == 9 ) buf << fTheta0XZ[e] << ", ";
1050  if(i == 10) buf << fTheta0YZ[e] << ", ";
1051  if(i == 11) buf << fSigmaThetaXZ[e] << ", ";
1052  if(i == 12) buf << fSigmaThetaYZ[e] << ", ";
1053  if(i == 13) buf << fT0[e] << ", ";
1054  if(i == 14) buf << fSigmaT[e] << ", ";
1055  values.append(buf.str());
1056  }
1057 
1058  values.erase(values.find_last_of(","));
1059  values.append(" ] \n");
1060 
1061  }// end loop over vector names in list
1062 
1063  mf::LogInfo("SingleGen") << values;
1064 
1065  return;
1066  }
1067 
1068 
1069  //____________________________________________________________________________
1070  double SingleGen::SelectFromHist(const TH1& h) // select from a 1D histogram
1071  {
1073  CLHEP::HepRandomEngine &engine = rng->getEngine();
1074  CLHEP::RandFlat flat(engine);
1075 
1076  double throw_value = h.Integral() * flat.fire();
1077  double cum_value(0);
1078  for (int i(0); i < h.GetNbinsX()+1; ++i){
1079  cum_value += h.GetBinContent(i);
1080  if (throw_value < cum_value){
1081  return flat.fire()*h.GetBinWidth(i) + h.GetBinLowEdge(i);
1082  }
1083  }
1084  return throw_value; // for some reason we've gone through all bins and failed?
1085  }
1086  //____________________________________________________________________________
1087  void SingleGen::SelectFromHist(const TH2& h, double &x, double &y) // select from a 2D histogram
1088  {
1090  CLHEP::HepRandomEngine &engine = rng->getEngine();
1091  CLHEP::RandFlat flat(engine);
1092 
1093  double throw_value = h.Integral() * flat.fire();
1094  double cum_value(0);
1095  for (int i(0); i < h.GetNbinsX()+1; ++i){
1096  for (int j(0); j < h.GetNbinsY()+1; ++j){
1097  cum_value += h.GetBinContent(i, j);
1098  if (throw_value < cum_value){
1099  x = flat.fire()*h.GetXaxis()->GetBinWidth(i) + h.GetXaxis()->GetBinLowEdge(i);
1100  y = flat.fire()*h.GetYaxis()->GetBinWidth(j) + h.GetYaxis()->GetBinLowEdge(j);
1101  return;
1102  }
1103  }
1104  }
1105  return; // for some reason we've gone through all bins and failed?
1106  }
1107  //____________________________________________________________________________
1108 
1109 
1110 }//end namespace evgen
1111 
1112 namespace evgen{
1113 
1115 
1116 }//end namespace evgen
1117 
1118 #endif
1119 
Float_t x
Definition: compare.C:6
int fPDist
How to distribute momenta (gaus or uniform)
static const std::map< int, std::string > DistributionNames
Names of all distribution modes.
Float_t s
Definition: plot.C:23
int fTDist
How to distribute t (gaus, or uniform)
module to produce single or multiple specified particles in the detector
bool PadVector(std::vector< double > &vec)
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:261
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:78
std::vector< std::unique_ptr< TH2 > > hThetaXzYzHist
actual TH1 for momentum distributions
std::vector< double > fSigmaT
Variation in t position (s)
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Float_t y
Definition: compare.C:6
std::vector< std::unique_ptr< TH1 > > hPHist
std::vector< std::string > fThetaXzYzHist
name of histogram for thetaxz/thetayz distribution
void SampleOne(unsigned int i, simb::MCTruth &mct)
bool fSingleVertex
if true - all particles produced at the same location
art::ProductID put(std::unique_ptr< PROD > &&)
Definition: Run.h:148
art::RandomNumberGenerator::seed_t seed_t
Particle class.
void Add(simb::MCParticle &part)
Definition: MCTruth.h:77
int fPosDist
How to distribute xyz (gaus, or uniform)
void printVecs(std::vector< std::string > const &list)
Definition: Run.h:30
std::vector< double > fSigmaZ
Variation in z position (cm)
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
std::vector< double > fT0
Central t position (s) in world coordinates.
base_engine_t & getEngine() const
base_engine_t & createEngine(seed_t seed)
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
long seed
Definition: chem4.cc:68
std::vector< double > fP0
Central momentum (GeV/c) to generate.
std::vector< double > fTheta0YZ
Angle in YZ plane (degrees)
single particles thrown at the detector
Definition: MCTruth.h:24
TString part[npart]
Definition: Style.C:32
double SelectFromHist(const TH1 &h)
void Sample(simb::MCTruth &mct)
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
std::vector< double > fY0
Central y position (cm) in world coordinates.
Framework includes.
std::vector< int > fPDG
PDG code of particles to generate.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::vector< double > fSigmaThetaYZ
Variation in angle in YZ plane.
std::vector< double > fSigmaP
Variation in momenta (GeV/c)
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
void produce(art::Event &evt)
std::vector< double > fZ0
Central z position (cm) in world coordinates.
std::string to_string(Flag_t< Storage > const flag)
Convert a flag into a stream (shows its index).
Definition: BitMask.h:187
std::vector< std::string > fPHist
name of histogram of momenta
std::vector< double > fSigmaX
Variation in x position (cm)
static constexpr int kGAUS
Gaussian distribution.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
void beginRun(art::Run &run)
std::vector< double > fX0
Central x position (cm) in world coordinates.
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.
#define LOG_DEBUG(id)
Char_t n[5]
static const std::map< int, std::string > ParticleSelectionModeNames
Names of all particle selection modes.
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
std::string fHistFileName
Filename containing histogram of momenta.
int fAngleDist
How to distribute angles (gaus, uniform)
static constexpr int kUNIF
Uniform distribution.
std::vector< double > fTheta0XZ
Angle in XZ plane (degrees)
Event generator information.
Definition: MCTruth.h:30
Float_t e
Definition: plot.C:34
static constexpr int kHIST
void setup()
Performs checks and initialization based on the current configuration.
Namespace collecting geometry-related classes utilities.
Event Generation using GENIE, cosmics or single particles.
std::vector< double > fSigmaThetaXZ
Variation in angle in XZ plane.
void SampleMany(simb::MCTruth &mct)
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
std::vector< double > fSigmaY
Variation in y position (cm)