GENIEGenerator
Loading...
Searching...
No Matches
gBeamHNLEvGen.cxx File Reference
Include dependency graph for gBeamHNLEvGen.cxx:

Go to the source code of this file.

Functions

void GetCommandLineArgs (int argc, char **argv)
void PrintSyntax (void)
int SelectDecayMode (std::vector< HNLDecayMode_t > *intChannels, SimpleHNL sh)
const EventRecordVisitorIHNLGenerator (void)
TLorentzVector GeneratePosition (GHepRecord *event)
int main (int argc, char **argv)

Variables

string kDefOptGeomLUnits = "mm"
string kDefOptGeomDUnits = "g_cm3"
NtpMCFormat_t kDefOptNtpFormat = kNFGHEP
string kDefOptEvFilePrefix = "gntp"
string kDefOptFluxFilePath = "./input-flux.root"
string kDefOptSName = "genie::EventGenerator"
string kDefOptSConfig = "BeamHNL"
string kDefOptSTune = "GHNL20_00a_00_000"
Long_t gOptRunNu = 1000
int gOptNev = 10
double gOptEnergyHNL = -1
double gOptMassHNL = -1
double gOptECoupling = -1
double gOptMCoupling = -1
double gOptTCoupling = -1
bool gOptIsMajorana = false
bool gOptIsMonoEnFlux = true
HNLDecayMode_t gOptDecayMode = kHNLDcyNull
std::vector< HNLDecayMode_tgOptIntChannels
string gOptEvFilePrefix = kDefOptEvFilePrefix
bool gOptUsingRootGeom = false
string gOptRootGeom
string gOptRootGeomTopVol = ""
double gOptGeomLUnits = 0
long int gOptRanSeed = -1
double fdx = 0
double fdy = 0
double fdz = 0
double fox = 0
double foy = 0
double foz = 0
double NTP_IS_E = 0.
double NTP_IS_PX = 0.
double NTP_IS_PY = 0.
double NTP_IS_PZ = 0.
double NTP_FS0_E = 0.
double NTP_FS0_PX = 0.
double NTP_FS0_PY = 0.
double NTP_FS0_PZ = 0.
double NTP_FS1_E = 0.
double NTP_FS1_PX = 0.
double NTP_FS1_PY = 0.
double NTP_FS1_PZ = 0.
double NTP_FS2_E = 0.
double NTP_FS2_PX = 0.
double NTP_FS2_PY = 0.
double NTP_FS2_PZ = 0.
int NTP_FS0_PDG = 0
int NTP_FS1_PDG = 0
int NTP_FS2_PDG = 0
double CoMLifetime = -1.0
double decayMod = 1.0
double evWeight = 1.0

Function Documentation

◆ GeneratePosition()

TLorentzVector GeneratePosition ( GHepRecord * event)

Definition at line 711 of file gBeamHNLEvGen.cxx.

712{
713 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
714 TRandom3 & rnd_geo = rnd->RndGeom();
715
716 double rndx = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
717 double rndy = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
718 double rndz = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
719
720 double t_gen = 0;
721 double x_gen = fox + rndx * fdx;
722 double y_gen = foy + rndy * fdy;
723 double z_gen = foz + rndz * fdz;
724
725 TLorentzVector pos(x_gen, y_gen, z_gen, t_gen);
726 return pos;
727}
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
TRandom3 & RndGeom(void) const
rnd number generator used by geometry drivers
Definition RandomGen.h:68
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
double foy
double fdz
double foz
double fox
double fdy
double fdx

References fdx, fdy, fdz, fox, foy, foz, genie::RandomGen::Instance(), and genie::RandomGen::RndGeom().

Referenced by main().

◆ GetCommandLineArgs()

void GetCommandLineArgs ( int argc,
char ** argv )

Definition at line 818 of file gBeamHNLEvGen.cxx.

