GENIEGenerator
Loading...
Searching...
No Matches
gBeamHNLParticleGun.cxx File Reference
Include dependency graph for gBeamHNLParticleGun.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 * GenerateOriginPosition (GHepRecord *event)
TLorentzVector * GenerateOriginMomentum (GHepRecord *event)
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 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
int gOptHNLKind = -1
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 evProdVtx [3] = {0.0, 0.0, 0.0}
double evWeight = 1.0

Function Documentation

◆ GenerateOriginMomentum()

TLorentzVector * GenerateOriginMomentum ( GHepRecord * event)

Definition at line 507 of file gBeamHNLParticleGun.cxx.

508{
509 double E = gOptEnergyHNL;
510 double M = gOptMassHNL;
511 double pMag = std::sqrt( E*E - M*M );
512
513 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
514 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
515
516 // first map directional cosines to unit sphere
517 /*
518 double cx = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cx" );
519 double cy = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cy" );
520 double cz = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cz" );
521 */
522 std::vector< double > PGDirection = hnlgen->GetPGunDirection();
523 double cx = PGDirection.at(0), cy = PGDirection.at(1), cz = PGDirection.at(2);
524 double c2 = std::sqrt( cx*cx + cy*cy + cz*cz );
525 assert( c2 > 0.0 );
526
527 cx *= 1.0 / c2; cy *= 1.0 / c2; cz *= 1.0 / c2;
528
529 double theta = TMath::ACos( cz / c2 );
530 double phi = ( TMath::Sin( theta ) != 0.0 ) ? TMath::ACos( cx / ( c2 * TMath::Sin( theta ) ) ) : 0.0;
531
532 // apply uniform random deviation
533 std::vector< double > PGunDeviation = hnlgen->GetPGunDeviation();
534 /*
535 double dthetaDeg = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-DTheta" );
536 double dphiDeg = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-DPhi" );
537 */
538 double dthetaDeg = PGunDeviation.at(0), dphiDeg = PGunDeviation.at(1);
539 double dThetaMax = dthetaDeg * constants::kPi / 180.0;
540 double dPhiMax = dphiDeg * constants::kPi / 180.0;
541
543 double dTheta = rng->RndGen().Uniform( -dThetaMax, dThetaMax );
544 double dPhi = rng->RndGen().Uniform( -dPhiMax, dPhiMax );
545
546 std::ostringstream asts;
547 asts
548 << "Output details for the momentum:"
549 << "\nDirectional cosines: ( " << cx << ", " << cy << ", " << cz << " )"
550 << "\nMapped to ( " << theta * 180.0 / constants::kPi
551 << ", " << phi * 180.0 / constants::kPi << " ) [deg] on unit sphere"
552 << "\nApplied deviation: dTheta = " << dTheta * 180.0 / constants::kPi
553 << ", dPhi = " << dPhi * 180.0 / constants::kPi << " [deg]";
554
555 // make momentum
556 theta += dTheta; phi += dPhi;
557 TLorentzVector * p4HNL = new TLorentzVector();
558 p4HNL->SetPxPyPzE( pMag * TMath::Sin( theta ) * TMath::Cos( phi ),
559 pMag * TMath::Sin( theta ) * TMath::Sin( phi ),
560 pMag * TMath::Cos( theta ), E );
561
562 asts << "\nMomentum = " << utils::print::P4AsString( p4HNL );
563
564 Interaction * interaction = event->Summary();
565 interaction->InitStatePtr()->SetProbeP4( *p4HNL );
566
567 LOG( "gevgen_pghnl", pDEBUG ) << asts.str();
568
569 //delete rng;
570 return p4HNL;
571}
#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
const Algorithm * GetAlgorithm(const AlgId &algid)
static AlgFactory * Instance()
Algorithm abstract base class.
Definition Algorithm.h:54
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition Interaction.h:56
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
TRandom3 & RndGen(void) const
rnd number generator for generic usage
Definition RandomGen.h:80
Heavy Neutral Lepton final-state product generator.
Definition HNLDecayer.h:41
std::vector< double > GetPGunDirection() const
std::vector< double > GetPGunDeviation() const
double gOptMassHNL
double gOptEnergyHNL
string P4AsString(const TLorentzVector *p)

References genie::AlgFactory::GetAlgorithm(), genie::hnl::Decayer::GetPGunDeviation(), genie::hnl::Decayer::GetPGunDirection(), gOptEnergyHNL, gOptMassHNL, genie::Interaction::InitStatePtr(), genie::AlgFactory::Instance(), genie::RandomGen::Instance(), genie::constants::kPi, LOG, genie::utils::print::P4AsString(), pDEBUG, genie::RandomGen::RndGen(), and genie::InitialState::SetProbeP4().

Referenced by main().

◆ GenerateOriginPosition()

TLorentzVector * GenerateOriginPosition ( GHepRecord * event)

Definition at line 573 of file gBeamHNLParticleGun.cxx.

