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

Go to the source code of this file.

Typedefs

typedef enum t_HNLValidation HNLValidation_t

Enumerations

enum  t_HNLValidation {
  kValNone = 0 , kValFluxFromDk2nu = 1 , kValDecay = 2 , kValGeom = 3 ,
  kValFullSim = 4
}

Functions

void GetCommandLineArgs (int argc, char **argv)
void ReadInConfig (void)
void PrintSyntax (void)
const EventRecordVisitorIHNLGenerator (void)
int TestDecay (void)
int main (int argc, char **argv)

Variables

string kDefOptGeomLUnits = "mm"
string kDefOptGeomAUnits = "rad"
string kDefOptGeomTUnits = "ns"
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
HNLValidation_t gOptValidationMode = kValNone
double gCfgMassHNL = -1
double gCfgECoupling = -1
double gCfgMCoupling = -1
double gCfgTCoupling = -1
bool gCfgIsMajorana = false
int gCfgHNLKind = 2
double CoMLifetime = -1
HNLProd_t gCfgProdMode = kHNLProdNull
HNLDecayMode_t gCfgDecayMode = kHNLDcyNull
std::vector< HNLDecayMode_tgCfgIntChannels
double gCfgUserOx = -1
double gCfgUserOy = -1
double gCfgUserOz = -1
double gCfgUserAx1 = -1
double gCfgUserAz = -1
double gCfgUserAx2 = -1
double gCfgBoxLx = -1
double gCfgBoxLy = -1
double gCfgBoxLz = -1
double gCfgParentEnergy = -1
double gCfgParentTheta = -1
double gCfgParentPhi = -1
double gCfgParentOx = -1
double gCfgParentOy = -1
double gCfgParentOz = -1
double gCfgHNLCx = -1
double gCfgHNLCy = -1
double gCfgHNLCz = -1
double gOptEnergyHNL = -1
bool gOptIsMonoEnFlux = true
string gOptEvFilePrefix = kDefOptEvFilePrefix
bool gOptUsingRootGeom = false
string gOptRootGeom
string geom
string lunits
string aunits
string tunits
double gOptGeomLUnits = 0
double gOptGeomAUnits = 0
double gOptGeomTUnits = 0
string gOptRootGeomTopVol = ""
double fdx = 0
double fdy = 0
double fdz = 0
double fox = 0
double foy = 0
double foz = 0
long int gOptRanSeed = -1

Typedef Documentation

◆ HNLValidation_t

Enumeration Type Documentation

◆ t_HNLValidation

Enumerator
kValNone 
kValFluxFromDk2nu 
kValDecay 
kValGeom 
kValFullSim 

Definition at line 134 of file gBeamHNLValidationApp.cxx.