819{
820 LOG("gevgen_hnl", pINFO) << "Parsing command line arguments";
821
822 // Parse run options for this app
823
824 CmdLnArgParser parser(argc,argv);
825
826 // force the app to look at HNL tune by default
827 // if user passes --tune argument, look at the user input tune instead
828 char * expargv[ argc + 2 ];
829 bool didFindUserInputTune = false;
830 std::string stExtraTuneBit = kDefOptSTune;
831
832 if( parser.OptionExists("tune") ){
833 didFindUserInputTune = true;
834 stExtraTuneBit = parser.ArgAsString("tune");
835 LOG( "gevgen_hnl", pWARN )
836 << "Using input HNL tune " << parser.ArgAsString("tune");
837 } else {
838 LOG( "gevgen_hnl", pWARN )
839 << "Using default HNL tune " << kDefOptSTune;
840 }
841 // append this to argv
842 for( int iArg = 0; iArg < argc; iArg++ ){
843 expargv[iArg] = argv[iArg];
844 }
845 if( !didFindUserInputTune ){
846 char * chBit = const_cast< char * >( stExtraTuneBit.c_str() ); // ugh. Ugly.
847 std::string stune("--tune"); char * tBit = const_cast< char * >( stune.c_str() );
848 expargv[argc] = tBit;
849 expargv[argc+1] = chBit;
850 }
851
852 // Common run options.
853 int expargc = ( didFindUserInputTune ) ? argc : argc+2;
854 std::string stnull(""); char * nBit = const_cast< char * >( stnull.c_str() );
855 expargv[expargc] = nBit;
856
857 RunOpt::Instance()->ReadFromCommandLine(expargc,expargv);
858
859 // help?
860 bool help = parser.OptionExists('h');
861 if(help) {
862 PrintSyntax();
863 exit(0);
864 }
865
866 // run number
867 if( parser.OptionExists('r') ) {
868 LOG("gevgen_hnl", pDEBUG) << "Reading MC run number";
869 gOptRunNu = parser.ArgAsLong('r');
870 } else {
871 LOG("gevgen_hnl", pDEBUG) << "Unspecified run number - Using default";
872 gOptRunNu = 1000;
873 } //-r
874
875 // number of events
876 if( parser.OptionExists('n') ) {
877 LOG("gevgen_hnl", pDEBUG)
878 << "Reading number of events to generate";
879 gOptNev = parser.ArgAsInt('n');
880 } else {
881 LOG("gevgen_hnl", pFATAL)
882 << "You need to specify the number of events";
883 PrintSyntax();
884 exit(0);
885 } //-n
886
887 // get HNL mass
888 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
889 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
890 gOptMassHNL = hnlgen->GetHNLMass();
891
892 bool isMonoEnergeticFlux = true;
893#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
894 if( parser.OptionExists('f') ) {
895 LOG("gevgen_hnl", pDEBUG)
896 << "A flux has been offered. Searching this path: " << parser.ArgAsString('f');
897 isMonoEnergeticFlux = false;
898 gOptFluxFilePath = parser.ArgAsString('f');
899
900 // check if this is valid path (assume these are dk2nu files)
901 if( gSystem->OpenDirectory( gOptFluxFilePath.c_str() ) != NULL ){
902 gOptIsUsingDk2nu = true;
903 LOG("gevgen_hnl", pDEBUG)
904 << "dk2nu flux files detected. Will create flux spectrum dynamically.";
905 } else {
906 LOG("gevgen_hnl", pFATAL)
907 << "Invalid flux file path " << gOptFluxFilePath;
908 exit(1);
909 }
910 } else {
911 // we need the 'E' option! Log it and pass below
912 LOG("gevgen_hnl", pINFO)
913 << "No flux file offered. Assuming monoenergetic flux.";
914 } //-f
915#endif // #ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
916
917 // HNL energy (only relevant if we do not have an input flux)
918 gOptEnergyHNL = -1;
919 if( isMonoEnergeticFlux ){
920 if( parser.OptionExists('E') ) {
921 LOG("gevgen_hnl", pDEBUG)
922 << "Reading HNL energy";
923 gOptEnergyHNL = parser.ArgAsDouble('E');
924 } else {
925 LOG("gevgen_hnl", pFATAL)
926 << "You need to specify the HNL energy";
927 PrintSyntax();
928 exit(0);
929 } //-E
930 assert(gOptEnergyHNL > gOptMassHNL);
931 }
932
933 gOptIsMonoEnFlux = isMonoEnergeticFlux;
934
935 // first flux entry to read
936 if( parser.OptionExists("firstEvent") ) {
937 gOptFirstEvent = parser.ArgAsInt("firstEvent");
938 LOG( "gevgen_hnl", pINFO )
939 << "Starting flux readin at first event = " << gOptFirstEvent;
940 } // --firstEvent
941
942 // HNL decay mode
943 int mode = -1;
944 if( parser.OptionExists('m') ) {
945 LOG("gevgen_hnl", pDEBUG)
946 << "Reading HNL decay mode";
947 mode = parser.ArgAsInt('m');
948 } else {
949 LOG("gevgen_hnl", pINFO)
950 << "No decay mode specified - will sample from allowed decay modes";
951 } //-m
953
955 if(!allowed) {
956 LOG("gevgen_hnl", pFATAL)
957 << "Specified decay is not allowed kinematically for the given HNL mass";
958 PrintSyntax();
959 exit(0);
960 }
961
962 //
963 // geometry
964 //
965
966 string geom = "";
967 string lunits;
968#ifdef __CAN_USE_ROOT_GEOM__
969 // string dunits;
970 if( parser.OptionExists('g') ) {
971 LOG("gevgen_hnl", pDEBUG) << "Getting input geometry";
972 geom = parser.ArgAsString('g');
973
974 // is it a ROOT file that contains a ROOT geometry?
975 bool accessible_geom_file =
976 ! (gSystem->AccessPathName(geom.c_str()));
977 if (accessible_geom_file) {
979 gOptUsingRootGeom = true;
980 } else {
981 LOG("gevgen_hnl", pFATAL)
982 << "Geometry option is not a ROOT file. Please use ROOT geom.";
983 PrintSyntax();
984 exit(1);
985 }
986 } else {
987 // LOG("gevgen_hnl", pFATAL)
988 // << "No geometry option specified - Exiting";
989 // PrintSyntax();
990 // exit(1);
991 } //-g
992
994 // using a ROOT geometry - get requested geometry units
995
996 // length units:
997 if( parser.OptionExists('L') ) {
998 LOG("gevgen_hnl", pDEBUG)
999 << "Checking for input geometry length units";
1000 lunits = parser.ArgAsString('L');
1001 } else {
1002 LOG("gevgen_hnl", pDEBUG) << "Using default geometry length units";
1004 } // -L
1005 // // density units:
1006 // if( parser.OptionExists('D') ) {
1007 // LOG("gevgen_hnl", pDEBUG)
1008 // << "Checking for input geometry density units";
1009 // dunits = parser.ArgAsString('D');
1010 // } else {
1011 // LOG("gevgen_hnl", pDEBUG) << "Using default geometry density units";
1012 // dunits = kDefOptGeomDUnits;
1013 // } // -D
1015 // gOptGeomDUnits = utils::units::UnitFromString(dunits);
1016
1017 } // using root geom?
1018#endif // #ifdef __CAN_USE_ROOT_GEOM__
1019
1020 // event file prefix
1021 if( parser.OptionExists('o') ) {
1022 LOG("gevgen_hnl", pDEBUG) << "Reading the event filename prefix";
1023 gOptEvFilePrefix = parser.ArgAsString('o');
1024 } else {
1025 LOG("gevgen_hnl", pDEBUG)
1026 << "Will set the default event filename prefix";
1028 } //-o
1029
1030 // random number seed
1031 if( parser.OptionExists("seed") ) {
1032 LOG("gevgen_hnl", pINFO) << "Reading random number seed";
1033 gOptRanSeed = parser.ArgAsLong("seed");
1034 } else {
1035 LOG("gevgen_hnl", pINFO) << "Unspecified random number seed - Using default";
1036 gOptRanSeed = -1;
1037 }
1038
1039 //
1040 // >>> print the command line options
1041 //
1042
1043 ostringstream gminfo;
1044 if (gOptUsingRootGeom) {
1045 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
1046 << ", top volume: "
1047 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
1048 << ", length units: " << lunits;
1049 // << ", density units: " << dunits;
1050 }
1051
1052 LOG("gevgen_hnl", pNOTICE)
1053 << "\n\n"
1054 << utils::print::PrintFramedMesg("gevgen_hnl job configuration");
1055
1056 LOG("gevgen_hnl", pNOTICE)
1057 << "\n @@ Run number : " << gOptRunNu
1058 << "\n @@ Random seed : " << gOptRanSeed
1059 << "\n @@ HNL mass : " << gOptMassHNL << " GeV"
1060 << "\n @@ Decay channel : " << utils::hnl::AsString(gOptDecayMode)
1061 << "\n @@ Geometry : " << gminfo.str()
1062 << "\n @@ Statistics : " << gOptNev << " events";
1063}
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pFATAL
Definition Messenger.h:56
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define pWARN
Definition Messenger.h:60
const Algorithm * GetAlgorithm(const AlgId &algid)
static AlgFactory * Instance()
Algorithm abstract base class.
Definition Algorithm.h:54
Command line argument parser.
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
Heavy Neutral Lepton final-state product generator.
Definition HNLDecayer.h:41
double GetHNLMass() const
long int gOptRanSeed
double gOptGeomLUnits
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
bool gOptUsingRootGeom
Long_t gOptRunNu
int gOptNev
string gOptRootGeom
string gOptEvFilePrefix
string gOptRootGeomTopVol
double gOptMassHNL
HNLDecayMode_t gOptDecayMode
bool gOptIsMonoEnFlux
double gOptEnergyHNL
void PrintSyntax(void)
string kDefOptSTune
string lunits
string geom
enum genie::hnl::EHNLDecayMode HNLDecayMode_t
string AsString(genie::hnl::HNLDecayMode_t hnldm)
bool IsKinematicallyAllowed(genie::hnl::HNLDecayMode_t hnldm, double Mhnl)
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f=' *')
double UnitFromString(string u)
Definition UnitUtils.cxx:18