574{
575 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
576 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
577
578 // get centre of box from config
579 std::vector< double > PGOrigin = hnlgen->GetPGunOrigin();
580 /*
581 double ox = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginX" );
582 double oy = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginY" );
583 double oz = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginZ" );
584 */
585 double ox = PGOrigin.at(0), oy = PGOrigin.at(1), oz = PGOrigin.at(2);
586
587 // allow uniform deviation
588 std::vector< double > PGDOrigin = hnlgen->GetPGunDOrigin();
589 /*
590 double Dxmax = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginDX" );
591 double Dymax = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginDY" );
592 double Dzmax = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginDZ" );
593 */
594 double Dxmax = PGDOrigin.at(0), Dymax = PGDOrigin.at(1), Dzmax = PGDOrigin.at(2);
595
597 double dx = rng->RndGen().Uniform( -Dxmax, Dxmax );
598 double dy = rng->RndGen().Uniform( -Dymax, Dymax );
599 double dz = rng->RndGen().Uniform( -Dzmax, Dzmax );
600
601 // make position
602 ox += dx; oy += dy; oz += dz;
603 TLorentzVector * x4HNL = new TLorentzVector();
604 x4HNL->SetXYZT( ox, oy, oz, 0.0 );
605
606 event->SetVertex( *x4HNL );
607 //delete rng;
608 return x4HNL;
609}
std::vector< double > GetPGunDOrigin() const
std::vector< double > GetPGunOrigin() const

References genie::AlgFactory::GetAlgorithm(), genie::hnl::Decayer::GetPGunDOrigin(), genie::hnl::Decayer::GetPGunOrigin(), genie::AlgFactory::Instance(), genie::RandomGen::Instance(), and genie::RandomGen::RndGen().

Referenced by main().

◆ GeneratePosition()

TLorentzVector GeneratePosition ( GHepRecord * event)

Definition at line 611 of file gBeamHNLParticleGun.cxx.

612{
613 if( gOptUsingRootGeom ){
614
615 const Algorithm * algVtxGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::VertexGenerator", "Default");
616
617 const VertexGenerator * vtxGen = dynamic_cast< const VertexGenerator * >( algVtxGen );
618 vtxGen->ProcessEventRecord( event );
619
620 TLorentzVector x4 = *(event->Vertex());
621 return x4;
622 }
623 else{
624 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
625 TRandom3 & rnd_geo = rnd->RndGeom();
626
627 double rndx = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
628 double rndy = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
629 double rndz = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
630
631 double t_gen = 0;
632 double x_gen = fox + rndx * fdx;
633 double y_gen = foy + rndy * fdy;
634 double z_gen = foz + rndz * fdz;
635
636 TLorentzVector pos(x_gen, y_gen, z_gen, t_gen);
637 return pos;
638 }
639
640 return TLorentzVector( 0, 0, 0, 0);
641}
TRandom3 & RndGeom(void) const
rnd number generator used by geometry drivers
Definition RandomGen.h:68
Heavy Neutral Lepton decay vertex generator given production vertex and momentum ***.
void ProcessEventRecord(GHepRecord *event_rec) const
bool gOptUsingRootGeom
double foy
double fdz
double foz
double fox
double fdy
double fdx

References fdx, fdy, fdz, fox, foy, foz, genie::AlgFactory::GetAlgorithm(), gOptUsingRootGeom, genie::AlgFactory::Instance(), genie::RandomGen::Instance(), genie::hnl::VertexGenerator::ProcessEventRecord(), and genie::RandomGen::RndGeom().

Referenced by main().

◆ GetCommandLineArgs()

void GetCommandLineArgs ( int argc,
char ** argv )

Definition at line 709 of file gBeamHNLParticleGun.cxx.