134 {
135
136 kValNone = 0,
138 kValDecay = 2,
139 kValGeom = 3,
140 kValFullSim = 4
141
enum t_HNLValidation HNLValidation_t

Function Documentation

◆ GetCommandLineArgs()

void GetCommandLineArgs ( int argc,
char ** argv )

Definition at line 1214 of file gBeamHNLValidationApp.cxx.

1215{
1216 LOG("gevald_hnl", pINFO) << "Parsing command line arguments";
1217
1218 // Parse run options for this app
1219
1220 CmdLnArgParser parser(argc,argv);
1221
1222 // help?
1223 bool help = parser.OptionExists('h');
1224 if(help) {
1225 PrintSyntax();
1226 exit(0);
1227 }
1228
1229 // force the app to look at HNL tune by default
1230 // if user passes --tune argument, look at the user input tune instead
1231 char * expargv[ argc + 2 ];
1232 bool didFindUserInputTune = false;
1233 std::string stExtraTuneBit = kDefOptSTune;
1234
1235 if( parser.OptionExists("tune") ){
1236 didFindUserInputTune = true;
1237 stExtraTuneBit = parser.ArgAsString("tune");
1238 LOG( "gevgen_hnl", pWARN )
1239 << "Using input HNL tune " << parser.ArgAsString("tune");
1240 } else {
1241 LOG( "gevgen_hnl", pWARN )
1242 << "Using default HNL tune " << kDefOptSTune;
1243 }
1244 // append this to argv
1245 for( int iArg = 0; iArg < argc; iArg++ ){
1246 expargv[iArg] = argv[iArg];
1247 }
1248 if( !didFindUserInputTune ){
1249 char * chBit = const_cast< char * >( stExtraTuneBit.c_str() ); // ugh. Ugly.
1250 std::string stune("--tune"); char * tBit = const_cast< char * >( stune.c_str() );
1251 expargv[argc] = tBit;
1252 expargv[argc+1] = chBit;
1253 }
1254
1255 // Common run options.
1256 int expargc = ( didFindUserInputTune ) ? argc : argc+2;
1257 std::string stnull(""); char * nBit = const_cast< char * >( stnull.c_str() );
1258 expargv[expargc] = nBit;
1259
1260 RunOpt::Instance()->ReadFromCommandLine(expargc,expargv);
1261
1262 // run number
1263 if( parser.OptionExists('r') ) {
1264 LOG("gevald_hnl", pDEBUG) << "Reading MC run number";
1265 gOptRunNu = parser.ArgAsLong('r');
1266 } else {
1267 LOG("gevald_hnl", pDEBUG) << "Unspecified run number - Using default";
1268 gOptRunNu = 1000;
1269 } //-r
1270
1271 // number of events
1272 if( parser.OptionExists('n') ) {
1273 LOG("gevald_hnl", pDEBUG)
1274 << "Reading number of events to generate";
1275 gOptNev = parser.ArgAsInt('n');
1276 } else {
1277 LOG("gevald_hnl", pFATAL)
1278 << "You need to specify the number of events";
1279 PrintSyntax();
1280 exit(0);
1281 } //-n
1282
1283 if( parser.OptionExists('M') ) {
1284 LOG("gevald_hnl", pDEBUG)
1285 << "Detecting mode. . .";
1286 gOptValidationMode = (HNLValidation_t) parser.ArgAsInt('M');
1287 } else {
1288 LOG("gevald_hnl", pFATAL)
1289 << "You must specify a validation mode.";
1290 PrintSyntax();
1291 exit(0);
1292 } // -M
1293
1294 bool isMonoEnergeticFlux = true;
1295#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
1296 if( parser.OptionExists('f') ) {
1297 LOG("gevald_hnl", pDEBUG)
1298 << "A flux has been offered. Searching this path: " << parser.ArgAsString('f');
1299 isMonoEnergeticFlux = false;
1300 gOptFluxFilePath = parser.ArgAsString('f');
1301
1302 // check if this is dk2nu
1303 //if( gOptFluxFilePath.find( "dk2nu" ) != string::npos ){
1304 if( gSystem->OpenDirectory( gOptFluxFilePath.c_str() ) != NULL ){
1305 gOptIsUsingDk2nu = true;
1306 LOG("gevald_hnl", pDEBUG)
1307 << "dk2nu flux files detected. Will create flux spectrum dynamically.";
1308 } else {
1309 LOG("gevald_hnl", pFATAL)
1310 << "Invalid flux file path " << gOptFluxFilePath;
1311 exit(1);
1312 }
1313 } else {
1314 // we need the 'E' option! Log it and pass below
1315 LOG("gevald_hnl", pINFO)
1316 << "No flux file offered. Assuming monoenergetic flux.";
1317 } //-f
1318 gOptIsMonoEnFlux = isMonoEnergeticFlux;
1319#endif // #ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
1320
1321#ifdef __CAN_USE_ROOT_GEOM__
1322 if( parser.OptionExists('g') ) {
1323 LOG("gevald_hnl", pDEBUG) << "Getting input geometry";
1324 geom = parser.ArgAsString('g');
1325
1326 // is it a ROOT file that contains a ROOT geometry?
1327 bool accessible_geom_file =
1328 ! (gSystem->AccessPathName(geom.c_str()));
1329 if (accessible_geom_file) {
1331 gOptUsingRootGeom = true;
1332 } else {
1333 LOG("gevald_hnl", pFATAL)
1334 << "Geometry option is not a ROOT file. Please use ROOT geom.";
1335 PrintSyntax();
1336 exit(1);
1337 }
1338 } else if( gOptValidationMode == kValGeom ) {
1339
1340 LOG("gevald_hnl", pFATAL)
1341 << "No geometry option specified - Exiting";
1342 PrintSyntax();
1343 exit(1);
1344 } //-g
1345
1346 if( parser.OptionExists('L') ) {
1347 lunits = parser.ArgAsString('L');
1348 LOG("gevald_hnl", pDEBUG) << "Setting length units to " << lunits.c_str();
1349 } else {
1350 LOG("gevald_hnl", pDEBUG) << "Using default geometry length units";
1352 } // -L
1354
1355 if( parser.OptionExists('A') ) {
1356 aunits = parser.ArgAsString('A');
1357 LOG("gevald_hnl", pDEBUG) << "Setting angle units to " << aunits.c_str();
1358 } else {
1359 LOG("gevald_hnl", pDEBUG) << "Using default angle length units";
1361 } // -A
1363
1364 if( parser.OptionExists('T') ) {
1365 tunits = parser.ArgAsString('T');
1366 LOG("gevald_hnl", pDEBUG) << "Setting time units to " << tunits.c_str();
1367 } else {
1368 LOG("gevald_hnl", pDEBUG) << "Using default geometry time units";
1370 } // -T
1372
1373#endif // #ifdef __CAN_USE_ROOT_GEOM__
1374
1375 // event file prefix
1376 if( parser.OptionExists('o') ) {
1377 LOG("gevald_hnl", pDEBUG) << "Reading the event filename prefix";
1378 gOptEvFilePrefix = parser.ArgAsString('o');
1379 } else {
1380 LOG("gevald_hnl", pDEBUG)
1381 << "Will set the default event filename prefix";
1383 } //-o
1384
1385 // random number seed
1386 if( parser.OptionExists("seed") ) {
1387 LOG("gevald_hnl", pINFO) << "Reading random number seed";
1388 gOptRanSeed = parser.ArgAsLong("seed");
1389 } else {
1390 LOG("gevald_hnl", pINFO) << "Unspecified random number seed - Using default";
1391 gOptRanSeed = -1;
1392 }
1393
1394 //
1395 // >>> print the command line options
1396 //
1397
1398 ostringstream gminfo;
1399 if (gOptUsingRootGeom) {
1400 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
1401 << ", top volume: "
1402 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
1403 << ", length units: " << lunits;
1404 // << ", density units: " << dunits;
1405 }
1406
1407 LOG("gevald_hnl", pNOTICE)
1408 << "\n\n"
1409 << utils::print::PrintFramedMesg("gevald_hnl job configuration");
1410
1411 LOG("gevald_hnl", pNOTICE)
1412 << "\n @@ Run number : " << gOptRunNu
1413 << "\n @@ Random seed : " << gOptRanSeed
1414 << "\n @@ HNL mass : " << gCfgMassHNL << " GeV"
1415 << "\n @@ Decay channel : " << utils::hnl::AsString(gCfgDecayMode)
1416 << "\n @@ Flux path : " << gOptFluxFilePath
1417 << "\n @@ Geometry : " << gminfo.str()
1418 << "\n @@ Statistics : " << gOptNev << " events";
1419}
#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
Command line argument parser.
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
long int gOptRanSeed
double gOptGeomLUnits
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
bool gOptUsingRootGeom
Long_t gOptRunNu
int gOptNev
string gOptRootGeom
string gOptEvFilePrefix
string gOptRootGeomTopVol
bool gOptIsMonoEnFlux
string kDefOptSTune
double gOptGeomTUnits
HNLValidation_t gOptValidationMode
HNLDecayMode_t gCfgDecayMode
string lunits
string aunits
string geom
void PrintSyntax(void)
string kDefOptGeomTUnits
double gOptGeomAUnits
string kDefOptGeomAUnits
double gCfgMassHNL
string tunits
string AsString(genie::hnl::HNLDecayMode_t hnldm)
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(), aunits, gCfgDecayMode, gCfgMassHNL, geom, gOptEvFilePrefix, gOptGeomAUnits, gOptGeomLUnits, gOptGeomTUnits, gOptIsMonoEnFlux, gOptNev, gOptRanSeed, gOptRootGeom, gOptRootGeomTopVol, gOptRunNu, gOptUsingRootGeom, gOptValidationMode, genie::RunOpt::Instance(), kDefOptEvFilePrefix, kDefOptGeomAUnits, kDefOptGeomLUnits, kDefOptGeomTUnits, kDefOptSTune, kValGeom, LOG, lunits, genie::CmdLnArgParser::OptionExists(), pDEBUG, pFATAL, pINFO, pNOTICE, genie::utils::print::PrintFramedMesg(), PrintSyntax(), pWARN, genie::RunOpt::ReadFromCommandLine(), tunits, and genie::utils::units::UnitFromString().

Referenced by main().

◆ HNLGenerator()

const EventRecordVisitorI * HNLGenerator ( void )

Definition at line 1187 of file gBeamHNLValidationApp.cxx.

1188{
1189 //string sname = "genie::EventGenerator";
1190 //string sconfig = "BeamHNL";
1192
1193 LOG("gevald_hnl", pINFO)
1194 << "Instantiating HNL generator.";
1195
1196 const Algorithm * algmcgen = algf->GetAlgorithm(kDefOptSName, kDefOptSConfig);
1197 LOG("gevald_hnl", pDEBUG)
1198 << "Got algorithm " << kDefOptSName.c_str() << "/" << kDefOptSConfig.c_str();;
1199
1200 const EventRecordVisitorI * mcgen =
1201 dynamic_cast< const EventRecordVisitorI * >( algmcgen );
1202 if(!mcgen) {
1203 LOG("gevald_hnl", pFATAL) << "Couldn't instantiate the HNL generator";
1204 gAbortingInErr = true;
1205 exit(1);
1206 }
1207
1208 LOG("gevald_hnl", pINFO)
1209 << "HNL generator instantiated successfully.";
1210
1211 return mcgen;
1212}
The GENIE Algorithm Factory.
Definition AlgFactory.h:39
const Algorithm * GetAlgorithm(const AlgId &algid)
static AlgFactory * Instance()
Algorithm abstract base class.
Definition Algorithm.h:54
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 TestDecay().

◆ main()

int main ( int argc,
char ** argv )

Definition at line 255 of file gBeamHNLValidationApp.cxx.

256{
257 // suppress ROOT Info messages
258 gErrorIgnoreLevel = kWarning;
259
260 // Parse command line arguments
261 GetCommandLineArgs(argc,argv);
262
263 // Get the validation configuration
264 ReadInConfig();
265
266 // Init messenger and random number seed
267 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
269
270 switch( gOptValidationMode ){
271
272 case kValFluxFromDk2nu: return TestFluxFromDk2nu(); break;
273 case kValDecay: return TestDecay(); break;
274 case kValGeom: return TestGeom(); break;
275 default: LOG( "gevald_hnl", pFATAL ) << "I didn't recognise this mode. Goodbye world!"; break;
276 }
277
278 return 0;
279}
void GetCommandLineArgs(int argc, char **argv)
int TestDecay(void)
void ReadInConfig(void)
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99

References GetCommandLineArgs(), gOptRanSeed, gOptValidationMode, genie::RunOpt::Instance(), kValDecay, kValFluxFromDk2nu, kValGeom, LOG, genie::utils::app_init::MesgThresholds(), pFATAL, genie::utils::app_init::RandGen(), ReadInConfig(), and TestDecay().

◆ PrintSyntax()

void PrintSyntax ( void )

Definition at line 1509 of file gBeamHNLValidationApp.cxx.

1510{
1511 LOG("gevald_hnl", pFATAL)
1512 << "\n **Syntax**"
1513 << "\n gevald_hnl [-h] "
1514 << "\n [-r run#]"
1515 << "\n -n n_of_events"
1516 << "\n [-f path_to_flux_files]"
1517 << "\n [-g geometry_file]"
1518 << "\n -M mode:"
1519 << "\n 1: Flux prediction from dk2nu files. Needs -f option"
1520 << "\n 2: HNL decay validation. Specify an origin point and 4-momentum"
1521 << "\n in the \"ParticleGun\" section in config. Needs -g option."
1522 << "\n 3: Custom geometry file validation. Needs -g option"
1523 << "\n Specify origin, momentum, and wiggle room for both of these in the"
1524 << "\n \"ParticleGun\" section in config"
1525 << "\n Regardless of how many events you ask for, this will evaluate 125x81"
1526 << "\n events: 5^3 from wiggling origin and 9^2 from wiggling momentum direction"
1527 << "\n 4: Full simulation (like gevgen_hnl but with lots of debug!)"
1528 << "\n"
1529 << "\n The configuration file lives at $GENIE/config/CommonHNL.xml - see"
1530 << " <param_set name=\"Validation\">"
1531 << "\n"
1532 << "\n Please also read the detailed documentation at http://www.genie-mc.org"
1533 << "\n or look at the source code: $GENIE/src/Apps/gBeamHNLValidationApp.cxx"
1534 << "\n";
1535}

References LOG, and pFATAL.

Referenced by GetCommandLineArgs().

◆ ReadInConfig()

void ReadInConfig ( void )

Definition at line 1421 of file gBeamHNLValidationApp.cxx.

1422{
1423 LOG("gevald_hnl", pFATAL)
1424 << "Reading in validation configuration. . .";
1425
1426 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
1427 __attribute__((unused)) const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
1428
1429 CoMLifetime = hnlgen->GetHNLLifetime();
1430
1431 gCfgMassHNL = hnlgen->GetHNLMass();
1432 std::vector<double> U4l2s = hnlgen->GetHNLCouplings();
1433 gCfgECoupling = U4l2s.at(0);
1434 gCfgMCoupling = U4l2s.at(1);
1435 gCfgTCoupling = U4l2s.at(2);
1436 gCfgIsMajorana = hnlgen->IsHNLMajorana();
1437
1438 //gOptEnergyHNL = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-Energy" );
1439 gOptEnergyHNL = hnlgen->GetPGunEnergy();
1440
1441 std::vector< double > PGDirection = hnlgen->GetPGunDirection();
1442 gCfgHNLCx = PGDirection.at(0), gCfgHNLCy = PGDirection.at(1), gCfgHNLCz = PGDirection.at(2);
1443 /*
1444 gCfgHNLCx = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cx" );
1445 gCfgHNLCy = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cy" );
1446 gCfgHNLCz = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-cz" );
1447 */
1448
1449 double dircosMag2 = std::pow( gCfgHNLCx, 2.0 ) +
1450 std::pow( gCfgHNLCy, 2.0 ) +
1451 std::pow( gCfgHNLCz, 2.0 );
1452 double invdircosmag = 1.0 / std::sqrt( dircosMag2 );
1453 gCfgHNLCx *= invdircosmag;
1454 gCfgHNLCy *= invdircosmag;
1455 gCfgHNLCz *= invdircosmag;
1456
1457 std::string stIntChannels = hnlgen->GetHNLInterestingChannels(); int iChan = -1;
1458 if( gCfgIntChannels.size() > 0 ) gCfgIntChannels.clear();
1459 while( stIntChannels.size() > 0 ){ // read channels from right (lowest mass) to left (highest mass)
1460 iChan++;
1461 HNLDecayMode_t md = static_cast< HNLDecayMode_t >( iChan );
1462 std::string tmpSt = stIntChannels.substr( stIntChannels.size()-1, stIntChannels.size() );
1463 if( std::strcmp( tmpSt.c_str(), "1" ) == 0 )
1464 gCfgIntChannels.emplace_back( md );
1465
1466 stIntChannels.erase( stIntChannels.end()-1, stIntChannels.end() );
1467 }
1468
1469 const Algorithm * algFluxCreator = AlgFactory::Instance()->GetAlgorithm("genie::hnl::FluxCreator", "Default");
1470 __attribute__((unused)) const FluxCreator * fluxCreator = dynamic_cast< const FluxCreator * >( algFluxCreator );
1471
1472 std::vector< double > b2uTranslation = fluxCreator->GetB2UTranslation();
1473 gCfgUserOx = b2uTranslation.at(0); gCfgUserOy = b2uTranslation.at(1); gCfgUserOz = b2uTranslation.at(2);
1474 std::vector< double > b2uRotation = fluxCreator->GetB2URotation();
1475 gCfgUserAx1 = b2uRotation.at(0); gCfgUserAz = b2uRotation.at(1); gCfgUserAx2 = b2uRotation.at(2);
1476
1477 // now transform the lengths and angles to the correct units
1481
1485
1486 ostringstream csts;
1487 csts << "Read out the following config:"
1488 << "\n"
1489 << "\nHNL mass = " << gCfgMassHNL << " [GeV]"
1490 << "\n|U_e4|^2 = " << gCfgECoupling
1491 << "\n|U_m4|^2 = " << gCfgMCoupling
1492 << "\n|U_t4|^2 = " << gCfgTCoupling
1493 << "\n"
1494 << "\nInteresting decay channels:";
1495 for( std::vector< HNLDecayMode_t >::iterator chit = gCfgIntChannels.begin();
1496 chit != gCfgIntChannels.end(); ++chit ){ csts << "\n\t" << utils::hnl::AsString(*chit); }
1497 csts << "\n"
1498 << "\nUser origin in beam coordinates = ( " << gCfgUserOx
1499 << ", " << gCfgUserOy << ", " << gCfgUserOz << " ) [" << lunits.c_str() << "]"
1500 << "\nEuler extrinsic x-z-x rotation = ( " << gCfgUserAx1
1501 << ", " << gCfgUserAz << ", " << gCfgUserAx2 << " ) [" << aunits.c_str() << "]"
1502 << "\nHNL particle-gun directional cosines: ( " << gCfgHNLCx << ", " << gCfgHNLCy
1503 << ", " << gCfgHNLCz << ") [ GeV / GeV ]";
1504
1505 LOG("gevald_hnl", pDEBUG) << csts.str();
1506
1507}
Heavy Neutral Lepton final-state product generator.
Definition HNLDecayer.h:41
std::vector< double > GetPGunDirection() const
std::string GetHNLInterestingChannels() const
double GetHNLMass() const
double GetHNLLifetime() const
bool IsHNLMajorana() const
double GetPGunEnergy() const
std::vector< double > GetHNLCouplings() const
Calculates HNL production kinematics & production vertex. Is a concrete implementation of the FluxRec...
std::vector< double > GetB2URotation() const
std::vector< double > GetB2UTranslation() const
double gOptEnergyHNL
double CoMLifetime
double gCfgUserOx
double gCfgUserOz
double gCfgTCoupling
double gCfgUserOy
double gCfgUserAx1
bool gCfgIsMajorana
double gCfgMCoupling
double gCfgHNLCz
double gCfgECoupling
std::vector< HNLDecayMode_t > gCfgIntChannels
double gCfgHNLCx
double gCfgUserAz
double gCfgUserAx2
double gCfgHNLCy
enum genie::hnl::EHNLDecayMode HNLDecayMode_t
static constexpr double m
Definition Units.h:71
static constexpr double rad
Definition Units.h:164

References genie::utils::hnl::AsString(), aunits, CoMLifetime, gCfgECoupling, gCfgHNLCx, gCfgHNLCy, gCfgHNLCz, gCfgIntChannels, gCfgIsMajorana, gCfgMassHNL, gCfgMCoupling, gCfgTCoupling, gCfgUserAx1, gCfgUserAx2, gCfgUserAz, gCfgUserOx, gCfgUserOy, gCfgUserOz, genie::AlgFactory::GetAlgorithm(), genie::hnl::FluxCreator::GetB2URotation(), genie::hnl::FluxCreator::GetB2UTranslation(), genie::hnl::Decayer::GetHNLCouplings(), genie::hnl::Decayer::GetHNLInterestingChannels(), genie::hnl::Decayer::GetHNLLifetime(), genie::hnl::Decayer::GetHNLMass(), genie::hnl::Decayer::GetPGunDirection(), genie::hnl::Decayer::GetPGunEnergy(), gOptEnergyHNL, gOptGeomAUnits, gOptGeomLUnits, genie::AlgFactory::Instance(), genie::hnl::Decayer::IsHNLMajorana(), LOG, lunits, genie::units::m, pDEBUG, pFATAL, and genie::units::rad.

Referenced by main().

◆ TestDecay()

int TestDecay ( void )

Definition at line 545 of file gBeamHNLValidationApp.cxx.

546{
547 string foutName("test_decay.root");
548
549 TFile * fout = TFile::Open( foutName.c_str(), "RECREATE" );
550
551 __attribute__((unused)) const EventRecordVisitorI * mcgen = HNLGenerator();
552 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
553 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
554
555 bool geom_is_accessible = ! (gSystem->AccessPathName(gOptRootGeom.c_str()));
556 if (!geom_is_accessible) {
557 LOG("gevald_hnl", pFATAL)
558 << "The specified ROOT geometry doesn't exist! Initialization failed!";
559 exit(1);
560 }
561
562 if( !gOptRootGeoManager ) gOptRootGeoManager = TGeoManager::Import(gOptRootGeom.c_str());
563
564 TGeoVolume * top_volume = gOptRootGeoManager->GetTopVolume();
565 assert( top_volume );
566 TGeoShape * ts = top_volume->GetShape();
567 __attribute__((unused)) TGeoBBox * box = (TGeoBBox *)ts;
568
569 LOG( "gevald_hnl", pDEBUG ) << "Imported box.";
570
571 const Algorithm * algVtxGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::VertexGenerator", "Default");
572
573 const VertexGenerator * vtxGen = dynamic_cast< const VertexGenerator * >( algVtxGen );
574 vtxGen->SetGeomFile( gOptRootGeom );
575
576 SimpleHNL sh = SimpleHNL( "HNLInstance", 0, kPdgHNL, kPdgKP,
578 std::map< HNLDecayMode_t, double > valMap = sh.GetValidChannels();
579 //const double CoMLifetime = sh.GetCoMLifetime();
580
581 assert( valMap.size() > 0 ); // must be able to decay to something!
582 assert( (*valMap.begin()).first == kHNLDcyNuNuNu );
583
584 LOG( "gevald_hnl", pINFO )
585 << "\n\nTesting decay modes for the HNL."
586 << "\nWill process gOptNev = " << gOptNev << " x ( N_channels = " << valMap.size()
587 << " ) = " << valMap.size() * gOptNev << " events, 1 for each valid channel."
588 << "\nWill produce 1 ROOT file ( " << foutName << " ) with:"
589 << "\n--> Energy spectrum for the decay products for each channel"
590 << "\n--> Rates of each decay channel";
591
592 // Set GHEP print level
593 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
594
595 // first set the 4-momentum of the HNL
596 //gOptEnergyHNL = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-Energy" );
597 gOptEnergyHNL = hnlgen->GetPGunEnergy();
598 double p3HNL = std::sqrt( gOptEnergyHNL * gOptEnergyHNL - gCfgMassHNL * gCfgMassHNL );
599 assert( p3HNL >= 0.0 );
600 TLorentzVector * p4HNL = new TLorentzVector( p3HNL * gCfgHNLCx,
601 p3HNL * gCfgHNLCy,
602 p3HNL * gCfgHNLCz, gOptEnergyHNL );
603
604 LOG( "gevald_hnl", pDEBUG )
605 << "\nUsing HNL with 4-momentum " << utils::print::P4AsString( p4HNL );
608
609 TLorentzVector * x4HNL = new TLorentzVector( 1.0, 2.0, 3.0, 0.0 ); // dummy
610
611 // now build array with indices of valid decay modes for speedy access
613 //double validRates[10] = { -1, -1, -1, -1, -1, -1, -1, -1, -1, -1 };
614 std::map< HNLDecayMode_t, double >::iterator vmit = valMap.begin(); int modeIdx = 0;
615 std::ostringstream msts;
616 for( ; vmit != valMap.end(); ++vmit ){
617 validModes[ modeIdx ] = (*vmit).first;
618 //validRates[ modeIdx ] = (*vmit).second;
619 msts << "\n" << utils::hnl::AsString( (*vmit).first );
620 modeIdx++;
621 }
622
623 LOG( "gevald_hnl", pDEBUG ) << "Here are the modes in order : " << msts.str();
624
625 // declare histos
626 // hSpectrum[i][j]: i iterates over HNLDecayMode_t, j over FS particle in same order as event record
627 TH1D hSpectrum[10][3], hCMSpectrum[10][3], hRates;
628 TH1D hParamSpace = TH1D( "hParamSpace", "Parameter space", 5, 0., 5. );
629
630 hParamSpace.SetBinContent( 1, 1000.0 * gCfgMassHNL ); // MeV
631 hParamSpace.SetBinContent( 2, gCfgECoupling );
632 hParamSpace.SetBinContent( 3, gCfgMCoupling );
633 hParamSpace.SetBinContent( 4, gCfgTCoupling );
634
635 hRates = TH1D( "hRates", "Rates of HNL decay channels", 10, 0, 10 );
636
637 std::string shortModes[10] = { "vvv", "vee", "vmue", "pi0v", "pie", "vmumu", "pimu", "pi0pi0v",
638 "pipi0e", "pipi0mu" };
639 std::string part0names[10] = { "v1", "v", "v", "pi0", "pi", "v", "pi", "pi01", "pi", "pi" };
640 std::string part1names[10] = { "v2", "e1", "mu", "v", "e", "mu1", "mu", "pi02", "pi0", "pi0" };
641 std::string part2names[10] = { "v3", "e2", "e", "None", "None", "mu2", "None", "v", "e", "mu" };
642 std::string partNames[3][10] = { part0names, part1names, part2names };
643 for( unsigned int iChan = 0; iChan < valMap.size(); iChan++ ){
644
645 std::string shortMode = shortModes[iChan];
646
647 for( Int_t iPart = 0; iPart < 3; iPart++ ){
648 std::string ParticleName = partNames[iPart][iChan];
649 if( strcmp( ParticleName.c_str(), "None" ) != 0 ){
650 hSpectrum[iChan][iPart] = TH1D( Form( "hSpectrum_%s_%s", shortMode.c_str(), ParticleName.c_str() ),
651 Form( "Fractional energy of particle: %s in decay: %s",
652 ParticleName.c_str(),
653 (utils::hnl::AsString( validModes[iChan] )).c_str() ),
654 100, 0., 1.0 );
655 hCMSpectrum[iChan][iPart] = TH1D( Form( "hCMSpectrum_%s_%s", shortMode.c_str(), ParticleName.c_str() ),
656 Form( "Rest frame energy of particle: %s in decay: %s",
657 ParticleName.c_str(),
658 (utils::hnl::AsString( validModes[iChan] )).c_str() ),
659 100, 0., gCfgMassHNL );
660 } // only declare histos of particles that exist in decay
661 // and are of allowed decays
662 }
663 }
664
665 // fill hRates now
666 double rnununu = ( (*valMap.find( kHNLDcyNuNuNu )) ).second;
667 double rnuee = ( valMap.find( kHNLDcyNuEE ) != valMap.end() ) ? (*(valMap.find( kHNLDcyNuEE ))).second : 0.0;
668 double rnumue = ( valMap.find( kHNLDcyNuMuE ) != valMap.end() ) ? (*(valMap.find( kHNLDcyNuMuE ))).second : 0.0;
669 double rpi0nu = ( valMap.find( kHNLDcyPi0Nu ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPi0Nu ))).second : 0.0;
670 double rpie = ( valMap.find( kHNLDcyPiE ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPiE ))).second : 0.0;
671 double rnumumu = ( valMap.find( kHNLDcyNuMuMu ) != valMap.end() ) ? (*(valMap.find( kHNLDcyNuMuMu ))).second : 0.0;
672 double rpimu = ( valMap.find( kHNLDcyPiMu ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPiMu ))).second : 0.0;
673 double rpi0pi0nu = ( valMap.find( kHNLDcyPi0Pi0Nu ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPi0Pi0Nu ))).second : 0.0;
674 double rpipi0e = ( valMap.find( kHNLDcyPiPi0E ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPiPi0E ))).second : 0.0;
675 double rpipi0mu = ( valMap.find( kHNLDcyPiPi0Mu ) != valMap.end() ) ? (*(valMap.find( kHNLDcyPiPi0Mu ))).second : 0.0;
676
677 hRates.SetBinContent( 1, rpimu );
678 hRates.SetBinContent( 2, rpie );
679 hRates.SetBinContent( 3, rpi0nu );
680 hRates.SetBinContent( 4, rnununu );
681 hRates.SetBinContent( 5, rnumumu );
682 hRates.SetBinContent( 6, rnuee );
683 hRates.SetBinContent( 7, rnumue );
684 hRates.SetBinContent( 8, rpipi0e );
685 hRates.SetBinContent( 9, rpipi0mu );
686 hRates.SetBinContent( 10, rpi0pi0nu );
687
688 int ievent = 0;
689 while( true ){
690
691 if( gOptNev >= 10000 ){
692 if( ievent % (gOptNev / 1000) == 0 ){
693 int irat = ievent / (gOptNev / 1000);
694 std::cerr << Form("%2.2f", 0.1 * irat) << " % ( " << ievent << " / "
695 << gOptNev << " ) \r" << std::flush;
696 }
697 } else if( gOptNev >= 100 ) {
698 if( ievent % (gOptNev / 10) == 0 ){
699 int irat = ievent / ( gOptNev / 10 );
700 std::cerr << 10.0 * irat << " % " << " ( " << ievent
701 << " / " << gOptNev << " ) \r" << std::flush;
702 }
703 }
704
705 if( ievent == gOptNev ){ std::cerr << " \n"; break; }
706
707 ostringstream asts;
708 for( unsigned int iMode = 0; iMode < valMap.size(); iMode++ ){
709 if( ievent == 0 ){
710 asts
711 << "\nDecay mode " << iMode << " is " << utils::hnl::AsString( validModes[ iMode ] );
712 }
713
714 // build an event
715 EventRecord * event = new EventRecord;
716 Interaction * interaction = Interaction::HNL( genie::kPdgHNL, gOptEnergyHNL, validModes[ iMode ] );
717 // set Particle(0) and a dummy vertex
718 GHepParticle ptHNL( kPdgHNL, kIStInitialState, -1, -1, -1, -1, *p4HNL, *x4HNL );
719 event->AddParticle( ptHNL );
720 interaction->InitStatePtr()->SetProbeP4( *p4HNL );
721 event->SetVertex( *x4HNL );
722
723 event->SetProbability( CoMLifetime );
724
725 event->AttachSummary( interaction );
726 LOG( "gevald_hnl", pDEBUG )
727 << "Simulating decay with mode "
728 << utils::hnl::AsString( (HNLDecayMode_t) interaction->ExclTag().DecayMode() );
729
730 // simulate the decay
731 hnlgen->ProcessEventRecord( event );
732
733 // now get a weight.
734 // = exp( - T_{box} / \tau_{HNL} ) = exp( - L_{box} / ( \beta_{HNL} \gamma_{HNL} c ) * h / \Gamma_{tot} )
735 // placing the HNL at a point configured by the user
736 std::vector< double > PGOrigin = hnlgen->GetPGunOrigin();
737 double ox = PGOrigin.at(0), oy = PGOrigin.at(1), oz = PGOrigin.at(2);
738 /*
739 double ox = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginX" );
740 double oy = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginY" );
741 double oz = utils::hnl::GetCfgDouble( "HNL", "ParticleGun", "PG-OriginZ" );
742 */
743 ox *= units::m / units::mm; oy *= units::m / units::mm; oz *= units::m / units::mm;
744 TVector3 startPoint( ox, oy, oz ); TVector3 entryPoint, exitPoint;
745 TLorentzVector tmpVtx( ox, oy, oz, 0.0 );
746 event->SetVertex( tmpVtx );
747
748 vtxGen->ProcessEventRecord(event);
749
750 LOG( "gevald_hnl", pDEBUG ) << *event;
751
752 // now fill the histos!
753 double wgt = event->Weight();
754 hSpectrum[iMode][0].Fill( (event->Particle(1))->E() / gOptEnergyHNL, wgt );
755 hSpectrum[iMode][1].Fill( (event->Particle(2))->E() / gOptEnergyHNL, wgt );
756 if( event->Particle(3) ) hSpectrum[iMode][2].Fill( (event->Particle(3))->E() / gOptEnergyHNL, wgt );
757
758 // let's also fill the CM spectra
759 // get particle 4-momenta and boost back to rest frame!
760 TLorentzVector * p4p1 = (event->Particle(1))->GetP4();
761 TLorentzVector * p4p2 = (event->Particle(2))->GetP4();
762 TLorentzVector * p4p3 = 0;
763 if( event->Particle(3) ) p4p3 = (event->Particle(3))->GetP4();
764
765 TVector3 boostVec = p4HNL->BoostVector();
766
767 p4p1->Boost( -boostVec );
768 p4p2->Boost( -boostVec );
769 if( p4p3 ) p4p3->Boost( -boostVec );
770
771 hCMSpectrum[iMode][0].Fill( p4p1->E(), wgt );
772 hCMSpectrum[iMode][1].Fill( p4p2->E(), wgt );
773 if( p4p3 ) hCMSpectrum[iMode][2].Fill( p4p3->E(), wgt );
774
775 // clean-up
776 delete event;
777
778 } // loop over valid decay channels
779
780 if( ievent == 0 ) LOG( "gevald_hnl", pDEBUG ) << asts.str();
781 LOG( "gevald_hnl", pDEBUG ) << "Finished event: " << ievent;
782
783 ievent++;
784
785 } // event loop
786
787 fout->cd();
788 hParamSpace.Write();
789 hRates.Write();
790 for( unsigned int i = 0; i < valMap.size(); i++ ){
791 for( Int_t j = 0; j < 3; j++ ){
792 std::string ParticleName = partNames[j][i];
793 if( strcmp( ParticleName.c_str(), "None" ) != 0 ){
794 hSpectrum[i][j].Write();
795 hCMSpectrum[i][j].Write();
796 }
797 }
798 }
799 fout->Write();
800 fout->Close();
801
802 delete p4HNL;
803 delete x4HNL;
804 return 0;
805}
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)
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition Interaction.h:56
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
static Interaction * HNL(int probe, double E=0, int decayed_mode=-1)
int DecayMode(void) const
Definition XclsTag.h:70
void ProcessEventRecord(GHepRecord *event) const
std::vector< double > GetPGunOrigin() const
void SetMomentumDirection(double ux, double uy, double uz)
Definition SimpleHNL.h:231
std::map< genie::hnl::HNLDecayMode_t, double > GetValidChannels()
Definition SimpleHNL.h:154
void SetEnergy(const double E)
Definition SimpleHNL.h:183
Heavy Neutral Lepton decay vertex generator given production vertex and momentum ***.
void SetGeomFile(std::string geomfile) const
void ProcessEventRecord(GHepRecord *event_rec) const
const EventRecordVisitorI * HNLGenerator(void)
static constexpr double mm
Definition Units.h:65
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::utils::hnl::AsString(), CoMLifetime, genie::XclsTag::DecayMode(), genie::Interaction::ExclTag(), gCfgECoupling, gCfgHNLCx, gCfgHNLCy, gCfgHNLCz, gCfgMassHNL, gCfgMCoupling, gCfgTCoupling, genie::AlgFactory::GetAlgorithm(), genie::hnl::Decayer::GetPGunEnergy(), genie::hnl::Decayer::GetPGunOrigin(), genie::hnl::SimpleHNL::GetValidChannels(), gOptEnergyHNL, gOptNev, gOptRootGeom, genie::Interaction::HNL(), HNLGenerator(), genie::Interaction::InitStatePtr(), genie::AlgFactory::Instance(), genie::RunOpt::Instance(), genie::hnl::kHNLDcyNuEE, genie::hnl::kHNLDcyNull, genie::hnl::kHNLDcyNuMuE, genie::hnl::kHNLDcyNuMuMu, genie::hnl::kHNLDcyNuNuNu, genie::hnl::kHNLDcyPi0Nu, genie::hnl::kHNLDcyPi0Pi0Nu, genie::hnl::kHNLDcyPiE, genie::hnl::kHNLDcyPiMu, genie::hnl::kHNLDcyPiPi0E, genie::hnl::kHNLDcyPiPi0Mu, genie::kIStInitialState, genie::kPdgHNL, genie::kPdgKP, LOG, genie::units::m, genie::units::mm, genie::utils::print::P4AsString(), pDEBUG, pFATAL, pINFO, genie::hnl::Decayer::ProcessEventRecord(), genie::hnl::VertexGenerator::ProcessEventRecord(), genie::hnl::SimpleHNL::SetEnergy(), genie::hnl::VertexGenerator::SetGeomFile(), genie::hnl::SimpleHNL::SetMomentumDirection(), genie::GHepRecord::SetPrintLevel(), and genie::InitialState::SetProbeP4().