References genie::CmdLnArgParser::ArgAsDouble(), genie::CmdLnArgParser::ArgAsInt(), genie::CmdLnArgParser::ArgAsLong(), genie::CmdLnArgParser::ArgAsString(), genie::utils::hnl::AsString(), geom, genie::AlgFactory::GetAlgorithm(), genie::hnl::Decayer::GetHNLMass(), gOptDecayMode, gOptEnergyHNL, gOptEvFilePrefix, gOptGeomLUnits, gOptIsMonoEnFlux, gOptMassHNL, gOptNev, gOptRanSeed, gOptRootGeom, gOptRootGeomTopVol, gOptRunNu, gOptUsingRootGeom, genie::AlgFactory::Instance(), genie::RunOpt::Instance(), genie::utils::hnl::IsKinematicallyAllowed(), kDefOptEvFilePrefix, kDefOptGeomLUnits, kDefOptSTune, LOG, lunits, genie::CmdLnArgParser::OptionExists(), pDEBUG, pFATAL, pINFO, pNOTICE, genie::utils::print::PrintFramedMesg(), PrintSyntax(), pWARN, genie::RunOpt::ReadFromCommandLine(), and genie::utils::units::UnitFromString().

Referenced by main().

◆ HNLGenerator()

const EventRecordVisitorI * HNLGenerator ( void )

Definition at line 729 of file gBeamHNLEvGen.cxx.

730{
731 //string sname = "genie::EventGenerator";
732 //string sconfig = "BeamHNL";
734
735 LOG("gevgen_hnl", pINFO)
736 << "Instantiating HNL generator.";
737
738 const Algorithm * algmcgen = algf->GetAlgorithm(kDefOptSName, kDefOptSConfig);
739 LOG("gevgen_hnl", pDEBUG)
740 << "Got algorithm " << kDefOptSName.c_str() << "/" << kDefOptSConfig.c_str();;
741
742 const EventRecordVisitorI * mcgen =
743 dynamic_cast< const EventRecordVisitorI * >( algmcgen );
744 if(!mcgen) {
745 LOG("gevgen_hnl", pFATAL) << "Couldn't instantiate the HNL generator";
746 gAbortingInErr = true;
747 exit(1);
748 }
749
750 LOG("gevgen_hnl", pINFO)
751 << "HNL generator instantiated successfully.";
752
753 return mcgen;
754}
The GENIE Algorithm Factory.
Definition AlgFactory.h:39
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the 'Visito...
string kDefOptSName
string kDefOptSConfig
bool gAbortingInErr
Definition Messenger.cxx:34

References genie::gAbortingInErr, genie::AlgFactory::GetAlgorithm(), genie::AlgFactory::Instance(), kDefOptSConfig, kDefOptSName, LOG, pDEBUG, pFATAL, and pINFO.

Referenced by main().

◆ main()

int main ( int argc,
char ** argv )

Definition at line 219 of file gBeamHNLEvGen.cxx.

