GENIEGenerator
Loading...
Searching...
No Matches
gNucleonDecayEvGen.cxx
Go to the documentation of this file.
1//________________________________________________________________________________________
2/*!
3
4\program gevgen_ndcy
5
6\brief A GENIE-based nucleon decay event generation application.
7
8 *** Synopsis :
9
10 gevgen_ndcy [-h]
11 [-r run#]
12 -n n_of_events
13 -m decay_mode
14 [-N decayed_nucleon_pdg]
15 -g geometry
16 [-L geometry_length_units]
17 [-D geometry_density_units]
18 [-t geometry_top_volume_name]
19 [-o output_event_file_prefix]
20 [--seed random_number_seed]
21 [--message-thresholds xml_file]
22 [--event-record-print-level level]
23 [--mc-job-status-refresh-rate rate]
24
25 *** Options :
26
27 [] Denotes an optional argument
28
29 -h
30 Prints out the gevgen_ndcy syntax and exits.
31 -r
32 Specifies the MC run number [default: 1000].
33 -n
34 Specifies how many events to generate.
35 -m
36 Nucleon decay mode ID:
37 see http://www-pdg.lbl.gov/2016/listings/rpp2016-list-p.pdf
38 for nucleon decay mode numbering convention
39 Example:
40 m = 1: N -> e+ pi
41 ...
42 m = 60: n -> 5nu
43
44
45 -N
46 Decayed nucleon PDG code.
47 Either N=2212 (proton) or N=2112 (neutron)
48 Example:
49 m = 1 and N = 2112: n -> e+ pi-
50 m = 1 and N = 2212: p -> e+ pi0
51
52 -g
53 Input 'geometry'.
54 This option can be used to specify any of:
55 1 > A ROOT file containing a ROOT/GEANT geometry description
56 [Examples]
57 - To use the master volume from the ROOT geometry stored
58 in the laguna-lbno.root file, type:
59 '-g /some/path/laguna-lbno.root'
60 2 > A mix of target materials, each with its corresponding weight,
61 typed as a comma-separated list of nuclear PDG codes (in the
62 std PDG2006 convention: 10LZZZAAAI) with the weight fractions
63 in brackets, eg code1[fraction1],code2[fraction2],...
64 If that option is used (no detailed input geometry description)
65 then the interaction vertices are distributed in the detector
66 by the detector MC.
67 [Examples]
68 - To use a target mix of 88.9% O16 and 11.1% Hydrogen type:
69 '-g 1000080160[0.889],1000010010[0.111]'
70 -L
71 Input geometry length units, eg 'm', 'cm', 'mm', ...
72 [default: 'mm']
73 -D
74 Input geometry density units, eg 'g_cm3', 'clhep_def_density_unit',...
75 [default: 'g_cm3']
76 -t
77 Input 'top volume' for event generation.
78 The option be used to force event generation in given sub-detector.
79 [default: the 'master volume' of the input geometry]
80 You can also use the -t option to switch generation on/off at
81 multiple volumes as, for example, in:
82 `-t +Vol1-Vol2+Vol3-Vol4',
83 `-t "+Vol1 -Vol2 +Vol3 -Vol4"',
84 `-t -Vol2-Vol4+Vol1+Vol3',
85 `-t "-Vol2 -Vol4 +Vol1 +Vol3"'m
86 where:
87 "+Vol1" and "+Vol3" tells GENIE to `switch on' Vol1 and Vol3, while
88 "-Vol2" and "-Vol4" tells GENIE to `switch off' Vol2 and Vol4.
89 If the very first character is a '+', GENIE will neglect all volumes
90 except the ones explicitly turned on. Vice versa, if the very first
91 character is a `-', GENIE will keep all volumes except the ones
92 explicitly turned off (feature contributed by J.Holeczek).
93 -o
94 Sets the prefix of the output event file.
95 The output filename is built as:
96 [prefix].[run_number].[event_tree_format].[file_format]
97 The default output filename is:
98 gntp.[run_number].ghep.root
99 This cmd line arguments lets you override 'gntp'
100 --seed
101 Random number seed.
102
103\author Costas Andreopoulos <c.andreopoulos \at cern.ch>
104 University of Liverpool
105
106\created November 03, 2011
107
108\cpright Copyright (c) 2003-2025, The GENIE Collaboration
109 For the full text of the license visit http://copyright.genie-mc.org
110
111
112*/
113//_________________________________________________________________________________________
114
115#include <cassert>
116#include <cstdlib>
117#include <string>
118#include <vector>
119#include <sstream>
120
121#include <TSystem.h>
122
142
143using std::string;
144using std::vector;
145using std::ostringstream;
146
147using namespace genie;
148
149// function prototypes
150void GetCommandLineArgs (int argc, char ** argv);
151void PrintSyntax (void);
152int SelectInitState (void);
154
155//
156string kDefOptGeomLUnits = "mm"; // default geometry length units
157string kDefOptGeomDUnits = "g_cm3"; // default geometry density units
158NtpMCFormat_t kDefOptNtpFormat = kNFGHEP; // default event tree format
159string kDefOptEvFilePrefix = "gntp";
160
161//
162Long_t gOptRunNu = 1000; // run number
163int gOptNev = 10; // number of events to generate
165int gOptDecayedNucleon = kNDNull; // decayed nucleon PDG
166string gOptEvFilePrefix = kDefOptEvFilePrefix; // event file prefix
167bool gOptUsingRootGeom = false; // using root geom or target mix?
168map<int,double> gOptTgtMix; // target mix (tgt pdg -> wght frac) / if not using detailed root geom
169string gOptRootGeom; // input ROOT file with realistic detector geometry
170string gOptRootGeomTopVol = ""; // input geometry top event generation volume
171double gOptGeomLUnits = 0; // input geometry length units
172double gOptGeomDUnits = 0; // input geometry density units
173long int gOptRanSeed = -1; // random number seed
174
175//_________________________________________________________________________________________
176int main(int argc, char ** argv)
177{
178 // Parse command line arguments
179 GetCommandLineArgs(argc,argv);
180
181 // Init messenger and random number seed
182 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
184
185 // Initialize an Ntuple Writer to save GHEP records into a TTree
188 ntpw.Initialize();
189
190 // Create a MC job monitor for a periodically updated status file
191 GMCJMonitor mcjmonitor(gOptRunNu);
192 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
193
194 // Set GHEP print level
195 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
196
197 // Get the nucleon decay generator
199
200 // Event loop
201 int ievent = 0;
202 int dpdg = 0;
203 if (gOptDecayedNucleon > 0) {
204 dpdg = gOptDecayedNucleon;
205 } else {
207 }
208
209 while (1)
210 {
211 if(ievent == gOptNev) break;
212
213 LOG("gevgen_ndcy", pNOTICE)
214 << " *** Generating event............ " << ievent;
215
216 EventRecord * event = new EventRecord;
217 int target = SelectInitState();
218 int decay = (int)gOptDecayMode;
219 Interaction * interaction = Interaction::NDecay(target,decay,dpdg);
220 event->AttachSummary(interaction);
221
222 // Simulate decay
223 mcgen->ProcessEventRecord(event);
224
225 LOG("gevgen_ndcy", pINFO)
226 << "Generated event: " << *event;
227
228 // Add event at the output ntuple, refresh the mc job monitor & clean-up
229 ntpw.AddEventRecord(ievent, event);
230 mcjmonitor.Update(ievent,event);
231 delete event;
232
233 ievent++;
234 } // event loop
235
236 // Save the generated event tree & close the output file
237 ntpw.Save();
238
239 LOG("gevgen_ndcy", pNOTICE) << "Done!";
240
241 return 0;
242}
243//_________________________________________________________________________________________
245{
247 LOG("gevgen_ndcy", pFATAL) << "Not a valid decay mode and/or decayed nucleon...";
248 gAbortingInErr = true;
249 exit(1);
250 }
251
252 int dpdg = 0;
253 if (gOptDecayedNucleon > 0) {
254 dpdg = gOptDecayedNucleon;
255 } else {
257 }
258
259 map<int,double> cprob; // cumulative probability
260 map<int,double>::const_iterator iter;
261
262 double sum_prob = 0;
263 for(iter = gOptTgtMix.begin(); iter != gOptTgtMix.end(); ++iter) {
264 int pdg_code = iter->first;
265 int A = pdg::IonPdgCodeToA(pdg_code);
266 int Z = pdg::IonPdgCodeToZ(pdg_code);
267
268 double nucleon_decay_fraction = 0.;
269 if (dpdg == kPdgProton ) { nucleon_decay_fraction = (double)Z / (double)A; }
270 else if (dpdg == kPdgNeutron ) { nucleon_decay_fraction = (double)(A-Z) / (double)A; }
271
272 double wgt = iter->second;
273 double prob = wgt*nucleon_decay_fraction;
274
275 sum_prob += prob;
276 cprob.insert(map<int, double>::value_type(pdg_code, sum_prob));
277 }
278
279 assert(sum_prob > 0.);
280
282 double r = sum_prob * rnd->RndEvg().Rndm();
283
284 for(iter = cprob.begin(); iter != cprob.end(); ++iter) {
285 int pdg_code = iter->first;
286 double prob = iter->second;
287 if(r < prob) {
288 LOG("gevgen_ndcy", pNOTICE) << "Selected initial state = " << pdg_code;
289 return pdg_code;
290 }
291 }
292
293 LOG("gevgen_ndcy", pFATAL) << "Couldn't select an initial state...";
294 gAbortingInErr = true;
295 exit(1);
296}
297//_________________________________________________________________________________________
299{
300 string sname = "genie::EventGenerator";
301 string sconfig = "NucleonDecay";
303 const EventRecordVisitorI * mcgen =
304 dynamic_cast<const EventRecordVisitorI *> (algf->GetAlgorithm(sname,sconfig));
305 if(!mcgen) {
306 LOG("gevgen_ndcy", pFATAL) << "Couldn't instantiate the nucleon decay generator";
307 gAbortingInErr = true;
308 exit(1);
309 }
310 return mcgen;
311}
312//_________________________________________________________________________________________
313void GetCommandLineArgs(int argc, char ** argv)
314{
315 LOG("gevgen_ndcy", pINFO) << "Parsing command line arguments";
316
317 // Common run options.
319
320 // Parse run options for this app
321
322 CmdLnArgParser parser(argc,argv);
323
324 // help?
325 bool help = parser.OptionExists('h');
326 if(help) {
327 PrintSyntax();
328 exit(0);
329 }
330
331 // run number
332 if( parser.OptionExists('r') ) {
333 LOG("gevgen_ndcy", pDEBUG) << "Reading MC run number";
334 gOptRunNu = parser.ArgAsLong('r');
335 } else {
336 LOG("gevgen_ndcy", pDEBUG) << "Unspecified run number - Using default";
337 gOptRunNu = 1000;
338 } //-r
339
340
341 // number of events
342 if( parser.OptionExists('n') ) {
343 LOG("gevgen_ndcy", pDEBUG)
344 << "Reading number of events to generate";
345 gOptNev = parser.ArgAsInt('n');
346 } else {
347 LOG("gevgen_ndcy", pFATAL)
348 << "You need to specify the number of events";
349 PrintSyntax();
350 exit(0);
351 } //-n
352
353 // decay mode
354 int mode = -1;
355 if( parser.OptionExists('m') ) {
356 LOG("gevgen_ndcy", pDEBUG)
357 << "Reading decay mode";
358 mode = parser.ArgAsInt('m');
359 } else {
360 LOG("gevgen_ndcy", pFATAL)
361 << "You need to specify the decay mode";
362 PrintSyntax();
363 exit(0);
364 } //-m
366
367 // decayed nucleon PDG
368 if( parser.OptionExists('N') ) {
369 LOG("gevgen_ndcy", pINFO) << "Reading decayed nucleon PDG";
370 gOptDecayedNucleon = parser.ArgAsInt('N');
371 } else {
372 LOG("gevgen_ndcy", pINFO) << "Unspecified decayed nucleon PDG - Using default";
374 }
375
377 if(!valid_mode) {
378 LOG("gevgen_ndcy", pFATAL)
379 << "You need to specify a valid decay mode / decayed nucleon PDG combination";
380 PrintSyntax();
381 exit(0);
382 }
383
384 //
385 // geometry
386 //
387
388 string geom = "";
389 string lunits, dunits;
390 if( parser.OptionExists('g') ) {
391 LOG("gevgen_ndcy", pDEBUG) << "Getting input geometry";
392 geom = parser.ArgAsString('g');
393
394 // is it a ROOT file that contains a ROOT geometry?
395 bool accessible_geom_file =
396 ! (gSystem->AccessPathName(geom.c_str()));
397 if (accessible_geom_file) {
399 gOptUsingRootGeom = true;
400 }
401 } else {
402 LOG("gevgen_ndcy", pFATAL)
403 << "No geometry option specified - Exiting";
404 PrintSyntax();
405 exit(1);
406 } //-g
407
409 // using a ROOT geometry - get requested geometry units
410
411 // legth units:
412 if( parser.OptionExists('L') ) {
413 LOG("gevgen_ndcy", pDEBUG)
414 << "Checking for input geometry length units";
415 lunits = parser.ArgAsString('L');
416 } else {
417 LOG("gevgen_ndcy", pDEBUG) << "Using default geometry length units";
419 } // -L
420 // density units:
421 if( parser.OptionExists('D') ) {
422 LOG("gevgen_ndcy", pDEBUG)
423 << "Checking for input geometry density units";
424 dunits = parser.ArgAsString('D');
425 } else {
426 LOG("gevgen_ndcy", pDEBUG) << "Using default geometry density units";
427 dunits = kDefOptGeomDUnits;
428 } // -D
431
432 // check whether an event generation volume name has been
433 // specified -- default is the 'top volume'
434 if( parser.OptionExists('t') ) {
435 LOG("gevgen_ndcy", pDEBUG) << "Checking for input volume name";
436 gOptRootGeomTopVol = parser.ArgAsString('t');
437 } else {
438 LOG("gevgen_ndcy", pDEBUG) << "Using the <master volume>";
439 } // -t
440
441 } // using root geom?
442
443 else {
444 // User has specified a target mix.
445 // Decode the list of target pdf codes & their corresponding weight fraction
446 // (specified as 'pdg_code_1[fraction_1],pdg_code_2[fraction_2],...')
447 // See documentation on top section of this file.
448 //
449 gOptTgtMix.clear();
450 vector<string> tgtmix = utils::str::Split(geom,",");
451 if(tgtmix.size()==1) {
452 int pdg = atoi(tgtmix[0].c_str());
453 double wgt = 1.0;
454 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
455 } else {
456 vector<string>::const_iterator tgtmix_iter = tgtmix.begin();
457 for( ; tgtmix_iter != tgtmix.end(); ++tgtmix_iter) {
458 string tgt_with_wgt = *tgtmix_iter;
459 string::size_type open_bracket = tgt_with_wgt.find("[");
460 string::size_type close_bracket = tgt_with_wgt.find("]");
461 if (open_bracket ==string::npos ||
462 close_bracket==string::npos)
463 {
464 LOG("gevgen_ndcy", pFATAL)
465 << "You made an error in specifying the target mix";
466 PrintSyntax();
467 exit(1);
468 }
469 string::size_type ibeg = 0;
470 string::size_type iend = open_bracket;
471 string::size_type jbeg = open_bracket+1;
472 string::size_type jend = close_bracket;
473 int pdg = atoi(tgt_with_wgt.substr(ibeg,iend-ibeg).c_str());
474 double wgt = atof(tgt_with_wgt.substr(jbeg,jend-jbeg).c_str());
475 LOG("gevgen_ndcy", pDEBUG)
476 << "Adding to target mix: pdg = " << pdg << ", wgt = " << wgt;
477 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
478
479 }// tgtmix_iter
480 } // >1 materials in mix
481 } // using tgt mix?
482
483 // event file prefix
484 if( parser.OptionExists('o') ) {
485 LOG("gevgen_ndcy", pDEBUG) << "Reading the event filename prefix";
486 gOptEvFilePrefix = parser.ArgAsString('o');
487 } else {
488 LOG("gevgen_ndcy", pDEBUG)
489 << "Will set the default event filename prefix";
491 } //-o
492
493
494 // random number seed
495 if( parser.OptionExists("seed") ) {
496 LOG("gevgen_ndcy", pINFO) << "Reading random number seed";
497 gOptRanSeed = parser.ArgAsLong("seed");
498 } else {
499 LOG("gevgen_ndcy", pINFO) << "Unspecified random number seed - Using default";
500 gOptRanSeed = -1;
501 }
502
503 //
504 // >>> print the command line options
505 //
506
507 PDGLibrary * pdglib = PDGLibrary::Instance();
508
509 ostringstream gminfo;
510 if (gOptUsingRootGeom) {
511 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
512 << ", top volume: "
513 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
514 << ", length units: " << lunits
515 << ", density units: " << dunits;
516 } else {
517 gminfo << "Using target mix - ";
518 map<int,double>::const_iterator iter;
519 for(iter = gOptTgtMix.begin(); iter != gOptTgtMix.end(); ++iter) {
520 int pdg_code = iter->first;
521 double wgt = iter->second;
522 TParticlePDG * p = pdglib->Find(pdg_code);
523 if(p) {
524 string name = p->GetName();
525 gminfo << "(" << name << ") -> " << 100*wgt << "% / ";
526 }//p?
527 }
528 }
529
530 LOG("gevgen_ndcy", pNOTICE)
531 << "\n\n"
532 << utils::print::PrintFramedMesg("gevgen_ndcy job configuration");
533
534 LOG("gevgen_ndcy", pNOTICE)
535 << "\n @@ Run number: " << gOptRunNu
536 << "\n @@ Random number seed: " << gOptRanSeed
538 << "\n @@ Geometry $ " << gminfo.str()
539 << "\n @@ Statistics $ " << gOptNev << " events";
540
541 //
542 // Temporary warnings...
543 //
545 LOG("gevgen_ndcy", pWARN)
546 << "\n ** ROOT geometries not supported yet in the nucleon decay mode"
547 << "\n ** (they will be in the very near future)"
548 << "\n ** Please specify a `target mix' instead.";
549 gAbortingInErr = true;
550 exit(1);
551 }
552}
553//_________________________________________________________________________________________
554void PrintSyntax(void)
555{
556 LOG("gevgen_ndcy", pFATAL)
557 << "\n **Syntax**"
558 << "\n gevgen_ndcy [-h] "
559 << "\n [-r run#]"
560 << "\n -m decay_mode"
561 << "\n [-N decayed_nucleon]"
562 << "\n -g geometry"
563 << "\n [-t top_volume_name_at_geom]"
564 << "\n [-L length_units_at_geom]"
565 << "\n [-D density_units_at_geom]"
566 << "\n -n n_of_events "
567 << "\n [-o output_event_file_prefix]"
568 << "\n [--seed random_number_seed]"
569 << "\n [--message-thresholds xml_file]"
570 << "\n [--event-record-print-level level]"
571 << "\n [--mc-job-status-refresh-rate rate]"
572 << "\n"
573 << " Please also read the detailed documentation at http://www.genie-mc.org"
574 << " or look at the source code: $GENIE/src/Apps/gNucleonDecayEvGen.cxx"
575 << "\n";
576}
577//_________________________________________________________________________________________
578
#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
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
int main()
The GENIE Algorithm Factory.
Definition AlgFactory.h:39
const Algorithm * GetAlgorithm(const AlgId &algid)
static AlgFactory * Instance()
Command line argument parser.
string ArgAsString(char opt)
bool OptionExists(char opt)
was option set?
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the 'Visito...
virtual void ProcessEventRecord(GHepRecord *event_rec) const =0
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition EventRecord.h:37
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 Update(int iev, const EventRecord *event)
void SetRefreshRate(int rate)
Summary information for an interaction.
Definition Interaction.h:56
static Interaction * NDecay(int tgt, int decay_mode=-1, int decayed_nucleon=0)
A utility class to facilitate creating the GENIE MC Ntuple from the output GENIE GHEP event records.
Definition NtpWriter.h:39
void Initialize(void)
add event
Definition NtpWriter.cxx:83
void Save(void)
get the even tree
void AddEventRecord(int ievent, const EventRecord *ev_rec)
save the event tree
Definition NtpWriter.cxx:57
void CustomizeFilenamePrefix(string prefix)
Singleton class to load & serve a TDatabasePDG.
Definition PDGLibrary.h:36
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
TRandom3 & RndEvg(void) const
rnd number generator used by the event generation drivers
Definition RandomGen.h:74
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
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
NtpMCFormat_t kDefOptNtpFormat
bool gOptUsingRootGeom
Long_t gOptRunNu
string kDefOptGeomDUnits
int gOptNev
map< int, double > gOptTgtMix
string gOptRootGeom
double gOptGeomDUnits
string gOptEvFilePrefix
string gOptRootGeomTopVol
HNLDecayMode_t gOptDecayMode
string lunits
string geom
int SelectInitState(void)
void GetCommandLineArgs(int argc, char **argv)
const EventRecordVisitorI * NucleonDecayGenerator(void)
void PrintSyntax(void)
int gOptDecayedNucleon
Utilities for improving the code readability when using PDG codes.
int IonPdgCodeToZ(int pdgc)
Definition PDGUtils.cxx:55
int IonPdgCodeToA(int pdgc)
Definition PDGUtils.cxx:63
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
bool IsValidMode(NucleonDecayMode_t ndm, int npdg=0)
string AsString(NucleonDecayMode_t ndm, int npdg=0)
int DecayedNucleonPdgCode(NucleonDecayMode_t ndm)
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f=' *')
vector< string > Split(string input, string delim)
double UnitFromString(string u)
Definition UnitUtils.cxx:18
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
bool gAbortingInErr
Definition Messenger.cxx:34
const int kPdgProton
Definition PDGCodes.h:81
const int kPdgNeutron
Definition PDGCodes.h:83
enum genie::ENucleonDecayMode NucleonDecayMode_t
@ kNFGHEP
Definition NtpMCFormat.h:30
enum genie::ENtpMCFormat NtpMCFormat_t