Referenced by main().

Variable Documentation

◆ aunits

string aunits

Definition at line 234 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), and ReadInConfig().

◆ CoMLifetime

double CoMLifetime = -1

Definition at line 190 of file gBeamHNLValidationApp.cxx.

◆ fdx

double fdx = 0

Definition at line 245 of file gBeamHNLValidationApp.cxx.

◆ fdy

double fdy = 0

Definition at line 246 of file gBeamHNLValidationApp.cxx.

◆ fdz

double fdz = 0

Definition at line 247 of file gBeamHNLValidationApp.cxx.

◆ fox

double fox = 0

Definition at line 248 of file gBeamHNLValidationApp.cxx.

◆ foy

double foy = 0

Definition at line 249 of file gBeamHNLValidationApp.cxx.

◆ foz

double foz = 0

Definition at line 250 of file gBeamHNLValidationApp.cxx.

◆ gCfgBoxLx

double gCfgBoxLx = -1

Definition at line 204 of file gBeamHNLValidationApp.cxx.

◆ gCfgBoxLy

double gCfgBoxLy = -1

Definition at line 205 of file gBeamHNLValidationApp.cxx.

◆ gCfgBoxLz

double gCfgBoxLz = -1

Definition at line 206 of file gBeamHNLValidationApp.cxx.