220{
221 // suppress ROOT Info messages
222 gErrorIgnoreLevel = kWarning;
223
224 // Parse command line arguments
225 GetCommandLineArgs(argc,argv);
226
227 // Init messenger and random number seed
228 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
230
231 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
232
233 // Get the HNL generator first to load config
234 // config loaded upon instantiation of HNLGenerator algorithm
235 // calls LoadConfig() of each sub-alg
236 __attribute__((unused)) const EventRecordVisitorI * mcgen = HNLGenerator();
237 const Algorithm * algFluxCreator = AlgFactory::Instance()->GetAlgorithm("genie::hnl::FluxCreator", "Default");
238 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
239 const Algorithm * algVtxGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::VertexGenerator", "Default");
240
241 const FluxCreator * fluxCreator = dynamic_cast< const FluxCreator * >( algFluxCreator );
242 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
243 const VertexGenerator * vtxGen = dynamic_cast< const VertexGenerator * >( algVtxGen );
244
245 //string confString = kDefOptSName + "/" + kDefOptSConfig;
246 string confString = kDefOptSConfig;
247
248 CoMLifetime = hnlgen->GetHNLLifetime(); // GeV^{-1}
249
250 std::vector< double > U4l2s = hnlgen->GetHNLCouplings();
251 gOptECoupling = U4l2s.at(0);
252 gOptMCoupling = U4l2s.at(1);
253 gOptTCoupling = U4l2s.at(2);
254 gOptIsMajorana = hnlgen->IsHNLMajorana();
255
256 std::string stIntChannels = hnlgen->GetHNLInterestingChannels(); int iChan = -1;
257 if( gOptIntChannels.size() > 0 ) gOptIntChannels.clear();
258 while( stIntChannels.size() > 0 ){ // read channels from right (lowest mass) to left (highest mass)
259 iChan++;
260 HNLDecayMode_t md = static_cast< HNLDecayMode_t >( iChan );
261 std::string tmpSt = stIntChannels.substr( stIntChannels.size()-1, stIntChannels.size() );
262 if( std::strcmp( tmpSt.c_str(), "1" ) == 0 )
263 gOptIntChannels.emplace_back( md );
264
265 stIntChannels.erase( stIntChannels.end()-1, stIntChannels.end() );
266 }
267
268 //gOptIntChannels = confIntChan;
269
270 // Initialize an Ntuple Writer to save GHEP records into a TTree
272 ntpw.CustomizeFilenamePrefix(gOptEvFilePrefix);
273 ntpw.Initialize();
274
275 // add flux info to the tree
276 FluxContainer gnmf;
277 if( gOptIsUsingDk2nu ) {
278 // fill the flux object with nonsense to start with
279 FluxContainer * ptGnmf = new FluxContainer();
280 gnmf = *ptGnmf;
281 delete ptGnmf;
282 FillFluxNonsense( gnmf );
283 TBranch * flux = ntpw.EventTree()->Branch( "flux",
284 "genie::hnl::FluxContainer",
285 &gnmf, 32000, 1 );
286 flux->SetAutoDelete(kFALSE);
287 }
288
289 LOG("gevgen_hnl", pNOTICE)
290 << "Initialised Ntuple Writer";
291
292 // add another few branches to tree.
293 ntpw.EventTree()->Branch("hnl_mass", &gOptMassHNL, "gOptMassHNL/D");
294 ntpw.EventTree()->Branch("hnl_coup_e", &gOptECoupling, "gOptECoupling/D");
295 ntpw.EventTree()->Branch("hnl_coup_m", &gOptMCoupling, "gOptMCoupling/D");
296 ntpw.EventTree()->Branch("hnl_coup_t", &gOptTCoupling, "gOptTCoupling/D");
297 ntpw.EventTree()->Branch("hnl_ismaj", &gOptIsMajorana, "gOptIsMajorana/I");
298
299 // Create a MC job monitor for a periodically updated status file
300 GMCJMonitor mcjmonitor(gOptRunNu);
301 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
302
303 LOG("gevgen_hnl", pNOTICE)
304 << "Initialised MC job monitor";
305
306 // Set GHEP print level
307 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
308
309#ifdef __CAN_USE_ROOT_GEOM__
310 // Read geometry bounding box - for vertex position generation
311 if( gOptUsingRootGeom ){
312 InitBoundingBox();
313 }
314#endif // #ifdef __CAN_USE_ROOT_GEOM__
315
316 if( !gOptIsMonoEnFlux ){
317 LOG( "gevgen_hnl", pWARN )
318 << "Using input flux files. These are *flat dk2nu-like ROOT trees, so far...*";
319
320 }
321 // Event loop
322 int iflux = (gOptFirstEvent < 0) ? 0 : gOptFirstEvent; int ievent = iflux;
323 int maxFluxEntries = -1;
324 fluxCreator->SetInputFluxPath( gOptFluxFilePath );
325 fluxCreator->SetGeomFile( gOptRootGeom );
326 fluxCreator->SetFirstFluxEntry( iflux );
327 vtxGen->SetGeomFile( gOptRootGeom );
328
329 bool tooManyEntries = false;
330 while (1)
331 {
332 if( tooManyEntries ){
333 if( gOptNev >= 10000 ){
334 if( (ievent-gOptFirstEvent) % (gOptNev / 1000) == 0 ){
335 int irat = (iflux-gOptFirstEvent) / ( gOptNev / 1000 );
336 std::cerr << 0.1 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
337 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
338 }
339 } else if( gOptNev >= 100 ) {
340 if( (ievent-gOptFirstEvent) % (gOptNev / 10) == 0 ){
341 int irat = (iflux-gOptFirstEvent) / ( gOptNev / 10 );
342 std::cerr << 10.0 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
343 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
344 }
345 }
346 } else {
347 if( gOptNev >= 10000 ){
348 if( (ievent-gOptFirstEvent) % (gOptNev / 1000) == 0 ){
349 int irat = (ievent-gOptFirstEvent) / ( gOptNev / 1000 );
350 std::cerr << 0.1 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
351 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
352 }
353 } else if( gOptNev >= 100 ) {
354 if( (ievent-gOptFirstEvent) % (gOptNev / 10) == 0 ){
355 int irat = (ievent-gOptFirstEvent) / ( gOptNev / 10 );
356 std::cerr << 10.0 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
357 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
358 }
359 }
360 }
361
362 if( tooManyEntries && ((iflux-gOptFirstEvent) == gOptNev) ) break;
363 else if( (ievent-gOptFirstEvent) == gOptNev ) break;
364
365 if( ievent < gOptFirstEvent ){ ievent++; continue; }
366
367 assert( ievent >= gOptFirstEvent && gOptFirstEvent >= 0 );
368
369 LOG("gevgen_hnl", pNOTICE)
370 << " *** Generating event............ " << (ievent-gOptFirstEvent);
371
372 EventRecord * event = new EventRecord;
373 event->SetWeight(1.0);
374 event->SetProbability( CoMLifetime );
375 event->SetXSec( iflux ); // will be overridden, use as handy container
376 evWeight = 1.0;
377
378 if( !gOptIsMonoEnFlux ){
379 fluxCreator->ProcessEventRecord( event );
380
381 // fluxCreator->ProcessEventRecord now tells us how many entries there are
382 maxFluxEntries = fluxCreator->GetNFluxEntries();
383 if( gOptNev > maxFluxEntries - gOptFirstEvent ){
384 LOG( "gevgen_hnl", pWARN )
385 << "You have asked for " << gOptNev << " events, but only provided "
386 << maxFluxEntries - gOptFirstEvent << " flux entries. Truncating events to " << maxFluxEntries - gOptFirstEvent << ".";
387 gOptNev = maxFluxEntries - gOptFirstEvent;
388 tooManyEntries = true;
389 }
390 if( iflux >= maxFluxEntries - 1 ){
391 LOG( "gevgen_hnl", pWARN )
392 << "Reached end of flux input (iflux = " << iflux << "), about to stop.";
393 break;
394 }
395
396 FluxContainer retGnmf = fluxCreator->RetrieveFluxInfo();
397 FillFlux( gnmf, retGnmf );
398
399 // check to see if this was nonsense
400 if( ! event->Particle(0) ){ iflux++; delete event; continue; }
401
402 gOptEnergyHNL = event->Particle(0)->GetP4()->E();
403 iflux++;
404 } else { // monoenergetic HNL. Add it with energy and momentum pointing on z axis
405
406 assert( gOptEnergyHNL > gOptMassHNL );
407 double HNLP = std::sqrt( gOptEnergyHNL*gOptEnergyHNL - gOptMassHNL*gOptMassHNL );
408 TLorentzVector probeP4( 0.0, 0.0, HNLP, gOptEnergyHNL );
409 TLorentzVector v4( 0.0, 0.0, 0.0, 0.0 );
410 GHepParticle ptHNL( kPdgHNL, kIStInitialState, -1, -1, -1, -1, probeP4, v4 );
411 event->AddParticle( ptHNL );
412
413 }
414 assert( gOptEnergyHNL > gOptMassHNL );
415
416 int hpdg = genie::kPdgHNL;
417 int typeMod = 1;
418 if( gOptIsMajorana ) typeMod = ( rnd->RndGen().Uniform( 0.0, 1.0 ) >= 0.5 ) ? -1 : 1;
419 else if( event->Particle(0)->Pdg() > 0 ) typeMod = 1;
420 else if( event->Particle(0)->Pdg() < 0 ) typeMod = -1;
421
422 // int target = SelectInitState();
423 int decay = (int) gOptDecayMode;
424
425 assert( gOptECoupling >= 0.0 && gOptMCoupling >= 0.0 && gOptTCoupling >= 0.0 );
426
427 // RETHERE assuming all these HNL have K+- parent. This is wrong
428 // (but not very wrong for interesting masses)
429 SimpleHNL sh( "HNL", ievent, hpdg, genie::kPdgKP,
431
432 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
433 CoMLifetime = sh.GetCoMLifetime();
434
435 if( gOptDecayMode == kHNLDcyNull ){ // select from available modes
436 LOG("gevgen_hnl", pNOTICE)
437 << "No decay mode specified - sampling from all available modes.";
438
439 std::vector< HNLDecayMode_t > * intChannels = &gOptIntChannels;
440
441 decay = SelectDecayMode( intChannels, sh );
442 }
443
444 Interaction * interaction = Interaction::HNL(typeMod * genie::kPdgHNL, gOptEnergyHNL, decay);
445
446 if( event->Particle(0) ){ // we have an HNL with definite momentum, so let's set it now
447 interaction->InitStatePtr()->SetProbeP4( *(event->Particle(0)->P4()) );
448 interaction->InitStatePtr()->SetProbePdg( event->Particle(0)->Pdg() );
449 LOG( "gevgen_hnl", pDEBUG )
450 << "\nsetting probe p4 = " << utils::print::P4AsString( event->Particle(0)->P4() );
451 }
452
453 event->AttachSummary(interaction);
454
455 // Simulate decay
456 hnlgen->ProcessEventRecord(event);
457
458 // add the FS 4-momenta to special branches
459 // Quite inelegant. Gets the job done, though
460 NTP_FS0_PDG = (event->Particle(1))->Pdg();
461 NTP_FS0_E = ((event->Particle(1))->P4())->E();
462 NTP_FS0_PX = ((event->Particle(1))->P4())->Px();
463 NTP_FS0_PY = ((event->Particle(1))->P4())->Py();
464 NTP_FS0_PZ = ((event->Particle(1))->P4())->Pz();
465 NTP_FS1_PDG = (event->Particle(2))->Pdg();
466 NTP_FS1_E = ((event->Particle(2))->P4())->E();
467 NTP_FS1_PX = ((event->Particle(2))->P4())->Px();
468 NTP_FS1_PY = ((event->Particle(2))->P4())->Py();
469 NTP_FS1_PZ = ((event->Particle(2))->P4())->Pz();
470 if( event->Particle(3) ){
471 NTP_FS2_PDG = (event->Particle(3))->Pdg();
472 NTP_FS2_E = ((event->Particle(3))->P4())->E();
473 NTP_FS2_PX = ((event->Particle(3))->P4())->Px();
474 NTP_FS2_PY = ((event->Particle(3))->P4())->Py();
475 NTP_FS2_PZ = ((event->Particle(3))->P4())->Pz();
476 }
477 else{
478 NTP_FS2_PDG = 0;
479 NTP_FS2_E = 0.0;
480 NTP_FS2_PX = 0.0;
481 NTP_FS2_PY = 0.0;
482 NTP_FS2_PZ = 0.0;
483 }
484
485 // Generate (or read) a position for the decay vertex
486 // also currently handles the geometrical weight
487 TLorentzVector x4mm;
488 if( gOptUsingRootGeom ){
489 TLorentzVector * p4HNL = event->Particle(0)->GetP4();
490 NTP_IS_E = p4HNL->E(); NTP_IS_PX = p4HNL->Px(); NTP_IS_PY = p4HNL->Py(); NTP_IS_PZ = p4HNL->Pz();
491 vtxGen->ProcessEventRecord( event );
492 x4mm = *(event->Vertex());
493 } else {
494 x4mm = GeneratePosition( event );
495 }
496
497 // update weight to scale for couplings, inhibited decays
498 // acceptance is already handled in FluxCreator
499 // geometry handled in VertexGenerator
500 evWeight = event->Weight();
502 evWeight *= 1.0 / decayMod;
503 event->SetWeight( evWeight / 1.0e+20 ); // in units of 1e+20 POT
504
505 // store the decayMod in the second daughter's polarisation
506 event->Particle(2)->SetPolarization( 1.0, decayMod );
507
508 // why does InitState show the wrong p4 here?
509 interaction->InitStatePtr()->SetProbeP4( *(event->Particle(0)->P4()) );
510
511 LOG("gevgen_hnl", pDEBUG) << "Weight = " << evWeight;
512
513 LOG("gevgen_hnl", pINFO)
514 << "Generated event: " << *event;
515
516 // Add event at the output ntuple, refresh the mc job monitor & clean-up
517 ntpw.AddEventRecord(ievent, event);
518 mcjmonitor.Update(ievent,event);
519
520 delete event;
521
522 ievent++;
523 } // event loop
524
525 // Save the generated event tree & close the output file
526 ntpw.Save();
527
528 LOG("gevgen_hnl", pNOTICE) << "Done!";
529
530 return 0;
531}
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition EventRecord.h:37
STDHEP-like event record entry that can fit a particle or a nucleus.
virtual void SetWeight(double wght)
Definition GHepRecord.h:130
static void SetPrintLevel(int print_level)
Simple class to create & update MC job status files and env. vars. This is used to be able to keep tr...
Definition GMCJMonitor.h:31
void SetProbePdg(int pdg_code)
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition Interaction.h:56
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
static Interaction * HNL(int probe, double E=0, int decayed_mode=-1)
A utility class to facilitate creating the GENIE MC Ntuple from the output GENIE GHEP event records.
Definition NtpWriter.h:39
TRandom3 & RndGen(void) const
rnd number generator for generic usage
Definition RandomGen.h:80
std::string GetHNLInterestingChannels() const
void ProcessEventRecord(GHepRecord *event) const
double GetHNLLifetime() const
bool IsHNLMajorana() const
std::vector< double > GetHNLCouplings() const
A GENIE flux container specific for HNL containers. Based on the dk2nu flux paradigm and genie::flux:...
Calculates HNL production kinematics & production vertex. Is a concrete implementation of the FluxRec...
void SetFirstFluxEntry(int iFirst) const
void ProcessEventRecord(GHepRecord *event_rec) const
FluxContainer RetrieveFluxInfo() const
void SetGeomFile(string geomfile) const
void SetInputFluxPath(std::string finpath) const
Heavy Neutral Lepton decay vertex generator given production vertex and momentum ***.
void SetGeomFile(std::string geomfile) const
void ProcessEventRecord(GHepRecord *event_rec) const
NtpMCFormat_t kDefOptNtpFormat
double gOptECoupling
TLorentzVector GeneratePosition(GHepRecord *event)
bool gOptIsMajorana
double NTP_FS1_E
double NTP_IS_E
double NTP_FS1_PX
double NTP_FS0_PZ
double NTP_FS1_PY
double gOptMCoupling
double NTP_FS0_PY
std::vector< HNLDecayMode_t > gOptIntChannels
double decayMod
const EventRecordVisitorI * HNLGenerator(void)
double NTP_FS0_PX
double evWeight
double NTP_IS_PX
double NTP_FS0_E
double NTP_FS2_PZ
int SelectDecayMode(std::vector< HNLDecayMode_t > *intChannels, SimpleHNL sh)
double NTP_FS2_PY
int NTP_FS1_PDG
int NTP_FS2_PDG
void GetCommandLineArgs(int argc, char **argv)
double NTP_FS2_PX
double NTP_IS_PY
double CoMLifetime
double NTP_FS2_E
int NTP_FS0_PDG
double gOptTCoupling
double NTP_FS1_PZ
double NTP_IS_PZ
GENIE flux drivers.
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
string P4AsString(const TLorentzVector *p)
@ kIStInitialState
Definition GHepStatus.h:29
const int kPdgKP
Definition PDGCodes.h:172
const int kPdgHNL
Definition PDGCodes.h:224

References genie::NtpWriter::AddEventRecord(), CoMLifetime, genie::NtpWriter::CustomizeFilenamePrefix(), decayMod, genie::NtpWriter::EventTree(), evWeight, GeneratePosition(), genie::AlgFactory::GetAlgorithm(), genie::hnl::SimpleHNL::GetCoMLifetime(), GetCommandLineArgs(), genie::hnl::Decayer::GetHNLCouplings(), genie::hnl::Decayer::GetHNLInterestingChannels(), genie::hnl::Decayer::GetHNLLifetime(), genie::hnl::FluxCreator::GetNFluxEntries(), genie::hnl::SimpleHNL::GetValidChannels(), gOptDecayMode, gOptECoupling, gOptEnergyHNL, gOptEvFilePrefix, gOptIntChannels, gOptIsMajorana, gOptIsMonoEnFlux, gOptMassHNL, gOptMCoupling, gOptNev, gOptRanSeed, gOptRootGeom, gOptRunNu, gOptTCoupling, gOptUsingRootGeom, genie::Interaction::HNL(), HNLGenerator(), genie::NtpWriter::Initialize(), genie::Interaction::InitStatePtr(), genie::AlgFactory::Instance(), genie::RandomGen::Instance(), genie::RunOpt::Instance(), genie::hnl::Decayer::IsHNLMajorana(), kDefOptNtpFormat, kDefOptSConfig, genie::hnl::kHNLDcyNull, genie::kIStInitialState, genie::kPdgHNL, genie::kPdgKP, LOG, genie::utils::app_init::MesgThresholds(), NTP_FS0_E, NTP_FS0_PDG, NTP_FS0_PX, NTP_FS0_PY, NTP_FS0_PZ, NTP_FS1_E, NTP_FS1_PDG, NTP_FS1_PX, NTP_FS1_PY, NTP_FS1_PZ, NTP_FS2_E, NTP_FS2_PDG, NTP_FS2_PX, NTP_FS2_PY, NTP_FS2_PZ, NTP_IS_E, NTP_IS_PX, NTP_IS_PY, NTP_IS_PZ, genie::utils::print::P4AsString(), pDEBUG, pINFO, pNOTICE, genie::hnl::Decayer::ProcessEventRecord(), genie::hnl::FluxCreator::ProcessEventRecord(), genie::hnl::VertexGenerator::ProcessEventRecord(), pWARN, genie::utils::app_init::RandGen(), genie::hnl::FluxCreator::RetrieveFluxInfo(), genie::RandomGen::RndGen(), genie::NtpWriter::Save(), SelectDecayMode(), genie::hnl::FluxCreator::SetFirstFluxEntry(), genie::hnl::FluxCreator::SetGeomFile(), genie::hnl::VertexGenerator::SetGeomFile(), genie::hnl::FluxCreator::SetInputFluxPath(), genie::GHepRecord::SetPrintLevel(), genie::InitialState::SetProbeP4(), genie::InitialState::SetProbePdg(), genie::GMCJMonitor::SetRefreshRate(), genie::GHepRecord::SetWeight(), and genie::GMCJMonitor::Update().

◆ PrintSyntax()

void PrintSyntax ( void )

Definition at line 1065 of file gBeamHNLEvGen.cxx.

1066{
1067 LOG("gevgen_hnl", pFATAL)
1068 << "\n **Syntax**"
1069 << "\n gevgen_hnl [-h] "
1070 << "\n [-r run#]"
1071 << "\n -n n_of_events"
1072 << "\n -f path/to/flux/files"
1073 << "\n [-E hnl_energy]"
1074 << "\n [--firstEvent first_event_for_dk2nu_readin]"
1075 << "\n [-m decay_mode]"
1076 << "\n [-g geometry (ROOT file)]"
1077 << "\n [-L length_units_at_geom]"
1078 << "\n [-o output_event_file_prefix]"
1079 << "\n [--seed random_number_seed]"
1080 << "\n [--message-thresholds xml_file]"
1081 << "\n [--event-record-print-level level]"
1082 << "\n [--mc-job-status-refresh-rate rate]"
1083 << "\n"
1084 << " Please also read the detailed documentation at http://www.genie-mc.org"
1085 << " or look at the source code: $GENIE/src/Apps/gBeamHNLEvGen.cxx"
1086 << "\n";
1087}

References LOG, and pFATAL.

Referenced by GetCommandLineArgs().

◆ SelectDecayMode()

int SelectDecayMode ( std::vector< HNLDecayMode_t > * intChannels,
SimpleHNL sh )

Definition at line 756 of file gBeamHNLEvGen.cxx.

757{
758 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
759
760 // set CoM lifetime now if unset
761 if( CoMLifetime < 0.0 ){
763 }
764
765 std::vector< HNLDecayMode_t > intAndValidChannels;
766 for( std::vector< HNLDecayMode_t >::iterator it = intChannels->begin(); it != intChannels->end(); ++it ){
767 HNLDecayMode_t mode = *it;
768
769 //check if this is a valid mode
770 if( !utils::hnl::IsKinematicallyAllowed( mode, gOptMassHNL ) ) continue;
771
772 intAndValidChannels.emplace_back( mode );
773 auto mapG = gammaMap.find( mode );
774 double theGamma = mapG->second;
775 LOG("gevgen_hnl", pDEBUG)
776 << "For mode " << utils::hnl::AsString( mode ) << " gamma = " << theGamma;
777 }
778
779 if( intAndValidChannels.size() == 0 ){ // all the modes picked by user are too heavy. Abort.
780 LOG( "gevgen_hnl", pFATAL )
781 << "None of the channels specified as interesting are kinematically allowed. Please either increase the HNL mass or change interesting channels in config.";
782 exit(1);
783 }
784
785 std::map< HNLDecayMode_t, double > intMap =
786 selector::SetInterestingChannels( intAndValidChannels, gammaMap );
787
788 sh.SetInterestingChannels( intMap );
789
790 // update fraction of total decay width that is not in inhibited channels
791 double gammaAll = 0.0, gammaInt = 0.0;
792 for( std::map< HNLDecayMode_t, double >::const_iterator itall = gammaMap.begin() ;
793 itall != gammaMap.end() ; ++itall ){
794 gammaAll += (*itall).second;
795 }
796 for( std::map< HNLDecayMode_t, double >::iterator itint = intMap.begin() ;
797 itint != intMap.end() ; ++itint ){
798 gammaInt += (*itint).second;
799 }
800 assert( gammaInt > 0.0 && gammaAll >= gammaInt );
801 decayMod = gammaInt / gammaAll;
802
803 // get probability that channels in intAndValidChannels will be selected
804 std::map< HNLDecayMode_t, double > PMap =
806
807 // now do a random throw
809 double ranThrow = rnd->RndGen().Uniform( 0., 1. ); // HNL's fate is sealed.
810
811 HNLDecayMode_t selectedDecayChan =
812 selector::SelectChannelInclusive( PMap, ranThrow );
813
814 int decay = ( int ) selectedDecayChan;
815 return decay;
816}
double GetCoMLifetime()
Definition SimpleHNL.h:96
void SetInterestingChannels(const std::map< genie::hnl::HNLDecayMode_t, double > gammaMap)
Definition SimpleHNL.h:276
std::map< genie::hnl::HNLDecayMode_t, double > GetValidChannels()
Definition SimpleHNL.h:154
genie::hnl::HNLDecayMode_t SelectChannelInclusive(std::map< genie::hnl::HNLDecayMode_t, double > Pmap, double ranThrow)
std::map< genie::hnl::HNLDecayMode_t, double > SetInterestingChannels(std::vector< genie::hnl::HNLDecayMode_t > intChannels, std::map< genie::hnl::HNLDecayMode_t, double > gammaMap)
std::map< genie::hnl::HNLDecayMode_t, double > GetProbabilities(std::map< genie::hnl::HNLDecayMode_t, double > gammaMap)

References genie::utils::hnl::AsString(), CoMLifetime, decayMod, genie::hnl::SimpleHNL::GetCoMLifetime(), genie::hnl::selector::GetProbabilities(), genie::hnl::SimpleHNL::GetValidChannels(), gOptMassHNL, genie::RandomGen::Instance(), genie::utils::hnl::IsKinematicallyAllowed(), LOG, pDEBUG, pFATAL, genie::RandomGen::RndGen(), genie::hnl::selector::SelectChannelInclusive(), genie::hnl::selector::SetInterestingChannels(), and genie::hnl::SimpleHNL::SetInterestingChannels().

Referenced by main().

Variable Documentation

◆ CoMLifetime

◆ decayMod

double decayMod = 1.0

Definition at line 214 of file gBeamHNLEvGen.cxx.

Referenced by main(), and SelectDecayMode().

◆ evWeight

double evWeight = 1.0

Definition at line 216 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ fdx

double fdx = 0

Definition at line 198 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ fdy

double fdy = 0

Definition at line 199 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ fdz

double fdz = 0

Definition at line 200 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ fox

double fox = 0

Definition at line 201 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ foy

double foy = 0

Definition at line 202 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ foz

double foz = 0

Definition at line 203 of file gBeamHNLEvGen.cxx.

Referenced by GeneratePosition().

◆ gOptDecayMode

HNLDecayMode_t gOptDecayMode = kHNLDcyNull

◆ gOptECoupling

double gOptECoupling = -1

Definition at line 167 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ gOptEnergyHNL

double gOptEnergyHNL = -1

◆ gOptEvFilePrefix

string gOptEvFilePrefix = kDefOptEvFilePrefix

Definition at line 184 of file gBeamHNLEvGen.cxx.

◆ gOptGeomLUnits

double gOptGeomLUnits = 0

Definition at line 194 of file gBeamHNLEvGen.cxx.

◆ gOptIntChannels

std::vector< HNLDecayMode_t > gOptIntChannels

Definition at line 182 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ gOptIsMajorana

bool gOptIsMajorana = false

Definition at line 171 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ gOptIsMonoEnFlux

bool gOptIsMonoEnFlux = true

Definition at line 179 of file gBeamHNLEvGen.cxx.

Referenced by GetCommandLineArgs(), and main().

◆ gOptMassHNL

double gOptMassHNL = -1

◆ gOptMCoupling

double gOptMCoupling = -1

Definition at line 168 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ gOptNev

int gOptNev = 10

Definition at line 163 of file gBeamHNLEvGen.cxx.

◆ gOptRanSeed

long int gOptRanSeed = -1

Definition at line 195 of file gBeamHNLEvGen.cxx.

◆ gOptRootGeom

string gOptRootGeom

Definition at line 186 of file gBeamHNLEvGen.cxx.

◆ gOptRootGeomTopVol

string gOptRootGeomTopVol = ""

Definition at line 193 of file gBeamHNLEvGen.cxx.

◆ gOptRunNu

Long_t gOptRunNu = 1000

Definition at line 162 of file gBeamHNLEvGen.cxx.

◆ gOptTCoupling

double gOptTCoupling = -1

Definition at line 169 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ gOptUsingRootGeom

bool gOptUsingRootGeom = false

Definition at line 185 of file gBeamHNLEvGen.cxx.

◆ kDefOptEvFilePrefix

string kDefOptEvFilePrefix = "gntp"

Definition at line 154 of file gBeamHNLEvGen.cxx.

◆ kDefOptFluxFilePath

string kDefOptFluxFilePath = "./input-flux.root"

Definition at line 155 of file gBeamHNLEvGen.cxx.

◆ kDefOptGeomDUnits

string kDefOptGeomDUnits = "g_cm3"

Definition at line 152 of file gBeamHNLEvGen.cxx.

◆ kDefOptGeomLUnits

string kDefOptGeomLUnits = "mm"

Definition at line 151 of file gBeamHNLEvGen.cxx.

◆ kDefOptNtpFormat

NtpMCFormat_t kDefOptNtpFormat = kNFGHEP

Definition at line 153 of file gBeamHNLEvGen.cxx.

◆ kDefOptSConfig

string kDefOptSConfig = "BeamHNL"

Definition at line 158 of file gBeamHNLEvGen.cxx.

Referenced by HNLGenerator(), and main().

◆ kDefOptSName

string kDefOptSName = "genie::EventGenerator"

Definition at line 157 of file gBeamHNLEvGen.cxx.

Referenced by HNLGenerator(), and main().

◆ kDefOptSTune

string kDefOptSTune = "GHNL20_00a_00_000"

Definition at line 159 of file gBeamHNLEvGen.cxx.

Referenced by GetCommandLineArgs().

◆ NTP_FS0_E

double NTP_FS0_E = 0.

Definition at line 206 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS0_PDG

int NTP_FS0_PDG = 0

Definition at line 209 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS0_PX

double NTP_FS0_PX = 0.

Definition at line 206 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS0_PY

double NTP_FS0_PY = 0.

Definition at line 206 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS0_PZ

double NTP_FS0_PZ = 0.

Definition at line 206 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS1_E

double NTP_FS1_E = 0.

Definition at line 207 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS1_PDG

int NTP_FS1_PDG = 0

Definition at line 209 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS1_PX

double NTP_FS1_PX = 0.

Definition at line 207 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS1_PY

double NTP_FS1_PY = 0.

Definition at line 207 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS1_PZ

double NTP_FS1_PZ = 0.

Definition at line 207 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS2_E

double NTP_FS2_E = 0.

Definition at line 208 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS2_PDG

int NTP_FS2_PDG = 0

Definition at line 209 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS2_PX

double NTP_FS2_PX = 0.

Definition at line 208 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS2_PY

double NTP_FS2_PY = 0.

Definition at line 208 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_FS2_PZ

double NTP_FS2_PZ = 0.

Definition at line 208 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_IS_E

double NTP_IS_E = 0.

Definition at line 205 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_IS_PX

double NTP_IS_PX = 0.

Definition at line 205 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_IS_PY

double NTP_IS_PY = 0.

Definition at line 205 of file gBeamHNLEvGen.cxx.

Referenced by main().

◆ NTP_IS_PZ

double NTP_IS_PZ = 0.

Definition at line 205 of file gBeamHNLEvGen.cxx.

Referenced by main().