710{
711 LOG("gevgen_pghnl", pINFO) << "Parsing command line arguments";
712
713 // Parse run options for this app
714
715 CmdLnArgParser parser(argc,argv);
716
717 // help?
718 bool help = parser.OptionExists('h');
719 if(help) {
720 PrintSyntax();
721 exit(0);
722 }
723
724 // force the app to look at HNL tune by default
725 // if user passes --tune argument, look at the user input tune instead
726 char * expargv[ argc + 2 ];
727 bool didFindUserInputTune = false;
728 std::string stExtraTuneBit = kDefOptSTune;
729
730 if( parser.OptionExists("tune") ){
731 didFindUserInputTune = true;
732 stExtraTuneBit = parser.ArgAsString("tune");
733 LOG( "gevgen_hnl", pWARN )
734 << "Using input HNL tune " << parser.ArgAsString("tune");
735 } else {
736 LOG( "gevgen_hnl", pWARN )
737 << "Using default HNL tune " << kDefOptSTune;
738 }
739 // append this to argv
740 for( int iArg = 0; iArg < argc; iArg++ ){
741 expargv[iArg] = argv[iArg];
742 }
743 if( !didFindUserInputTune ){
744 char * chBit = const_cast< char * >( stExtraTuneBit.c_str() ); // ugh. Ugly.
745 std::string stune("--tune"); char * tBit = const_cast< char * >( stune.c_str() );
746 expargv[argc] = tBit;
747 expargv[argc+1] = chBit;
748 }
749
750 // Common run options.
751 int expargc = ( didFindUserInputTune ) ? argc : argc+2;
752 std::string stnull(""); char * nBit = const_cast< char * >( stnull.c_str() );
753 expargv[expargc] = nBit;
754
755 RunOpt::Instance()->ReadFromCommandLine(expargc,expargv);
756
757 // run number
758 if( parser.OptionExists('r') ) {
759 LOG("gevgen_pghnl", pDEBUG) << "Reading MC run number";
760 gOptRunNu = parser.ArgAsLong('r');
761 } else {
762 LOG("gevgen_pghnl", pDEBUG) << "Unspecified run number - Using default";
763 gOptRunNu = 1000;
764 } //-r
765
766 // number of events
767 if( parser.OptionExists('n') ) {
768 LOG("gevgen_pghnl", pDEBUG)
769 << "Reading number of events to generate";
770 gOptNev = parser.ArgAsInt('n');
771 } else {
772 LOG("gevgen_pghnl", pFATAL)
773 << "You need to specify the number of events";
774 PrintSyntax();
775 exit(0);
776 } //-n
777
778 // get HNL mass
779 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
780 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
781 gOptMassHNL = hnlgen->GetHNLMass();
782
783 // HNL decay mode
784 int mode = -1;
785 if( parser.OptionExists('m') ) {
786 LOG("gevgen_pghnl", pDEBUG)
787 << "Reading HNL decay mode";
788 mode = parser.ArgAsInt('m');
789 } else {
790 LOG("gevgen_pghnl", pINFO)
791 << "No decay mode specified - will sample from allowed decay modes";
792 } //-m
794
796 if(!allowed) {
797 LOG("gevgen_pghnl", pFATAL)
798 << "Specified decay is not allowed kinematically for the given HNL mass";
799 PrintSyntax();
800 exit(0);
801 }
802
803 //
804 // geometry
805 //
806
807 string geom = "";
808 string lunits;
809#ifdef __CAN_USE_ROOT_GEOM__
810 // string dunits;
811 if( parser.OptionExists('g') ) {
812 LOG("gevgen_pghnl", pDEBUG) << "Getting input geometry";
813 geom = parser.ArgAsString('g');
814
815 // is it a ROOT file that contains a ROOT geometry?
816 bool accessible_geom_file =
817 ! (gSystem->AccessPathName(geom.c_str()));
818 if (accessible_geom_file) {
820 gOptUsingRootGeom = true;
821 } else {
822 LOG("gevgen_pghnl", pFATAL)
823 << "Geometry option is not a ROOT file. Please use ROOT geom.";
824 PrintSyntax();
825 exit(1);
826 }
827 } else {
828 // LOG("gevgen_pghnl", pFATAL)
829 // << "No geometry option specified - Exiting";
830 // PrintSyntax();
831 // exit(1);
832 } //-g
833
835 // using a ROOT geometry - get requested geometry units
836
837 // length units:
838 if( parser.OptionExists('L') ) {
839 LOG("gevgen_pghnl", pDEBUG)
840 << "Checking for input geometry length units";
841 lunits = parser.ArgAsString('L');
842 } else {
843 LOG("gevgen_pghnl", pDEBUG) << "Using default geometry length units";
845 } // -L
846 // // density units:
847 // if( parser.OptionExists('D') ) {
848 // LOG("gevgen_pghnl", pDEBUG)
849 // << "Checking for input geometry density units";
850 // dunits = parser.ArgAsString('D');
851 // } else {
852 // LOG("gevgen_pghnl", pDEBUG) << "Using default geometry density units";
853 // dunits = kDefOptGeomDUnits;
854 // } // -D
856 // gOptGeomDUnits = utils::units::UnitFromString(dunits);
857
858 } // using root geom?
859#endif // #ifdef __CAN_USE_ROOT_GEOM__
860
861 // event file prefix
862 if( parser.OptionExists('o') ) {
863 LOG("gevgen_pghnl", pDEBUG) << "Reading the event filename prefix";
864 gOptEvFilePrefix = parser.ArgAsString('o');
865 } else {
866 LOG("gevgen_pghnl", pDEBUG)
867 << "Will set the default event filename prefix";
869 } //-o
870
871 // random number seed
872 if( parser.OptionExists("seed") ) {
873 LOG("gevgen_pghnl", pINFO) << "Reading random number seed";
874 gOptRanSeed = parser.ArgAsLong("seed");
875 } else {
876 LOG("gevgen_pghnl", pINFO) << "Unspecified random number seed - Using default";
877 gOptRanSeed = -1;
878 }
879
880 //
881 // >>> print the command line options
882 //
883
884 ostringstream gminfo;
885 if (gOptUsingRootGeom) {
886 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
887 << ", top volume: "
888 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
889 << ", length units: " << lunits;
890 // << ", density units: " << dunits;
891 }
892
893 LOG("gevgen_pghnl", pNOTICE)
894 << "\n\n"
895 << utils::print::PrintFramedMesg("gevgen_pghnl job configuration");
896
897 LOG("gevgen_pghnl", pNOTICE)
898 << "\n @@ Run number : " << gOptRunNu
899 << "\n @@ Random seed : " << gOptRanSeed
900 << "\n @@ Decay channel : " << utils::hnl::AsString(gOptDecayMode)
901 << "\n @@ Geometry : " << gminfo.str()
902 << "\n @@ Statistics : " << gOptNev << " events";
903}
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pFATAL
Definition Messenger.h:56
#define pWARN
Definition Messenger.h:60
Command line argument parser.
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
double GetHNLMass() const
long int gOptRanSeed
double gOptGeomLUnits
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
Long_t gOptRunNu
int gOptNev
string gOptRootGeom
string gOptEvFilePrefix
string gOptRootGeomTopVol
HNLDecayMode_t gOptDecayMode
string kDefOptSTune
void PrintSyntax(void)
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::ArgAsInt(), genie::CmdLnArgParser::ArgAsLong(), genie::CmdLnArgParser::ArgAsString(), genie::utils::hnl::AsString(), geom, genie::AlgFactory::GetAlgorithm(), genie::hnl::Decayer::GetHNLMass(), gOptDecayMode, gOptEvFilePrefix, gOptGeomLUnits, 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 643 of file gBeamHNLParticleGun.cxx.