◆ gCfgDecayMode

HNLDecayMode_t gCfgDecayMode = kHNLDcyNull

Definition at line 194 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs().

◆ gCfgECoupling

double gCfgECoupling = -1

Definition at line 184 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgHNLCx

double gCfgHNLCx = -1

Definition at line 215 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgHNLCy

double gCfgHNLCy = -1

Definition at line 216 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgHNLCz

double gCfgHNLCz = -1

Definition at line 217 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgHNLKind

int gCfgHNLKind = 2

Definition at line 188 of file gBeamHNLValidationApp.cxx.

◆ gCfgIntChannels

std::vector< HNLDecayMode_t > gCfgIntChannels

Definition at line 195 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgIsMajorana

bool gCfgIsMajorana = false

Definition at line 187 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgMassHNL

double gCfgMassHNL = -1

Definition at line 183 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), ReadInConfig(), and TestDecay().

◆ gCfgMCoupling

double gCfgMCoupling = -1

Definition at line 185 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgParentEnergy

double gCfgParentEnergy = -1

Definition at line 208 of file gBeamHNLValidationApp.cxx.

◆ gCfgParentOx

double gCfgParentOx = -1

Definition at line 211 of file gBeamHNLValidationApp.cxx.

◆ gCfgParentOy

double gCfgParentOy = -1

