1 #ifndef EVGEN_SINGLEGEN 10 #define EVGEN_SINGLEGEN 21 #include <initializer_list> 41 #include "cetlib_except/exception.h" 42 #include "cetlib/filesystem.h" 43 #include "cetlib/search_path.h" 58 #include "TDatabasePDG.h" 64 #include "CLHEP/Random/RandFlat.h" 65 #include "CLHEP/Random/RandGaussQ.h" 67 namespace simb {
class MCTruth; }
81 Name(
"ParticleSelectionMode"),
82 Comment(
"generate one particle, or one particle per PDG ID: " + presentOptions(ParticleSelectionModeNames))
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")
92 Comment(
"PDG ID of the particles to be generated; this is the key for the other options marked as \"per PDG ID\"")
97 Comment(
"momentum distribution type: " + presentOptions(DistributionNames)),
98 optionName(kHIST, DistributionNames)
103 Comment(
"central momentum (GeV/c) to generate"),
104 [
this](){
return !fromHistogram(PDist()); }
109 Comment(
"variation in momenta (GeV/c)"),
110 [
this](){
return !fromHistogram(PDist()); }
115 Comment(
"central x position (cm) in world coordinates [per PDG ID]")
120 Comment(
"central y position (cm) in world coordinates [per PDG ID]")
125 Comment(
"central z position (cm) in world coordinates [per PDG ID]")
130 Comment(
"central time (s) [per PDG ID]")
135 Comment(
"variation (radius or RMS) in x position (cm) [per PDG ID]")
140 Comment(
"variation (radius or RMS) in y position (cm) [per PDG ID]")
145 Comment(
"variation (radius or RMS) in z position (cm) [per PDG ID]")
150 Comment(
"variation (semi-interval or RMS) in time (s) [per PDG ID]")
155 Comment(
"distribution of starting position: " + presentOptions(DistributionNames,
true, { kHIST }))
160 Comment(
"time distribution type: " + presentOptions(DistributionNames,
true, { kHIST }))
164 Name(
"SingleVertex"),
165 Comment(
"if true, all particles are produced at the same location"),
171 Comment(
"angle from Z axis on world X-Z plane (degrees)")
176 Comment(
"angle from Z axis on world Y-Z plane (degrees)")
180 Name(
"SigmaThetaXZ"),
181 Comment(
"variation in angle in X-Z plane (degrees)")
185 Name(
"SigmaThetaYZ"),
186 Comment(
"variation in angle in Y-Z plane (degrees)")
191 Comment(
"angular distribution type: " + presentOptions(DistributionNames)),
192 optionName(kHIST, DistributionNames)
196 Name(
"HistogramFile"),
197 Comment(
"ROOT file containing the required distributions for the generation"),
198 [
this](){
return fromHistogram(AngleDist()) || fromHistogram(PDist()); }
203 Comment(
"name of the histograms of momentum distributions"),
204 [
this](){
return fromHistogram(PDist()); }
215 Name(
"ThetaXzYzHist"),
216 Comment(
"name of the histograms of angular (X-Z and Y-Z) distribution"),
217 [
this](){
return fromHistogram(AngleDist()); }
222 Comment(
"override the random number generator seed")
229 bool fromHistogram(std::string
const& key)
const;
250 void SampleOne(
unsigned int i,
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);
262 static constexpr
int kSelectAllParts = 0;
263 static constexpr
int kSelectOneRandPart = 1;
269 static constexpr
int kUNIF = 0;
270 static constexpr
int kGAUS = 1;
271 static constexpr
int kHIST = 2;
278 std::vector<int> fPDG;
306 std::vector<std::unique_ptr<TH1>>
hPHist ;
315 static std::map<int, std::string> makeParticleSelectionModeNames();
318 static std::map<int, std::string> makeDistributionNames();
346 template <
typename OptionList>
347 static auto selectOption
348 (std::string Option, OptionList
const& allowedOptions) -> decltype(
auto);
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
369 template <
typename OptionList>
370 static std::string presentOptions
371 (OptionList
const& allowedOptions,
bool printKey =
true)
372 {
return presentOptions(allowedOptions, printKey, {}); }
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>" 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";
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";
403 const std::map<int, std::string> SingleGen::ParticleSelectionModeNames
404 = SingleGen::makeParticleSelectionModeNames();
405 const std::map<int, std::string> SingleGen::DistributionNames
406 = SingleGen::makeDistributionNames();
409 template <
typename OptionList>
410 auto SingleGen::selectOption
411 (std::string Option, OptionList
const& allowedOptions) -> decltype(
auto)
413 using key_type =
typename OptionList::value_type::first_type;
414 using tolower_type = int(*)(int);
415 auto toLower = [](
auto const& S)
419 std::transform(S.cbegin(), S.cend(), std::back_inserter(s),
420 (tolower_type) &std::tolower);
423 auto option = toLower(Option);
424 for (
auto const& candidate: allowedOptions) {
425 if (toLower(candidate.second) == option)
return candidate.first;
429 key_type num = std::stoi(Option, &end);
430 if (allowedOptions.count(num) && (end == Option.length()))
return num;
432 catch (std::invalid_argument
const&) {}
433 throw std::runtime_error(
"Option '" + Option +
"' not supported.");
437 template <
typename OptionList>
438 std::string SingleGen::presentOptions(
439 OptionList
const& allowedOptions,
bool printKey ,
440 std::initializer_list<typename OptionList::value_type::first_type> exclude
445 for (
auto const& option: allowedOptions) {
446 auto const& key = option.first;
447 if (std::find(exclude.begin(), exclude.end(), key) != exclude.end())
449 if (n++ > 0) msg +=
", ";
450 msg +=
'\"' + std::string(option.second) +
'\"';
458 template <
typename OptionList>
459 std::string SingleGen::optionName(
460 typename OptionList::value_type::first_type optionKey,
461 OptionList
const& allowedOptions,
464 auto iOption = allowedOptions.find(optionKey);
465 return (iOption != allowedOptions.end())? iOption->second: defName;
470 bool SingleGen::Config::fromHistogram(std::string
const& key)
const {
471 return selectOption(PDist(), DistributionNames) == kHIST;
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())
500 , fThetaXzYzHist(config().ThetaXzYzHist())
508 if (config().Seed(seed)) {
513 produces< std::vector<simb::MCTruth> >();
524 std::vector<std::string> vlist(15);
534 vlist[9] =
"Theta0XZ";
535 vlist[10] =
"Theta0YZ";
536 vlist[11] =
"SigmaThetaXZ";
537 vlist[12] =
"SigmaThetaYZ";
539 vlist[14] =
"SigmaT";
548 if( !this->
PadVector(
fP0 ) ){ list.append(vlist[1].append(
", \n")); }
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")); }
561 if( !this->
PadVector(
fT0 ) ){ list.append(vlist[13].append(
", \n")); }
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";
577 TFile* histFile =
nullptr;
583 if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
595 std::string relative_filename{
"./"};
597 if (cet::file_exists(relative_filename)) {
598 histFile =
new TFile(relative_filename.c_str());
599 if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
602 throw art::Exception(
art::errors::NotFound) <<
"Cannot open ROOT file found using relative path and originally specified in parameter HistogramFile: \"" << relative_filename <<
'"';
606 cet::search_path sp{
"FW_SEARCH_PATH"};
607 std::string found_filename;
612 histFile =
new TFile(found_filename.c_str());
613 if (!histFile || histFile->IsZombie() || !histFile->IsOpen()) {
616 throw art::Exception(
art::errors::NotFound) <<
"Cannot open ROOT file found on FW_SEARCH_PATH and originally specified in parameter HistogramFile: \"" << found_filename <<
'"';
629 <<
"Position distribution of type '" 641 <<
"Time distribution of type '" 653 <<
fPHist.size() <<
" momentum histograms to describe " <<
fPDG.size() <<
" particle types...";
656 for (
auto const& histName:
fPHist) {
657 TH1* pHist =
dynamic_cast<TH1*
>(histFile->Get(histName.c_str()));
660 <<
"Failed to read momentum histogram '" << histName <<
"' from '" << histFile->GetPath() <<
"\'";
662 pHist->SetDirectory(
nullptr);
663 hPHist.emplace_back(pHist);
674 <<
fThetaXzYzHist.size() <<
" direction histograms to describe " <<
fPDG.size() <<
" particle types...";
678 TH2* pHist =
dynamic_cast<TH2*
>(histFile->Get(histName.c_str()));
681 <<
"Failed to read direction histogram '" << histName <<
"' from '" << histFile->GetPath() <<
"\'";
683 pHist->SetDirectory(
nullptr);
698 if( vec.size() !=
fPDG.size() ){
710 if(vec.size() != 1)
return false;
713 vec.resize(
fPDG.size(), vec[0]);
739 std::unique_ptr< std::vector<simb::MCTruth> > truthcol(
new std::vector<simb::MCTruth>);
747 truthcol->push_back(truth);
749 evt.
put(std::move(truthcol));
761 CLHEP::HepRandomEngine &engine = rng->
getEngine();
762 CLHEP::RandFlat flat(engine);
763 CLHEP::RandGaussQ gauss(engine);
775 p =
fP0[i] +
fSigmaP[i]*(2.0*flat.fire()-1.0);
779 static TDatabasePDG pdgt;
780 TParticlePDG* pdgp = pdgt.GetParticle(
fPDG[i]);
781 if (pdgp) m = pdgp->Mass();
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);
801 t =
fT0[i] +
fSigmaT[i]*(2.0*flat.fire()-1.0);
804 TLorentzVector pos(x[0], x[1], x[2], t);
811 double thyzradsplussigma = 0;
812 double thyzradsminussigma = 0;
822 thxz = (180./M_PI)*thetaxz;
823 thyz = (180./M_PI)*thetayz;
832 thyzrads = std::asin(std::sin((M_PI/180.)*(
fTheta0YZ[i])));
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.);
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);
845 double thxzrad=thxz*M_PI/180.0;
846 double thyzrad=thyz*M_PI/180.0;
848 TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
850 p*std::cos(thxzrad)*std::cos(thyzrad),
854 int trackid = -1*(i+1);
855 std::string primary(
"primary");
875 CLHEP::HepRandomEngine &engine = rng->
getEngine();
876 CLHEP::RandFlat flat(engine);
877 CLHEP::RandGaussQ gauss(engine);
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);
897 t =
fT0[0] +
fSigmaT[0]*(2.0*flat.fire()-1.0);
900 TLorentzVector pos(x[0], x[1], x[2], t);
903 for (
unsigned int i(0); i<
fPDG.size(); ++i){
914 p =
fP0[i] +
fSigmaP[i]*(2.0*flat.fire()-1.0);
917 static TDatabasePDG pdgt;
918 TParticlePDG* pdgp = pdgt.GetParticle(
fPDG[i]);
919 if (pdgp) m = pdgp->Mass();
927 double thyzradsplussigma = 0;
928 double thyzradsminussigma = 0;
938 thxz = (180./M_PI)*thetaxz;
939 thyz = (180./M_PI)*thetayz;
948 thyzrads = std::asin(std::sin((M_PI/180.)*(
fTheta0YZ[i])));
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.);
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);
961 double thxzrad=thxz*M_PI/180.0;
962 double thyzrad=thyz*M_PI/180.0;
964 TLorentzVector pvec(p*std::cos(thyzrad)*std::sin(thxzrad),
966 p*std::cos(thxzrad)*std::cos(thyzrad),
970 int trackid = -1*(i+1);
971 std::string primary(
"primary");
995 for (
unsigned int i=0; i<
fPDG.size(); ++i) {
1004 CLHEP::HepRandomEngine &engine = rng->
getEngine();
1005 CLHEP::RandFlat flat(engine);
1007 unsigned int i=flat.fireInt(
fPDG.size());
1012 mf::LogWarning(
"UnrecognizeOption") <<
"SingleGen does not recognize ParticleSelectionMode " 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" ;
1031 for(
size_t i = 0; i <=1; ++i){
1033 values.append(list[i]);
1034 values.append(
": [ ");
1036 for(
size_t e = 0;
e <
fPDG.size(); ++
e){
1037 std::stringstream buf;
1039 if(i == 0 ) buf <<
fPDG[
e] <<
", ";
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] <<
", ";
1053 if(i == 13) buf <<
fT0[
e] <<
", ";
1054 if(i == 14) buf <<
fSigmaT[
e] <<
", ";
1055 values.append(buf.str());
1058 values.erase(values.find_last_of(
","));
1059 values.append(
" ] \n");
1073 CLHEP::HepRandomEngine &engine = rng->
getEngine();
1074 CLHEP::RandFlat flat(engine);
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);
1090 CLHEP::HepRandomEngine &engine = rng->
getEngine();
1091 CLHEP::RandFlat flat(engine);
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);
1119
int fPDist
How to distribute momenta (gaus or uniform)
static const std::map< int, std::string > DistributionNames
Names of all distribution modes.
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)
void SetOrigin(simb::Origin_t origin)
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
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::RandomNumberGenerator::seed_t seed_t
void Add(simb::MCParticle &part)
int fPosDist
How to distribute xyz (gaus, or uniform)
void printVecs(std::vector< std::string > const &list)
std::vector< double > fSigmaZ
Variation in z position (cm)
ProductID put(std::unique_ptr< PROD > &&product)
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)
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
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.
std::vector< int > fPDG
PDG code of particles to generate.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
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).
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.
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.
static constexpr int kHIST
void setup()
Performs checks and initialization based on the current configuration.
Event Generation using GENIE, cosmics or single particles.
std::vector< double > fSigmaThetaXZ
Variation in angle in XZ plane.
void SampleMany(simb::MCTruth &mct)
cet::coded_exception< error, detail::translate > exception
std::vector< double > fSigmaY
Variation in y position (cm)