644{
645 //string sname = "genie::EventGenerator";
646 //string sconfig = "BeamHNL";
648
649 LOG("gevgen_pghnl", pINFO)
650 << "Instantiating HNL generator.";
651
652 const Algorithm * algmcgen = algf->GetAlgorithm(kDefOptSName, kDefOptSConfig);
653 LOG("gevgen_pghnl", pDEBUG)
654 << "Got algorithm " << kDefOptSName.c_str() << "/" << kDefOptSConfig.c_str();;
655
656 const EventRecordVisitorI * mcgen =
657 dynamic_cast< const EventRecordVisitorI * >( algmcgen );
658 if(!mcgen) {
659 LOG("gevgen_pghnl", pFATAL) << "Couldn't instantiate the HNL generator";
660 gAbortingInErr = true;
661 exit(1);
662 }
663
664 LOG("gevgen_pghnl", pINFO)
665 << "HNL generator instantiated successfully.";
666
667 return mcgen;
668}
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 194 of file gBeamHNLParticleGun.cxx.

195{
196 // suppress ROOT Info messages
197 gErrorIgnoreLevel = kWarning;
198
199 // Parse command line arguments
200 GetCommandLineArgs(argc,argv);
201
202 // Init messenger and random number seed
203 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
205
206 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
207
208 // Get the HNL generator first to load config
209 // config loaded upon instantiation of HNLGenerator algorithm
210 // ==> Decayer::LoadConfig()
211 __attribute__((unused)) const EventRecordVisitorI * mcgen = HNLGenerator();
212 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
213 const Algorithm * algVtxGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::VertexGenerator", "Default");
214
215 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
216 const VertexGenerator * vtxGen = dynamic_cast< const VertexGenerator * >( algVtxGen );
217
218 bool geom_is_accessible = ! (gSystem->AccessPathName(gOptRootGeom.c_str()));
219 if (!geom_is_accessible) {
220 LOG("gevgen_pghnl", pFATAL)
221 << "The specified ROOT geometry doesn't exist! Initialization failed!";
222 exit(1);
223 }
224 vtxGen->SetGeomFile( gOptRootGeom );
225
226 if( !gOptRootGeoManager ) gOptRootGeoManager = TGeoManager::Import(gOptRootGeom.c_str());
227
228 TGeoVolume * top_volume = gOptRootGeoManager->GetTopVolume();
229 assert( top_volume );
230 __attribute__((unused)) TGeoShape * ts = top_volume->GetShape();
231
232 //TGeoBBox * box = (TGeoBBox *)ts;
233
234 string confString = kDefOptSName + kDefOptSConfig;
235
236 CoMLifetime = hnlgen->GetHNLLifetime();
237
238 gOptMassHNL = hnlgen->GetHNLMass();
239 std::vector< double > U4l2s = hnlgen->GetHNLCouplings();
240 gOptECoupling = U4l2s.at(0);
241 gOptMCoupling = U4l2s.at(1);
242 gOptTCoupling = U4l2s.at(2);
243 gOptIsMajorana = hnlgen->IsHNLMajorana();
244
245 std::string stIntChannels = hnlgen->GetHNLInterestingChannels(); int iChan = -1;
246 if( gOptIntChannels.size() > 0 ) gOptIntChannels.clear();
247 while( stIntChannels.size() > 0 ){ // read channels from right (lowest mass) to left (highest mass)
248 iChan++;
249 HNLDecayMode_t md = static_cast< HNLDecayMode_t >( iChan );
250 std::string tmpSt = stIntChannels.substr( stIntChannels.size()-1, stIntChannels.size() );
251 if( std::strcmp( tmpSt.c_str(), "1" ) == 0 )
252 gOptIntChannels.emplace_back( md );
253
254 stIntChannels.erase( stIntChannels.end()-1, stIntChannels.end() );
255 }
256
257 assert( gOptECoupling >= 0.0 && gOptMCoupling >= 0.0 && gOptTCoupling >= 0.0 );
258
259 // Initialize an Ntuple Writer to save GHEP records into a TTree
261 ntpw.CustomizeFilenamePrefix(gOptEvFilePrefix);
262 ntpw.Initialize();
263
264 LOG("gevgen_pghnl", pNOTICE)
265 << "Initialised Ntuple Writer";
266
267 // add another few branches to tree.
268 ntpw.EventTree()->Branch("hnl_mass", &gOptMassHNL, "gOptMassHNL/D");
269 ntpw.EventTree()->Branch("hnl_coup_e", &gOptECoupling, "gOptECoupling/D");
270 ntpw.EventTree()->Branch("hnl_coup_m", &gOptMCoupling, "gOptMCoupling/D");
271 ntpw.EventTree()->Branch("hnl_coup_t", &gOptTCoupling, "gOptTCoupling/D");
272 ntpw.EventTree()->Branch("hnl_ismaj", &gOptIsMajorana, "gOptIsMajorana/I");
273
274 // let's make HNL-specific FS branches until we get gntpc sorted out
275 ntpw.EventTree()->Branch("hnl_IS_E", &NTP_IS_E, "NTP_IS_E/D");
276 ntpw.EventTree()->Branch("hnl_IS_PX", &NTP_IS_PX, "NTP_IS_PX/D");
277 ntpw.EventTree()->Branch("hnl_IS_PY", &NTP_IS_PY, "NTP_IS_PY/D");
278 ntpw.EventTree()->Branch("hnl_IS_PZ", &NTP_IS_PZ, "NTP_IS_PZ/D");
279 ntpw.EventTree()->Branch("hnl_FS0_PDG", &NTP_FS0_PDG, "NTP_FS0_PDG/I");
280 ntpw.EventTree()->Branch("hnl_FS0_E", &NTP_FS0_E, "NTP_FS0_E/D");
281 ntpw.EventTree()->Branch("hnl_FS0_PX", &NTP_FS0_PX, "NTP_FS0_PX/D");
282 ntpw.EventTree()->Branch("hnl_FS0_PY", &NTP_FS0_PY, "NTP_FS0_PY/D");
283 ntpw.EventTree()->Branch("hnl_FS0_PZ", &NTP_FS0_PZ, "NTP_FS0_PZ/D");
284 ntpw.EventTree()->Branch("hnl_FS1_PDG", &NTP_FS1_PDG, "NTP_FS1_PDG/I");
285 ntpw.EventTree()->Branch("hnl_FS1_E", &NTP_FS1_E, "NTP_FS1_E/D");
286 ntpw.EventTree()->Branch("hnl_FS1_PX", &NTP_FS1_PX, "NTP_FS1_PX/D");
287 ntpw.EventTree()->Branch("hnl_FS1_PY", &NTP_FS1_PY, "NTP_FS1_PY/D");
288 ntpw.EventTree()->Branch("hnl_FS1_PZ", &NTP_FS1_PZ, "NTP_FS1_PZ/D");
289 ntpw.EventTree()->Branch("hnl_FS2_PDG", &NTP_FS2_PDG, "NTP_FS2_PDG/I");
290 ntpw.EventTree()->Branch("hnl_FS2_E", &NTP_FS2_E, "NTP_FS2_E/D");
291 ntpw.EventTree()->Branch("hnl_FS2_PX", &NTP_FS2_PX, "NTP_FS2_PX/D");
292 ntpw.EventTree()->Branch("hnl_FS2_PY", &NTP_FS2_PY, "NTP_FS2_PY/D");
293 ntpw.EventTree()->Branch("hnl_FS2_PZ", &NTP_FS2_PZ, "NTP_FS2_PZ/D");
294
295 // Create a MC job monitor for a periodically updated status file
296 GMCJMonitor mcjmonitor(gOptRunNu);
297 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
298
299 LOG("gevgen_pghnl", pNOTICE)
300 << "Initialised MC job monitor";
301
302 // Set GHEP print level
303 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
304
305#ifdef __CAN_USE_ROOT_GEOM__
306 // Read geometry bounding box - for vertex position generation
307 if( gOptUsingRootGeom ){
308 InitBoundingBox();
309 }
310#endif // #ifdef __CAN_USE_ROOT_GEOM__
311
312 // read in energy. This is always constant
313 //gOptEnergyHNL = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-Energy" );
314 gOptEnergyHNL = hnlgen->GetPGunEnergy();
315 if( gOptEnergyHNL <= gOptMassHNL ){
316 LOG( "gevgen_pghnl", pFATAL )
317 << "\nInsufficient energy " << gOptEnergyHNL << " GeV for HNL of mass " << gOptMassHNL << " GeV. Exiting."
318 << "\nPlease check ${GENIE}/config/CommonHNL.xml sections \"ParamSpace\" and \"ParticleGun\"";
319 exit(1);
320 }
321 assert( gOptEnergyHNL > gOptMassHNL );
322
323 // Event loop
324 int ievent = 0;
325
326 while (1)
327 {
328 if( gOptNev >= 10000 ){
329 if( ievent % (gOptNev / 1000 ) == 0 ){
330 int irat = ievent / ( gOptNev / 1000 );
331 std::cerr << 0.1 * irat << " % " << " ( " << ievent
332 << " / " << gOptNev << " ) \r" << std::flush;
333 }
334 } else if( gOptNev >= 100 ){
335 if( ievent % (gOptNev / 10 ) == 0 ){
336 int irat = ievent / ( gOptNev / 10 );
337 std::cerr << 10.0 * irat << " % " << " ( " << ievent
338 << " / " << gOptNev << " ) \r" << std::flush;
339 }
340 }
341
342 if((ievent) == gOptNev) break;
343
344 LOG("gevgen_pghnl", pNOTICE)
345 << " *** Generating event............ " << ievent;
346
347 int hpdg = genie::kPdgHNL;
348 EventRecord * event = new EventRecord;
349
350 event->SetProbability( CoMLifetime );
351
352 int decay = (int) gOptDecayMode;
353
354 SimpleHNL sh( "HNL", ievent, hpdg, genie::kPdgKP,
356
357 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
358 CoMLifetime = sh.GetCoMLifetime();
359
360 if( gOptDecayMode == kHNLDcyNull ){ // select from available modes
361 LOG("gevgen_pghnl", pNOTICE)
362 << "No decay mode specified - sampling from all available modes.";
363
364 std::vector< HNLDecayMode_t > * intChannels = &gOptIntChannels;
365
366 decay = SelectDecayMode( intChannels, sh );
367 }
368
370 event->AttachSummary(interaction);
371
372 // generate origin position and momentum direction for this HNL
373 TLorentzVector * x4orig = GenerateOriginPosition( event );
374 TLorentzVector * p4HNL = GenerateOriginMomentum( event );
375
376 // and add as Particle(0)
377 GHepParticle ptHNL( kPdgHNL, kIStInitialState, -1, -1, -1, -1, *p4HNL, *x4orig );
378 event->AddParticle( ptHNL );
379
380 LOG("gevgen_pghnl", pDEBUG)
381 << "Note decay mode is " << utils::hnl::AsString(gOptDecayMode);
382
383 // Simulate decay
384 hnlgen->ProcessEventRecord(event);
385 vtxGen->ProcessEventRecord(event);
386
387 // add the FS 4-momenta to special branches
388 // Quite inelegant. Gets the job done, though
389 NTP_FS0_PDG = (event->Particle(1))->Pdg();
390 NTP_FS0_E = ((event->Particle(1))->P4())->E();
391 NTP_FS0_PX = ((event->Particle(1))->P4())->Px();
392 NTP_FS0_PY = ((event->Particle(1))->P4())->Py();
393 NTP_FS0_PZ = ((event->Particle(1))->P4())->Pz();
394 NTP_FS1_PDG = (event->Particle(2))->Pdg();
395 NTP_FS1_E = ((event->Particle(2))->P4())->E();
396 NTP_FS1_PX = ((event->Particle(2))->P4())->Px();
397 NTP_FS1_PY = ((event->Particle(2))->P4())->Py();
398 NTP_FS1_PZ = ((event->Particle(2))->P4())->Pz();
399 if( event->Particle(3) ){
400 NTP_FS2_PDG = (event->Particle(3))->Pdg();
401 NTP_FS2_E = ((event->Particle(3))->P4())->E();
402 NTP_FS2_PX = ((event->Particle(3))->P4())->Px();
403 NTP_FS2_PY = ((event->Particle(3))->P4())->Py();
404 NTP_FS2_PZ = ((event->Particle(3))->P4())->Pz();
405 }
406 else{
407 NTP_FS2_PDG = 0;
408 NTP_FS2_E = 0.0;
409 NTP_FS2_PX = 0.0;
410 NTP_FS2_PY = 0.0;
411 NTP_FS2_PZ = 0.0;
412 }
413
414 // Generate (or read) a position for the decay vertex
415 // also currently handles the event weight
416 TLorentzVector x4mm = GeneratePosition( event );
417 event->SetWeight( evWeight );
418
419 // why does InitState show the wrong p4 here?
420 interaction->InitStatePtr()->SetProbeP4( *(event->Particle(0)->P4()) );
421
422 LOG("gevgen_pghnl", pINFO)
423 << "Generated event: " << *event;
424
425 // Add event at the output ntuple, refresh the mc job monitor & clean-up
426 ntpw.AddEventRecord(ievent, event);
427 mcjmonitor.Update(ievent,event);
428 delete event;
429
430 ievent++;
431 } // event loop
432
433 // Save the generated event tree & close the output file
434 ntpw.Save();
435
436 LOG("gevgen_pghnl", pNOTICE) << "Done!";
437
438 return 0;
439}
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.
static void SetPrintLevel(int print_level)
virtual void SetProbability(double prob)
Definition GHepRecord.h:131
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
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
std::string GetHNLInterestingChannels() const
void ProcessEventRecord(GHepRecord *event) const
double GetHNLLifetime() const
bool IsHNLMajorana() const
double GetPGunEnergy() const
std::vector< double > GetHNLCouplings() const
void SetGeomFile(std::string geomfile) const
NtpMCFormat_t kDefOptNtpFormat
double gOptECoupling
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 NTP_FS0_PX
double evWeight
double NTP_IS_PX
double NTP_FS0_E
double NTP_FS2_PZ
double NTP_FS2_PY
int NTP_FS1_PDG
int NTP_FS2_PDG
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
TLorentzVector GeneratePosition(GHepRecord *event)
TLorentzVector * GenerateOriginPosition(GHepRecord *event)
const EventRecordVisitorI * HNLGenerator(void)
int SelectDecayMode(std::vector< HNLDecayMode_t > *intChannels, SimpleHNL sh)
void GetCommandLineArgs(int argc, char **argv)
TLorentzVector * GenerateOriginMomentum(GHepRecord *event)
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
@ kIStInitialState
Definition GHepStatus.h:29
const int kPdgKP
Definition PDGCodes.h:172
const int kPdgHNL
Definition PDGCodes.h:224

References genie::NtpWriter::AddEventRecord(), genie::utils::hnl::AsString(), CoMLifetime, genie::NtpWriter::CustomizeFilenamePrefix(), genie::NtpWriter::EventTree(), evWeight, GenerateOriginMomentum(), GenerateOriginPosition(), GeneratePosition(), genie::AlgFactory::GetAlgorithm(), genie::hnl::SimpleHNL::GetCoMLifetime(), GetCommandLineArgs(), genie::hnl::Decayer::GetHNLCouplings(), genie::hnl::Decayer::GetHNLInterestingChannels(), genie::hnl::Decayer::GetHNLLifetime(), genie::hnl::Decayer::GetHNLMass(), genie::hnl::Decayer::GetPGunEnergy(), genie::hnl::SimpleHNL::GetValidChannels(), gOptDecayMode, gOptECoupling, gOptEnergyHNL, gOptEvFilePrefix, gOptIntChannels, gOptIsMajorana, 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, kDefOptSName, 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, pDEBUG, pFATAL, pINFO, pNOTICE, genie::hnl::Decayer::ProcessEventRecord(), genie::hnl::VertexGenerator::ProcessEventRecord(), genie::utils::app_init::RandGen(), genie::NtpWriter::Save(), SelectDecayMode(), genie::hnl::VertexGenerator::SetGeomFile(), genie::GHepRecord::SetPrintLevel(), genie::GHepRecord::SetProbability(), genie::InitialState::SetProbeP4(), genie::GMCJMonitor::SetRefreshRate(), and genie::GMCJMonitor::Update().

◆ PrintSyntax()

void PrintSyntax ( void )

Definition at line 905 of file gBeamHNLParticleGun.cxx.

906{
907 LOG("gevgen_pghnl", pFATAL)
908 << "\n **Syntax**"
909 << "\n gevgen_pghnl [-h] "
910 << "\n [-r run#]"
911 << "\n -n n_of_events"
912 << "\n [-m decay_mode]"
913 << "\n [-g geometry (ROOT file)]"
914 << "\n [-L length_units_at_geom]"
915 << "\n [-o output_event_file_prefix]"
916 << "\n [--seed random_number_seed]"
917 << "\n [--message-thresholds xml_file]"
918 << "\n [--event-record-print-level level]"
919 << "\n [--mc-job-status-refresh-rate rate]"
920 << "\n"
921 << " Please also read the detailed documentation at http://www.genie-mc.org"
922 << " or look at the source code: $GENIE/src/Apps/gBeamHNLParticleGun.cxx"
923 << "\n";
924}

References LOG, and pFATAL.

Referenced by GetCommandLineArgs().

◆ SelectDecayMode()

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

Definition at line 670 of file gBeamHNLParticleGun.cxx.

671{
672 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
673
674 // set CoM lifetime now if unset
675 if( CoMLifetime < 0.0 ){
677 LOG( "gevgen_pghnl", pDEBUG )
678 << "Rest frame CoMLifetime = " << CoMLifetime << " [GeV^{-1}]";
679 }
680
681 for( std::vector< HNLDecayMode_t >::iterator it = intChannels->begin(); it != intChannels->end(); ++it ){
682 HNLDecayMode_t mode = *it;
683 auto mapG = gammaMap.find( mode );
684 double theGamma = mapG->second;
685 LOG("gevgen_pghnl", pDEBUG)
686 << "For mode " << utils::hnl::AsString( mode ) << " gamma = " << theGamma;
687 }
688
689 std::map< HNLDecayMode_t, double > intMap =
690 selector::SetInterestingChannels( (*intChannels), gammaMap );
691
692 sh.SetInterestingChannels( intMap );
693
694 // get probability that channels in intChannels will be selected
695 std::map< HNLDecayMode_t, double > PMap =
697
698 // now do a random throw
700 double ranThrow = rnd->RndGen().Uniform( 0., 1. ); // HNL's fate is sealed.
701
702 HNLDecayMode_t selectedDecayChan =
703 selector::SelectChannelInclusive( PMap, ranThrow );
704
705 int decay = ( int ) selectedDecayChan;
706 return decay;
707}
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, genie::hnl::SimpleHNL::GetCoMLifetime(), genie::hnl::selector::GetProbabilities(), genie::hnl::SimpleHNL::GetValidChannels(), genie::RandomGen::Instance(), LOG, pDEBUG, genie::RandomGen::RndGen(), genie::hnl::selector::SelectChannelInclusive(), genie::hnl::selector::SetInterestingChannels(), and genie::hnl::SimpleHNL::SetInterestingChannels().

Referenced by main().

Variable Documentation

◆ CoMLifetime

double CoMLifetime = -1.0

Definition at line 187 of file gBeamHNLParticleGun.cxx.

◆ evProdVtx

double evProdVtx[3] = {0.0, 0.0, 0.0}

Definition at line 189 of file gBeamHNLParticleGun.cxx.

189{0.0, 0.0, 0.0}; // mm

◆ evWeight

double evWeight = 1.0

Definition at line 191 of file gBeamHNLParticleGun.cxx.

◆ fdx

double fdx = 0

Definition at line 172 of file gBeamHNLParticleGun.cxx.

◆ fdy

double fdy = 0

Definition at line 173 of file gBeamHNLParticleGun.cxx.

◆ fdz

double fdz = 0

Definition at line 174 of file gBeamHNLParticleGun.cxx.

◆ fox

double fox = 0

Definition at line 175 of file gBeamHNLParticleGun.cxx.

◆ foy

double foy = 0

Definition at line 176 of file gBeamHNLParticleGun.cxx.

◆ foz

double foz = 0

Definition at line 177 of file gBeamHNLParticleGun.cxx.

◆ gOptDecayMode

HNLDecayMode_t gOptDecayMode = kHNLDcyNull

Definition at line 155 of file gBeamHNLParticleGun.cxx.

◆ gOptECoupling

double gOptECoupling = -1

Definition at line 148 of file gBeamHNLParticleGun.cxx.

◆ gOptEnergyHNL

double gOptEnergyHNL = -1

Definition at line 146 of file gBeamHNLParticleGun.cxx.

◆ gOptEvFilePrefix

string gOptEvFilePrefix = kDefOptEvFilePrefix

Definition at line 158 of file gBeamHNLParticleGun.cxx.

◆ gOptGeomLUnits

double gOptGeomLUnits = 0

Definition at line 168 of file gBeamHNLParticleGun.cxx.

◆ gOptHNLKind

int gOptHNLKind = -1

Definition at line 153 of file gBeamHNLParticleGun.cxx.

◆ gOptIntChannels

std::vector< HNLDecayMode_t > gOptIntChannels

Definition at line 156 of file gBeamHNLParticleGun.cxx.

◆ gOptIsMajorana

bool gOptIsMajorana = false

Definition at line 152 of file gBeamHNLParticleGun.cxx.

◆ gOptMassHNL

double gOptMassHNL = -1

Definition at line 147 of file gBeamHNLParticleGun.cxx.

◆ gOptMCoupling

double gOptMCoupling = -1

Definition at line 149 of file gBeamHNLParticleGun.cxx.

◆ gOptNev

int gOptNev = 10

Definition at line 144 of file gBeamHNLParticleGun.cxx.

◆ gOptRanSeed

long int gOptRanSeed = -1

Definition at line 169 of file gBeamHNLParticleGun.cxx.

◆ gOptRootGeom

string gOptRootGeom

Definition at line 160 of file gBeamHNLParticleGun.cxx.

◆ gOptRootGeomTopVol

string gOptRootGeomTopVol = ""

Definition at line 167 of file gBeamHNLParticleGun.cxx.

◆ gOptRunNu

Long_t gOptRunNu = 1000

Definition at line 143 of file gBeamHNLParticleGun.cxx.

◆ gOptTCoupling

double gOptTCoupling = -1

Definition at line 150 of file gBeamHNLParticleGun.cxx.

◆ gOptUsingRootGeom

bool gOptUsingRootGeom = false

Definition at line 159 of file gBeamHNLParticleGun.cxx.

◆ kDefOptEvFilePrefix

string kDefOptEvFilePrefix = "gntp"

Definition at line 136 of file gBeamHNLParticleGun.cxx.

◆ kDefOptGeomDUnits

string kDefOptGeomDUnits = "g_cm3"

Definition at line 134 of file gBeamHNLParticleGun.cxx.

◆ kDefOptGeomLUnits

string kDefOptGeomLUnits = "mm"

Definition at line 133 of file gBeamHNLParticleGun.cxx.

◆ kDefOptNtpFormat

NtpMCFormat_t kDefOptNtpFormat = kNFGHEP

Definition at line 135 of file gBeamHNLParticleGun.cxx.

◆ kDefOptSConfig

string kDefOptSConfig = "BeamHNL"

Definition at line 139 of file gBeamHNLParticleGun.cxx.

◆ kDefOptSName

string kDefOptSName = "genie::EventGenerator"

Definition at line 138 of file gBeamHNLParticleGun.cxx.

◆ kDefOptSTune

string kDefOptSTune = "GHNL20_00a_00_000"

Definition at line 140 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS0_E

double NTP_FS0_E = 0.

Definition at line 180 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS0_PDG

int NTP_FS0_PDG = 0

Definition at line 183 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS0_PX

double NTP_FS0_PX = 0.

Definition at line 180 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS0_PY

double NTP_FS0_PY = 0.

Definition at line 180 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS0_PZ

double NTP_FS0_PZ = 0.

Definition at line 180 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS1_E

double NTP_FS1_E = 0.

Definition at line 181 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS1_PDG

int NTP_FS1_PDG = 0

Definition at line 183 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS1_PX

double NTP_FS1_PX = 0.

Definition at line 181 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS1_PY

double NTP_FS1_PY = 0.

Definition at line 181 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS1_PZ

double NTP_FS1_PZ = 0.

Definition at line 181 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS2_E

double NTP_FS2_E = 0.

Definition at line 182 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS2_PDG

int NTP_FS2_PDG = 0

Definition at line 183 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS2_PX

double NTP_FS2_PX = 0.

Definition at line 182 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS2_PY

double NTP_FS2_PY = 0.

Definition at line 182 of file gBeamHNLParticleGun.cxx.

◆ NTP_FS2_PZ

double NTP_FS2_PZ = 0.

Definition at line 182 of file gBeamHNLParticleGun.cxx.

◆ NTP_IS_E

double NTP_IS_E = 0.

Definition at line 179 of file gBeamHNLParticleGun.cxx.

◆ NTP_IS_PX

double NTP_IS_PX = 0.

Definition at line 179 of file gBeamHNLParticleGun.cxx.

◆ NTP_IS_PY

double NTP_IS_PY = 0.

Definition at line 179 of file gBeamHNLParticleGun.cxx.

◆ NTP_IS_PZ

double NTP_IS_PZ = 0.

Definition at line 179 of file gBeamHNLParticleGun.cxx.