Definition at line 212 of file gBeamHNLValidationApp.cxx.

◆ gCfgParentOz

double gCfgParentOz = -1

Definition at line 213 of file gBeamHNLValidationApp.cxx.

◆ gCfgParentPhi

double gCfgParentPhi = -1

Definition at line 210 of file gBeamHNLValidationApp.cxx.

◆ gCfgParentTheta

double gCfgParentTheta = -1

Definition at line 209 of file gBeamHNLValidationApp.cxx.

◆ gCfgProdMode

HNLProd_t gCfgProdMode = kHNLProdNull

Definition at line 193 of file gBeamHNLValidationApp.cxx.

◆ gCfgTCoupling

double gCfgTCoupling = -1

Definition at line 186 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig(), and TestDecay().

◆ gCfgUserAx1

double gCfgUserAx1 = -1

Definition at line 200 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgUserAx2

double gCfgUserAx2 = -1

Definition at line 202 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgUserAz

double gCfgUserAz = -1

Definition at line 201 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgUserOx

double gCfgUserOx = -1

Definition at line 197 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgUserOy

double gCfgUserOy = -1

Definition at line 198 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ gCfgUserOz

double gCfgUserOz = -1

Definition at line 199 of file gBeamHNLValidationApp.cxx.

Referenced by ReadInConfig().

◆ geom

string geom

Definition at line 233 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), GetTargetCodes(), and main().

◆ gOptEnergyHNL

double gOptEnergyHNL = -1

Definition at line 221 of file gBeamHNLValidationApp.cxx.

◆ gOptEvFilePrefix

string gOptEvFilePrefix = kDefOptEvFilePrefix

Definition at line 229 of file gBeamHNLValidationApp.cxx.

◆ gOptGeomAUnits

double gOptGeomAUnits = 0

Definition at line 236 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), and ReadInConfig().

◆ gOptGeomLUnits

double gOptGeomLUnits = 0

Definition at line 235 of file gBeamHNLValidationApp.cxx.

◆ gOptGeomTUnits

double gOptGeomTUnits = 0

Definition at line 237 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs().

◆ gOptIsMonoEnFlux

bool gOptIsMonoEnFlux = true

Definition at line 227 of file gBeamHNLValidationApp.cxx.

◆ gOptNev

int gOptNev = 10

Definition at line 177 of file gBeamHNLValidationApp.cxx.

◆ gOptRanSeed

long int gOptRanSeed = -1

Definition at line 252 of file gBeamHNLValidationApp.cxx.

◆ gOptRootGeom

string gOptRootGeom

Definition at line 231 of file gBeamHNLValidationApp.cxx.

◆ gOptRootGeomTopVol

string gOptRootGeomTopVol = ""

Definition at line 242 of file gBeamHNLValidationApp.cxx.

◆ gOptRunNu

Long_t gOptRunNu = 1000

Definition at line 176 of file gBeamHNLValidationApp.cxx.

◆ gOptUsingRootGeom

bool gOptUsingRootGeom = false

Definition at line 230 of file gBeamHNLValidationApp.cxx.

◆ gOptValidationMode

HNLValidation_t gOptValidationMode = kValNone

Definition at line 179 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), and main().

◆ kDefOptEvFilePrefix

string kDefOptEvFilePrefix = "gntp"

Definition at line 168 of file gBeamHNLValidationApp.cxx.

◆ kDefOptFluxFilePath

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

Definition at line 169 of file gBeamHNLValidationApp.cxx.

◆ kDefOptGeomAUnits

string kDefOptGeomAUnits = "rad"

Definition at line 163 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs().

◆ kDefOptGeomDUnits

string kDefOptGeomDUnits = "g_cm3"

Definition at line 165 of file gBeamHNLValidationApp.cxx.

◆ kDefOptGeomLUnits

string kDefOptGeomLUnits = "mm"

Definition at line 162 of file gBeamHNLValidationApp.cxx.

◆ kDefOptGeomTUnits

string kDefOptGeomTUnits = "ns"

Definition at line 164 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs().

◆ kDefOptNtpFormat

NtpMCFormat_t kDefOptNtpFormat = kNFGHEP

Definition at line 167 of file gBeamHNLValidationApp.cxx.

◆ kDefOptSConfig

string kDefOptSConfig = "BeamHNL"

Definition at line 172 of file gBeamHNLValidationApp.cxx.

◆ kDefOptSName

string kDefOptSName = "genie::EventGenerator"

Definition at line 171 of file gBeamHNLValidationApp.cxx.

◆ kDefOptSTune

string kDefOptSTune = "GHNL20_00a_00_000"

Definition at line 173 of file gBeamHNLValidationApp.cxx.

◆ lunits

string lunits

Definition at line 234 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs(), and ReadInConfig().

◆ tunits

string tunits

Definition at line 234 of file gBeamHNLValidationApp.cxx.

Referenced by GetCommandLineArgs().