GENIEGenerator
Loading...
Searching...
No Matches
gBeamHNLEvGen.cxx
Go to the documentation of this file.
1//________________________________________________________________________________________
2/*!
3
4\program gevgen_hnl
5
6\brief A GENIE-based neutral heavy lepton event generation application.
7
8 *** Synopsis :
9
10 gevgen_hnl [-h]
11 [-r run#]
12 -n n_of_events
13 -f path/to/flux/files
14 [-E hnl_energy]
15 [--firstEvent first event for dk2nu flux readin]
16 [-m decay_mode]
17 [-g geometry (ROOT file)]
18 [-L geometry_length_units]
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_hnl 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 HNL decay mode ID:
37 -f
38 Input HNL flux.
39 --firstEvent
40 If using dk2nu fluxes, start reading at this entry
41 -g
42 Input detector geometry.
43 If a geometry is specified, HNL decay vertices will be distributed
44 in the desired detector volume.
45 Using this argument, you can pass a ROOT file containing your
46 detector geometry description.
47 -L
48 Input geometry length units, eg 'm', 'cm', 'mm', ...
49 [default: 'mm']
50 -o
51 Sets the prefix of the output event file.
52 The output filename is built as:
53 [prefix].[run_number].[event_tree_format].[file_format]
54 The default output filename is:
55 gntp.[run_number].ghep.root
56 This cmd line arguments lets you override 'gntp'
57 --seed
58 Random number seed.
59
60\author John Plows <komninos-john.plows \at cern.ch>
61 University of Oxford
62
63 Costas Andreopoulos <c.andreopoulos \at cern.ch>
64 University of Liverpool
65
66\created February 11, 2020
67
68\cpright Copyright (c) 2003-2025, The GENIE Collaboration
69 For the full text of the license visit http://copyright.genie-mc.org
70
71*/
72//_________________________________________________________________________________________
73
74#include <cassert>
75#include <cstdlib>
76#include <string>
77#include <vector>
78#include <sstream>
79
80#include <TSystem.h>
81
109
110using std::string;
111using std::vector;
112using std::ostringstream;
113
114using namespace genie;
115using namespace genie::hnl;
116using namespace genie::hnl::enums;
117
118#ifdef __GENIE_FLUX_DRIVERS_ENABLED__
119#define __CAN_GENERATE_EVENTS_USING_A_FLUX__
120//#include "Tools/Flux/GNuMIFlux.h"
121#include <TH1.h>
122#endif // #ifdef __GENIE_FLUX_DRIVERS_ENABLED__
123
124#ifdef __GENIE_GEOM_DRIVERS_ENABLED__
125#define __CAN_USE_ROOT_GEOM__
126#include <TGeoVolume.h>
127#include <TGeoManager.h>
128#include <TGeoShape.h>
129#include <TGeoBBox.h>
130#endif // #ifdef __GENIE_GEOM_DRIVERS_ENABLED__
131
132// function prototypes
133void GetCommandLineArgs (int argc, char ** argv);
134void PrintSyntax (void);
135
136int SelectDecayMode (std::vector<HNLDecayMode_t> *intChannels, SimpleHNL sh);
138
139#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
140void GenerateEventsUsingFlux (void);
141void FillFluxNonsense (FluxContainer &ggn);
142void FillFlux (FluxContainer &ggn, FluxContainer &tgn);
143#endif // #ifdef __GENIE_FLUX_DRIVERS_ENABLED__
144
145TLorentzVector GeneratePosition( GHepRecord * event );
146#ifdef __CAN_USE_ROOT_GEOM__
147void InitBoundingBox (void);
148#endif // #ifdef __CAN_USE_ROOT_GEOM__
149
150//
151string kDefOptGeomLUnits = "mm"; // default geometry length units
152string kDefOptGeomDUnits = "g_cm3"; // default geometry density units
153NtpMCFormat_t kDefOptNtpFormat = kNFGHEP; // default event tree format
154string kDefOptEvFilePrefix = "gntp";
155string kDefOptFluxFilePath = "./input-flux.root";
156
157string kDefOptSName = "genie::EventGenerator";
158string kDefOptSConfig = "BeamHNL";
159string kDefOptSTune = "GHNL20_00a_00_000";
160
161//
162Long_t gOptRunNu = 1000; // run number
163int gOptNev = 10; // number of events to generate
164
165double gOptEnergyHNL = -1; // HNL energy
166double gOptMassHNL = -1; // HNL mass
167double gOptECoupling = -1; // |U_e4|^2
168double gOptMCoupling = -1; // |U_m4|^2
169double gOptTCoupling = -1; // |U_t4|^2
170
171bool gOptIsMajorana = false; // is this Majorana?
172
173#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
174string gOptFluxFilePath = kDefOptFluxFilePath; // where flux files live
175map<string,string> gOptFluxShortNames;
176bool gOptIsUsingDk2nu = false; // using flat dk2nu files?
177int gOptFirstEvent = 0; // skip to this entry in dk2nu
178#endif // #ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
179bool gOptIsMonoEnFlux = true; // do we have monoenergetic flux?
180
182std::vector< HNLDecayMode_t > gOptIntChannels; // decays to un-inhibit
183
184string gOptEvFilePrefix = kDefOptEvFilePrefix; // event file prefix
185bool gOptUsingRootGeom = false; // using root geom or target mix?
186string gOptRootGeom; // input ROOT file with realistic detector geometry
187
188#ifdef __CAN_USE_ROOT_GEOM__
189TGeoManager * gOptRootGeoManager = 0; // the workhorse geometry manager
190TGeoVolume * gOptRootGeoVolume = 0;
191#endif // #ifdef __CAN_USE_ROOT_GEOM__
192
193string gOptRootGeomTopVol = ""; // input geometry top event generation volume
194double gOptGeomLUnits = 0; // input geometry length units
195long int gOptRanSeed = -1; // random number seed
196
197// Geometry bounding box and origin - read from the input geometry file (if any)
198double fdx = 0; // half-length - x
199double fdy = 0; // half-length - y
200double fdz = 0; // half-length - z
201double fox = 0; // origin - x
202double foy = 0; // origin - y
203double foz = 0; // origin - z
204
205double NTP_IS_E = 0., NTP_IS_PX = 0., NTP_IS_PY = 0., NTP_IS_PZ = 0.;
206double NTP_FS0_E = 0., NTP_FS0_PX = 0., NTP_FS0_PY = 0., NTP_FS0_PZ = 0.;
207double NTP_FS1_E = 0., NTP_FS1_PX = 0., NTP_FS1_PY = 0., NTP_FS1_PZ = 0.;
208double NTP_FS2_E = 0., NTP_FS2_PX = 0., NTP_FS2_PY = 0., NTP_FS2_PZ = 0.;
210
211// HNL lifetime in rest frame
212double CoMLifetime = -1.0; // GeV^{-1}
213// == Gamma( all valid channels ) / Gamma( all interesting channels )
214double decayMod = 1.0;
215// event weight
216double evWeight = 1.0;
217
218//_________________________________________________________________________________________
219int main(int argc, char ** argv)
220{
221 // suppress ROOT Info messages
222 gErrorIgnoreLevel = kWarning;
223
224 // Parse command line arguments
225 GetCommandLineArgs(argc,argv);
226
227 // Init messenger and random number seed
228 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
230
231 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
232
233 // Get the HNL generator first to load config
234 // config loaded upon instantiation of HNLGenerator algorithm
235 // calls LoadConfig() of each sub-alg
236 __attribute__((unused)) const EventRecordVisitorI * mcgen = HNLGenerator();
237 const Algorithm * algFluxCreator = AlgFactory::Instance()->GetAlgorithm("genie::hnl::FluxCreator", "Default");
238 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
239 const Algorithm * algVtxGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::VertexGenerator", "Default");
240
241 const FluxCreator * fluxCreator = dynamic_cast< const FluxCreator * >( algFluxCreator );
242 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
243 const VertexGenerator * vtxGen = dynamic_cast< const VertexGenerator * >( algVtxGen );
244
245 //string confString = kDefOptSName + "/" + kDefOptSConfig;
246 string confString = kDefOptSConfig;
247
248 CoMLifetime = hnlgen->GetHNLLifetime(); // GeV^{-1}
249
250 std::vector< double > U4l2s = hnlgen->GetHNLCouplings();
251 gOptECoupling = U4l2s.at(0);
252 gOptMCoupling = U4l2s.at(1);
253 gOptTCoupling = U4l2s.at(2);
254 gOptIsMajorana = hnlgen->IsHNLMajorana();
255
256 std::string stIntChannels = hnlgen->GetHNLInterestingChannels(); int iChan = -1;
257 if( gOptIntChannels.size() > 0 ) gOptIntChannels.clear();
258 while( stIntChannels.size() > 0 ){ // read channels from right (lowest mass) to left (highest mass)
259 iChan++;
260 HNLDecayMode_t md = static_cast< HNLDecayMode_t >( iChan );
261 std::string tmpSt = stIntChannels.substr( stIntChannels.size()-1, stIntChannels.size() );
262 if( std::strcmp( tmpSt.c_str(), "1" ) == 0 )
263 gOptIntChannels.emplace_back( md );
264
265 stIntChannels.erase( stIntChannels.end()-1, stIntChannels.end() );
266 }
267
268 //gOptIntChannels = confIntChan;
269
270 // Initialize an Ntuple Writer to save GHEP records into a TTree
273 ntpw.Initialize();
274
275 // add flux info to the tree
276 FluxContainer gnmf;
277 if( gOptIsUsingDk2nu ) {
278 // fill the flux object with nonsense to start with
279 FluxContainer * ptGnmf = new FluxContainer();
280 gnmf = *ptGnmf;
281 delete ptGnmf;
282 FillFluxNonsense( gnmf );
283 TBranch * flux = ntpw.EventTree()->Branch( "flux",
284 "genie::hnl::FluxContainer",
285 &gnmf, 32000, 1 );
286 flux->SetAutoDelete(kFALSE);
287 }
288
289 LOG("gevgen_hnl", pNOTICE)
290 << "Initialised Ntuple Writer";
291
292 // add another few branches to tree.
293 ntpw.EventTree()->Branch("hnl_mass", &gOptMassHNL, "gOptMassHNL/D");
294 ntpw.EventTree()->Branch("hnl_coup_e", &gOptECoupling, "gOptECoupling/D");
295 ntpw.EventTree()->Branch("hnl_coup_m", &gOptMCoupling, "gOptMCoupling/D");
296 ntpw.EventTree()->Branch("hnl_coup_t", &gOptTCoupling, "gOptTCoupling/D");
297 ntpw.EventTree()->Branch("hnl_ismaj", &gOptIsMajorana, "gOptIsMajorana/I");
298
299 // Create a MC job monitor for a periodically updated status file
300 GMCJMonitor mcjmonitor(gOptRunNu);
301 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
302
303 LOG("gevgen_hnl", pNOTICE)
304 << "Initialised MC job monitor";
305
306 // Set GHEP print level
307 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
308
309#ifdef __CAN_USE_ROOT_GEOM__
310 // Read geometry bounding box - for vertex position generation
311 if( gOptUsingRootGeom ){
312 InitBoundingBox();
313 }
314#endif // #ifdef __CAN_USE_ROOT_GEOM__
315
316 if( !gOptIsMonoEnFlux ){
317 LOG( "gevgen_hnl", pWARN )
318 << "Using input flux files. These are *flat dk2nu-like ROOT trees, so far...*";
319
320 }
321 // Event loop
322 int iflux = (gOptFirstEvent < 0) ? 0 : gOptFirstEvent; int ievent = iflux;
323 int maxFluxEntries = -1;
324 fluxCreator->SetInputFluxPath( gOptFluxFilePath );
325 fluxCreator->SetGeomFile( gOptRootGeom );
326 fluxCreator->SetFirstFluxEntry( iflux );
327 vtxGen->SetGeomFile( gOptRootGeom );
328
329 bool tooManyEntries = false;
330 while (1)
331 {
332 if( tooManyEntries ){
333 if( gOptNev >= 10000 ){
334 if( (ievent-gOptFirstEvent) % (gOptNev / 1000) == 0 ){
335 int irat = (iflux-gOptFirstEvent) / ( gOptNev / 1000 );
336 std::cerr << 0.1 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
337 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
338 }
339 } else if( gOptNev >= 100 ) {
340 if( (ievent-gOptFirstEvent) % (gOptNev / 10) == 0 ){
341 int irat = (iflux-gOptFirstEvent) / ( gOptNev / 10 );
342 std::cerr << 10.0 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
343 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
344 }
345 }
346 } else {
347 if( gOptNev >= 10000 ){
348 if( (ievent-gOptFirstEvent) % (gOptNev / 1000) == 0 ){
349 int irat = (ievent-gOptFirstEvent) / ( gOptNev / 1000 );
350 std::cerr << 0.1 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
351 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
352 }
353 } else if( gOptNev >= 100 ) {
354 if( (ievent-gOptFirstEvent) % (gOptNev / 10) == 0 ){
355 int irat = (ievent-gOptFirstEvent) / ( gOptNev / 10 );
356 std::cerr << 10.0 * irat << " % " << " ( " << (iflux-gOptFirstEvent)
357 << " seen ), ( " << (ievent-gOptFirstEvent) << " / " << gOptNev << " processed ) \r" << std::flush;
358 }
359 }
360 }
361
362 if( tooManyEntries && ((iflux-gOptFirstEvent) == gOptNev) ) break;
363 else if( (ievent-gOptFirstEvent) == gOptNev ) break;
364
365 if( ievent < gOptFirstEvent ){ ievent++; continue; }
366
367 assert( ievent >= gOptFirstEvent && gOptFirstEvent >= 0 );
368
369 LOG("gevgen_hnl", pNOTICE)
370 << " *** Generating event............ " << (ievent-gOptFirstEvent);
371
372 EventRecord * event = new EventRecord;
373 event->SetWeight(1.0);
374 event->SetProbability( CoMLifetime );
375 event->SetXSec( iflux ); // will be overridden, use as handy container
376 evWeight = 1.0;
377
378 if( !gOptIsMonoEnFlux ){
379 fluxCreator->ProcessEventRecord( event );
380
381 // fluxCreator->ProcessEventRecord now tells us how many entries there are
382 maxFluxEntries = fluxCreator->GetNFluxEntries();
383 if( gOptNev > maxFluxEntries - gOptFirstEvent ){
384 LOG( "gevgen_hnl", pWARN )
385 << "You have asked for " << gOptNev << " events, but only provided "
386 << maxFluxEntries - gOptFirstEvent << " flux entries. Truncating events to " << maxFluxEntries - gOptFirstEvent << ".";
387 gOptNev = maxFluxEntries - gOptFirstEvent;
388 tooManyEntries = true;
389 }
390 if( iflux >= maxFluxEntries - 1 ){
391 LOG( "gevgen_hnl", pWARN )
392 << "Reached end of flux input (iflux = " << iflux << "), about to stop.";
393 break;
394 }
395
396 FluxContainer retGnmf = fluxCreator->RetrieveFluxInfo();
397 FillFlux( gnmf, retGnmf );
398
399 // check to see if this was nonsense
400 if( ! event->Particle(0) ){ iflux++; delete event; continue; }
401
402 gOptEnergyHNL = event->Particle(0)->GetP4()->E();
403 iflux++;
404 } else { // monoenergetic HNL. Add it with energy and momentum pointing on z axis
405
406 assert( gOptEnergyHNL > gOptMassHNL );
407 double HNLP = std::sqrt( gOptEnergyHNL*gOptEnergyHNL - gOptMassHNL*gOptMassHNL );
408 TLorentzVector probeP4( 0.0, 0.0, HNLP, gOptEnergyHNL );
409 TLorentzVector v4( 0.0, 0.0, 0.0, 0.0 );
410 GHepParticle ptHNL( kPdgHNL, kIStInitialState, -1, -1, -1, -1, probeP4, v4 );
411 event->AddParticle( ptHNL );
412
413 }
414 assert( gOptEnergyHNL > gOptMassHNL );
415
416 int hpdg = genie::kPdgHNL;
417 int typeMod = 1;
418 if( gOptIsMajorana ) typeMod = ( rnd->RndGen().Uniform( 0.0, 1.0 ) >= 0.5 ) ? -1 : 1;
419 else if( event->Particle(0)->Pdg() > 0 ) typeMod = 1;
420 else if( event->Particle(0)->Pdg() < 0 ) typeMod = -1;
421
422 // int target = SelectInitState();
423 int decay = (int) gOptDecayMode;
424
425 assert( gOptECoupling >= 0.0 && gOptMCoupling >= 0.0 && gOptTCoupling >= 0.0 );
426
427 // RETHERE assuming all these HNL have K+- parent. This is wrong
428 // (but not very wrong for interesting masses)
429 SimpleHNL sh( "HNL", ievent, hpdg, genie::kPdgKP,
431
432 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
434
435 if( gOptDecayMode == kHNLDcyNull ){ // select from available modes
436 LOG("gevgen_hnl", pNOTICE)
437 << "No decay mode specified - sampling from all available modes.";
438
439 std::vector< HNLDecayMode_t > * intChannels = &gOptIntChannels;
440
441 decay = SelectDecayMode( intChannels, sh );
442 }
443
444 Interaction * interaction = Interaction::HNL(typeMod * genie::kPdgHNL, gOptEnergyHNL, decay);
445
446 if( event->Particle(0) ){ // we have an HNL with definite momentum, so let's set it now
447 interaction->InitStatePtr()->SetProbeP4( *(event->Particle(0)->P4()) );
448 interaction->InitStatePtr()->SetProbePdg( event->Particle(0)->Pdg() );
449 LOG( "gevgen_hnl", pDEBUG )
450 << "\nsetting probe p4 = " << utils::print::P4AsString( event->Particle(0)->P4() );
451 }
452
453 event->AttachSummary(interaction);
454
455 // Simulate decay
456 hnlgen->ProcessEventRecord(event);
457
458 // add the FS 4-momenta to special branches
459 // Quite inelegant. Gets the job done, though
460 NTP_FS0_PDG = (event->Particle(1))->Pdg();
461 NTP_FS0_E = ((event->Particle(1))->P4())->E();
462 NTP_FS0_PX = ((event->Particle(1))->P4())->Px();
463 NTP_FS0_PY = ((event->Particle(1))->P4())->Py();
464 NTP_FS0_PZ = ((event->Particle(1))->P4())->Pz();
465 NTP_FS1_PDG = (event->Particle(2))->Pdg();
466 NTP_FS1_E = ((event->Particle(2))->P4())->E();
467 NTP_FS1_PX = ((event->Particle(2))->P4())->Px();
468 NTP_FS1_PY = ((event->Particle(2))->P4())->Py();
469 NTP_FS1_PZ = ((event->Particle(2))->P4())->Pz();
470 if( event->Particle(3) ){
471 NTP_FS2_PDG = (event->Particle(3))->Pdg();
472 NTP_FS2_E = ((event->Particle(3))->P4())->E();
473 NTP_FS2_PX = ((event->Particle(3))->P4())->Px();
474 NTP_FS2_PY = ((event->Particle(3))->P4())->Py();
475 NTP_FS2_PZ = ((event->Particle(3))->P4())->Pz();
476 }
477 else{
478 NTP_FS2_PDG = 0;
479 NTP_FS2_E = 0.0;
480 NTP_FS2_PX = 0.0;
481 NTP_FS2_PY = 0.0;
482 NTP_FS2_PZ = 0.0;
483 }
484
485 // Generate (or read) a position for the decay vertex
486 // also currently handles the geometrical weight
487 TLorentzVector x4mm;
488 if( gOptUsingRootGeom ){
489 TLorentzVector * p4HNL = event->Particle(0)->GetP4();
490 NTP_IS_E = p4HNL->E(); NTP_IS_PX = p4HNL->Px(); NTP_IS_PY = p4HNL->Py(); NTP_IS_PZ = p4HNL->Pz();
491 vtxGen->ProcessEventRecord( event );
492 x4mm = *(event->Vertex());
493 } else {
494 x4mm = GeneratePosition( event );
495 }
496
497 // update weight to scale for couplings, inhibited decays
498 // acceptance is already handled in FluxCreator
499 // geometry handled in VertexGenerator
500 evWeight = event->Weight();
502 evWeight *= 1.0 / decayMod;
503 event->SetWeight( evWeight / 1.0e+20 ); // in units of 1e+20 POT
504
505 // store the decayMod in the second daughter's polarisation
506 event->Particle(2)->SetPolarization( 1.0, decayMod );
507
508 // why does InitState show the wrong p4 here?
509 interaction->InitStatePtr()->SetProbeP4( *(event->Particle(0)->P4()) );
510
511 LOG("gevgen_hnl", pDEBUG) << "Weight = " << evWeight;
512
513 LOG("gevgen_hnl", pINFO)
514 << "Generated event: " << *event;
515
516 // Add event at the output ntuple, refresh the mc job monitor & clean-up
517 ntpw.AddEventRecord(ievent, event);
518 mcjmonitor.Update(ievent,event);
519
520 delete event;
521
522 ievent++;
523 } // event loop
524
525 // Save the generated event tree & close the output file
526 ntpw.Save();
527
528 LOG("gevgen_hnl", pNOTICE) << "Done!";
529
530 return 0;
531}
532//_________________________________________________________________________________________
533//............................................................................
534#ifdef __CAN_USE_ROOT_GEOM__
535void InitBoundingBox(void)
536{
537// Initialise geometry bounding box, used for generating HNL vertex positions
538
539 LOG("gevgen_hnl", pINFO)
540 << "Initialising geometry bounding box.";
541
542 fdx = 0; // half-length - x
543 fdy = 0; // half-length - y
544 fdz = 0; // half-length - z
545 fox = 0; // origin - x
546 foy = 0; // origin - y
547 foz = 0; // origin - z
548
549 if(!gOptUsingRootGeom){ // make a unit-m sided box
550 LOG("gevgen_hnl", pINFO)
551 << "No geometry file input detected, making a unit-m side box volume.";
552
553 TGeoManager * geom = new TGeoManager( "box1", "A simple box detector" );
554
555 //--- define some materials
556 TGeoMaterial *matVacuum = new TGeoMaterial("Vacuum", 0,0,0);
557 TGeoMaterial *matAl = new TGeoMaterial("Al", 26.98,13,2.7);
558 //--- define some media
559 TGeoMedium *Vacuum = new TGeoMedium("Vacuum",1, matVacuum);
560 TGeoMedium *Al = new TGeoMedium("Root Material",2, matAl);
561
562 //--- make the top container volume
563 //const double boxSideX = 2.5, boxSideY = 2.5, boxSideZ = 2.5; // m
564 //const double bigBoxSide = 2.0 * std::max( boxSideX, std::max( boxSideY, boxSideZ ) ); // m
565 //const double worldLen = 1.01 * bigBoxSide; // m
566
567 TGeoVolume * topvol = geom->MakeBox( "TOP", Vacuum, 101.0, 101.0, 101.0 );
568 geom->SetTopVolume( topvol );
569
570 //--- make the detector box container
571 TGeoVolume * boxvol = geom->MakeBox( "VOL", Vacuum, 100.5, 100.5, 100.5 );
572 boxvol->SetVisibility(kFALSE);
573
574 //--- origin is at centre of the box
575 TGeoVolume * box = geom->MakeBox( "BOX", Al, 100.0, 100.0, 100.0 );
576 //TGeoTranslation * tr0 = new TGeoTranslation( 0.0, 0.0, 0.0 );
577 TGeoRotation * rot0 = new TGeoRotation( "rot0", 90.0, 0.0, 90.0, 90.0, 0.0, 0.0 );
578
579 //--- add directly to top volume
580 topvol->AddNode( box, 1, rot0 );
581
582 gOptRootGeoManager = geom;
583
584 return;
585 }
586
587 bool geom_is_accessible = ! (gSystem->AccessPathName(gOptRootGeom.c_str()));
588 if (!geom_is_accessible) {
589 LOG("gevgen_hnl", pFATAL)
590 << "The specified ROOT geometry doesn't exist! Initialization failed!";
591 exit(1);
592 }
593
594 if( !gOptRootGeoManager ) gOptRootGeoManager = TGeoManager::Import(gOptRootGeom.c_str());
595
596 TGeoVolume * top_volume = gOptRootGeoManager->GetTopVolume();
597 assert( top_volume );
598 TGeoShape * ts = top_volume->GetShape();
599
600 TGeoBBox * box = (TGeoBBox *)ts;
601
602 //get box origin and dimensions (in the same units as the geometry)
603 fdx = box->GetDX();
604 fdy = box->GetDY();
605 fdz = box->GetDZ();
606 fox = (box->GetOrigin())[0];
607 foy = (box->GetOrigin())[1];
608 foz = (box->GetOrigin())[2];
609
610 LOG("gevgen_hnl", pDEBUG)
611 << "Before conversion the bounding box has:"
612 << "\nOrigin = ( " << fox << " , " << foy << " , " << foz << " )"
613 << "\nDimensions = " << fdx << " x " << fdy << " x " << fdz
614 << "\n1cm = 1.0 unit";
615
616 // Convert from local to SI units
623
624 LOG("gevgen_hnl", pINFO)
625 << "Initialised bounding box successfully.";
626
627}
628#endif // #ifdef __CAN_USE_ROOT_GEOM__
629//............................................................................
630//_________________________________________________________________________________________
631//............................................................................
632#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
633void FillFluxNonsense( FluxContainer &ggn )
634{
635 ggn.evtno = -9999;
636
637 ggn.pdg = -9999;
638 ggn.parPdg = -9999;
639 ggn.lepPdg = -9999;
640 ggn.nuPdg = -9999;
641
642 ggn.prodChan = -9999;
643 ggn.nuProdChan = -9999;
644
645 ggn.startPoint.SetXYZ(-9999.9, -9999.9, -9999.9);
646 ggn.targetPoint.SetXYZ(-9999.9, -9999.9, -9999.9);
647 ggn.startPointUser.SetXYZ(-9999.9, -9999.9, -9999.9);
648 ggn.targetPointUser.SetXYZ(-9999.9, -9999.9, -9999.9);
649 ggn.delay = -9999.9;
650
651 ggn.polz.SetXYZ(-9999.9, -9999.9, -9999.9);
652
653 ggn.p4.SetPxPyPzE(-9999.9, -9999.9, -9999.9, -9999.9);
654 ggn.parp4.SetPxPyPzE(-9999.9, -9999.9, -9999.9, -9999.9);
655 ggn.p4User.SetPxPyPzE(-9999.9, -9999.9, -9999.9, -9999.9);
656 ggn.parp4User.SetPxPyPzE(-9999.9, -9999.9, -9999.9, -9999.9);
657
658 ggn.Ecm = -9999.9;
659 ggn.nuEcm = -9999.9;
660 ggn.XYWgt = -9999.9;
661 ggn.boostCorr = -9999.9;
662 ggn.accCorr = -9999.9;
663 ggn.zetaMinus = -9999.9;
664 ggn.zetaPlus = -9999.9;
665 ggn.acceptance = -9999.9;
666
667 ggn.nimpwt = -9999.9;
668
669 return;
670}
671//_________________________________________________________________________________________
672void FillFlux( FluxContainer &ggn, FluxContainer &tgn )
673{
674 ggn.evtno = tgn.evtno;
675
676 ggn.pdg = tgn.pdg;
677 ggn.parPdg = tgn.parPdg;
678 ggn.lepPdg = tgn.lepPdg;
679 ggn.nuPdg = tgn.nuPdg;
680
681 ggn.prodChan = tgn.prodChan;
682 ggn.nuProdChan = tgn.nuProdChan;
683
684 ggn.startPoint = tgn.startPoint;
685 ggn.targetPoint = tgn.targetPoint;
688 ggn.delay = tgn.delay;
689
690 ggn.polz = tgn.polz;
691
692 ggn.p4 = tgn.p4;
693 ggn.parp4 = tgn.parp4;
694 ggn.p4User = tgn.p4User;
695 ggn.parp4User = tgn.parp4User;
696
697 ggn.Ecm = tgn.Ecm;
698 ggn.nuEcm = tgn.nuEcm;
699 ggn.XYWgt = tgn.XYWgt;
700 ggn.boostCorr = tgn.boostCorr;
701 ggn.accCorr = tgn.accCorr;
702 ggn.zetaMinus = tgn.zetaMinus;
703 ggn.zetaPlus = tgn.zetaPlus;
704 ggn.acceptance = tgn.acceptance;
705
706 ggn.nimpwt = tgn.nimpwt;
707}
708//............................................................................
709#endif // #ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
710//_________________________________________________________________________________________
711TLorentzVector GeneratePosition( GHepRecord * event )
712{
713 __attribute__((unused)) RandomGen * rnd = RandomGen::Instance();
714 TRandom3 & rnd_geo = rnd->RndGeom();
715
716 double rndx = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
717 double rndy = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
718 double rndz = 2 * (rnd_geo.Rndm() - 0.5); // [-1,1]
719
720 double t_gen = 0;
721 double x_gen = fox + rndx * fdx;
722 double y_gen = foy + rndy * fdy;
723 double z_gen = foz + rndz * fdz;
724
725 TLorentzVector pos(x_gen, y_gen, z_gen, t_gen);
726 return pos;
727}
728//_________________________________________________________________________________________
730{
731 //string sname = "genie::EventGenerator";
732 //string sconfig = "BeamHNL";
734
735 LOG("gevgen_hnl", pINFO)
736 << "Instantiating HNL generator.";
737
738 const Algorithm * algmcgen = algf->GetAlgorithm(kDefOptSName, kDefOptSConfig);
739 LOG("gevgen_hnl", pDEBUG)
740 << "Got algorithm " << kDefOptSName.c_str() << "/" << kDefOptSConfig.c_str();;
741
742 const EventRecordVisitorI * mcgen =
743 dynamic_cast< const EventRecordVisitorI * >( algmcgen );
744 if(!mcgen) {
745 LOG("gevgen_hnl", pFATAL) << "Couldn't instantiate the HNL generator";
746 gAbortingInErr = true;
747 exit(1);
748 }
749
750 LOG("gevgen_hnl", pINFO)
751 << "HNL generator instantiated successfully.";
752
753 return mcgen;
754}
755//_________________________________________________________________________________________
756int SelectDecayMode( std::vector< HNLDecayMode_t > * intChannels, SimpleHNL sh )
757{
758 const std::map< HNLDecayMode_t, double > gammaMap = sh.GetValidChannels();
759
760 // set CoM lifetime now if unset
761 if( CoMLifetime < 0.0 ){
763 }
764
765 std::vector< HNLDecayMode_t > intAndValidChannels;
766 for( std::vector< HNLDecayMode_t >::iterator it = intChannels->begin(); it != intChannels->end(); ++it ){
767 HNLDecayMode_t mode = *it;
768
769 //check if this is a valid mode
770 if( !utils::hnl::IsKinematicallyAllowed( mode, gOptMassHNL ) ) continue;
771
772 intAndValidChannels.emplace_back( mode );
773 auto mapG = gammaMap.find( mode );
774 double theGamma = mapG->second;
775 LOG("gevgen_hnl", pDEBUG)
776 << "For mode " << utils::hnl::AsString( mode ) << " gamma = " << theGamma;
777 }
778
779 if( intAndValidChannels.size() == 0 ){ // all the modes picked by user are too heavy. Abort.
780 LOG( "gevgen_hnl", pFATAL )
781 << "None of the channels specified as interesting are kinematically allowed. Please either increase the HNL mass or change interesting channels in config.";
782 exit(1);
783 }
784
785 std::map< HNLDecayMode_t, double > intMap =
786 selector::SetInterestingChannels( intAndValidChannels, gammaMap );
787
788 sh.SetInterestingChannels( intMap );
789
790 // update fraction of total decay width that is not in inhibited channels
791 double gammaAll = 0.0, gammaInt = 0.0;
792 for( std::map< HNLDecayMode_t, double >::const_iterator itall = gammaMap.begin() ;
793 itall != gammaMap.end() ; ++itall ){
794 gammaAll += (*itall).second;
795 }
796 for( std::map< HNLDecayMode_t, double >::iterator itint = intMap.begin() ;
797 itint != intMap.end() ; ++itint ){
798 gammaInt += (*itint).second;
799 }
800 assert( gammaInt > 0.0 && gammaAll >= gammaInt );
801 decayMod = gammaInt / gammaAll;
802
803 // get probability that channels in intAndValidChannels will be selected
804 std::map< HNLDecayMode_t, double > PMap =
806
807 // now do a random throw
809 double ranThrow = rnd->RndGen().Uniform( 0., 1. ); // HNL's fate is sealed.
810
811 HNLDecayMode_t selectedDecayChan =
812 selector::SelectChannelInclusive( PMap, ranThrow );
813
814 int decay = ( int ) selectedDecayChan;
815 return decay;
816}
817//_________________________________________________________________________________________
818void GetCommandLineArgs(int argc, char ** argv)
819{
820 LOG("gevgen_hnl", pINFO) << "Parsing command line arguments";
821
822 // Parse run options for this app
823
824 CmdLnArgParser parser(argc,argv);
825
826 // force the app to look at HNL tune by default
827 // if user passes --tune argument, look at the user input tune instead
828 char * expargv[ argc + 2 ];
829 bool didFindUserInputTune = false;
830 std::string stExtraTuneBit = kDefOptSTune;
831
832 if( parser.OptionExists("tune") ){
833 didFindUserInputTune = true;
834 stExtraTuneBit = parser.ArgAsString("tune");
835 LOG( "gevgen_hnl", pWARN )
836 << "Using input HNL tune " << parser.ArgAsString("tune");
837 } else {
838 LOG( "gevgen_hnl", pWARN )
839 << "Using default HNL tune " << kDefOptSTune;
840 }
841 // append this to argv
842 for( int iArg = 0; iArg < argc; iArg++ ){
843 expargv[iArg] = argv[iArg];
844 }
845 if( !didFindUserInputTune ){
846 char * chBit = const_cast< char * >( stExtraTuneBit.c_str() ); // ugh. Ugly.
847 std::string stune("--tune"); char * tBit = const_cast< char * >( stune.c_str() );
848 expargv[argc] = tBit;
849 expargv[argc+1] = chBit;
850 }
851
852 // Common run options.
853 int expargc = ( didFindUserInputTune ) ? argc : argc+2;
854 std::string stnull(""); char * nBit = const_cast< char * >( stnull.c_str() );
855 expargv[expargc] = nBit;
856
857 RunOpt::Instance()->ReadFromCommandLine(expargc,expargv);
858
859 // help?
860 bool help = parser.OptionExists('h');
861 if(help) {
862 PrintSyntax();
863 exit(0);
864 }
865
866 // run number
867 if( parser.OptionExists('r') ) {
868 LOG("gevgen_hnl", pDEBUG) << "Reading MC run number";
869 gOptRunNu = parser.ArgAsLong('r');
870 } else {
871 LOG("gevgen_hnl", pDEBUG) << "Unspecified run number - Using default";
872 gOptRunNu = 1000;
873 } //-r
874
875 // number of events
876 if( parser.OptionExists('n') ) {
877 LOG("gevgen_hnl", pDEBUG)
878 << "Reading number of events to generate";
879 gOptNev = parser.ArgAsInt('n');
880 } else {
881 LOG("gevgen_hnl", pFATAL)
882 << "You need to specify the number of events";
883 PrintSyntax();
884 exit(0);
885 } //-n
886
887 // get HNL mass
888 const Algorithm * algHNLGen = AlgFactory::Instance()->GetAlgorithm("genie::hnl::Decayer", "Default");
889 const Decayer * hnlgen = dynamic_cast< const Decayer * >( algHNLGen );
890 gOptMassHNL = hnlgen->GetHNLMass();
891
892 bool isMonoEnergeticFlux = true;
893#ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
894 if( parser.OptionExists('f') ) {
895 LOG("gevgen_hnl", pDEBUG)
896 << "A flux has been offered. Searching this path: " << parser.ArgAsString('f');
897 isMonoEnergeticFlux = false;
898 gOptFluxFilePath = parser.ArgAsString('f');
899
900 // check if this is valid path (assume these are dk2nu files)
901 if( gSystem->OpenDirectory( gOptFluxFilePath.c_str() ) != NULL ){
902 gOptIsUsingDk2nu = true;
903 LOG("gevgen_hnl", pDEBUG)
904 << "dk2nu flux files detected. Will create flux spectrum dynamically.";
905 } else {
906 LOG("gevgen_hnl", pFATAL)
907 << "Invalid flux file path " << gOptFluxFilePath;
908 exit(1);
909 }
910 } else {
911 // we need the 'E' option! Log it and pass below
912 LOG("gevgen_hnl", pINFO)
913 << "No flux file offered. Assuming monoenergetic flux.";
914 } //-f
915#endif // #ifdef __CAN_GENERATE_EVENTS_USING_A_FLUX__
916
917 // HNL energy (only relevant if we do not have an input flux)
918 gOptEnergyHNL = -1;
919 if( isMonoEnergeticFlux ){
920 if( parser.OptionExists('E') ) {
921 LOG("gevgen_hnl", pDEBUG)
922 << "Reading HNL energy";
923 gOptEnergyHNL = parser.ArgAsDouble('E');
924 } else {
925 LOG("gevgen_hnl", pFATAL)
926 << "You need to specify the HNL energy";
927 PrintSyntax();
928 exit(0);
929 } //-E
930 assert(gOptEnergyHNL > gOptMassHNL);
931 }
932
933 gOptIsMonoEnFlux = isMonoEnergeticFlux;
934
935 // first flux entry to read
936 if( parser.OptionExists("firstEvent") ) {
937 gOptFirstEvent = parser.ArgAsInt("firstEvent");
938 LOG( "gevgen_hnl", pINFO )
939 << "Starting flux readin at first event = " << gOptFirstEvent;
940 } // --firstEvent
941
942 // HNL decay mode
943 int mode = -1;
944 if( parser.OptionExists('m') ) {
945 LOG("gevgen_hnl", pDEBUG)
946 << "Reading HNL decay mode";
947 mode = parser.ArgAsInt('m');
948 } else {
949 LOG("gevgen_hnl", pINFO)
950 << "No decay mode specified - will sample from allowed decay modes";
951 } //-m
953
955 if(!allowed) {
956 LOG("gevgen_hnl", pFATAL)
957 << "Specified decay is not allowed kinematically for the given HNL mass";
958 PrintSyntax();
959 exit(0);
960 }
961
962 //
963 // geometry
964 //
965
966 string geom = "";
967 string lunits;
968#ifdef __CAN_USE_ROOT_GEOM__
969 // string dunits;
970 if( parser.OptionExists('g') ) {
971 LOG("gevgen_hnl", pDEBUG) << "Getting input geometry";
972 geom = parser.ArgAsString('g');
973
974 // is it a ROOT file that contains a ROOT geometry?
975 bool accessible_geom_file =
976 ! (gSystem->AccessPathName(geom.c_str()));
977 if (accessible_geom_file) {
979 gOptUsingRootGeom = true;
980 } else {
981 LOG("gevgen_hnl", pFATAL)
982 << "Geometry option is not a ROOT file. Please use ROOT geom.";
983 PrintSyntax();
984 exit(1);
985 }
986 } else {
987 // LOG("gevgen_hnl", pFATAL)
988 // << "No geometry option specified - Exiting";
989 // PrintSyntax();
990 // exit(1);
991 } //-g
992
994 // using a ROOT geometry - get requested geometry units
995
996 // length units:
997 if( parser.OptionExists('L') ) {
998 LOG("gevgen_hnl", pDEBUG)
999 << "Checking for input geometry length units";
1000 lunits = parser.ArgAsString('L');
1001 } else {
1002 LOG("gevgen_hnl", pDEBUG) << "Using default geometry length units";
1004 } // -L
1005 // // density units:
1006 // if( parser.OptionExists('D') ) {
1007 // LOG("gevgen_hnl", pDEBUG)
1008 // << "Checking for input geometry density units";
1009 // dunits = parser.ArgAsString('D');
1010 // } else {
1011 // LOG("gevgen_hnl", pDEBUG) << "Using default geometry density units";
1012 // dunits = kDefOptGeomDUnits;
1013 // } // -D
1015 // gOptGeomDUnits = utils::units::UnitFromString(dunits);
1016
1017 } // using root geom?
1018#endif // #ifdef __CAN_USE_ROOT_GEOM__
1019
1020 // event file prefix
1021 if( parser.OptionExists('o') ) {
1022 LOG("gevgen_hnl", pDEBUG) << "Reading the event filename prefix";
1023 gOptEvFilePrefix = parser.ArgAsString('o');
1024 } else {
1025 LOG("gevgen_hnl", pDEBUG)
1026 << "Will set the default event filename prefix";
1028 } //-o
1029
1030 // random number seed
1031 if( parser.OptionExists("seed") ) {
1032 LOG("gevgen_hnl", pINFO) << "Reading random number seed";
1033 gOptRanSeed = parser.ArgAsLong("seed");
1034 } else {
1035 LOG("gevgen_hnl", pINFO) << "Unspecified random number seed - Using default";
1036 gOptRanSeed = -1;
1037 }
1038
1039 //
1040 // >>> print the command line options
1041 //
1042
1043 ostringstream gminfo;
1044 if (gOptUsingRootGeom) {
1045 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
1046 << ", top volume: "
1047 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
1048 << ", length units: " << lunits;
1049 // << ", density units: " << dunits;
1050 }
1051
1052 LOG("gevgen_hnl", pNOTICE)
1053 << "\n\n"
1054 << utils::print::PrintFramedMesg("gevgen_hnl job configuration");
1055
1056 LOG("gevgen_hnl", pNOTICE)
1057 << "\n @@ Run number : " << gOptRunNu
1058 << "\n @@ Random seed : " << gOptRanSeed
1059 << "\n @@ HNL mass : " << gOptMassHNL << " GeV"
1060 << "\n @@ Decay channel : " << utils::hnl::AsString(gOptDecayMode)
1061 << "\n @@ Geometry : " << gminfo.str()
1062 << "\n @@ Statistics : " << gOptNev << " events";
1063}
1064//_________________________________________________________________________________________
1065void PrintSyntax(void)
1066{
1067 LOG("gevgen_hnl", pFATAL)
1068 << "\n **Syntax**"
1069 << "\n gevgen_hnl [-h] "
1070 << "\n [-r run#]"
1071 << "\n -n n_of_events"
1072 << "\n -f path/to/flux/files"
1073 << "\n [-E hnl_energy]"
1074 << "\n [--firstEvent first_event_for_dk2nu_readin]"
1075 << "\n [-m decay_mode]"
1076 << "\n [-g geometry (ROOT file)]"
1077 << "\n [-L length_units_at_geom]"
1078 << "\n [-o output_event_file_prefix]"
1079 << "\n [--seed random_number_seed]"
1080 << "\n [--message-thresholds xml_file]"
1081 << "\n [--event-record-print-level level]"
1082 << "\n [--mc-job-status-refresh-rate rate]"
1083 << "\n"
1084 << " Please also read the detailed documentation at http://www.genie-mc.org"
1085 << " or look at the source code: $GENIE/src/Apps/gBeamHNLEvGen.cxx"
1086 << "\n";
1087}
1088//_________________________________________________________________________________________
#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()
Algorithm abstract base class.
Definition Algorithm.h:54
Command line argument parser.
string ArgAsString(char opt)
double ArgAsDouble(char opt)
bool OptionExists(char opt)
was option set?
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the 'Visito...
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.
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
virtual void SetWeight(double wght)
Definition GHepRecord.h:130
static void SetPrintLevel(int print_level)
Simple class to create & update MC job status files and env. vars. This is used to be able to keep tr...
Definition GMCJMonitor.h:31
void Update(int iev, const EventRecord *event)
void SetRefreshRate(int rate)
void SetProbePdg(int pdg_code)
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition Interaction.h:56
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
static Interaction * HNL(int probe, double E=0, int decayed_mode=-1)
A utility class to facilitate creating the GENIE MC Ntuple from the output GENIE GHEP event records.
Definition NtpWriter.h:39
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)
TTree * EventTree(void)
Definition NtpWriter.h:55
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
TRandom3 & RndGeom(void) const
rnd number generator used by geometry drivers
Definition RandomGen.h:68
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
TRandom3 & RndGen(void) const
rnd number generator for generic usage
Definition RandomGen.h:80
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
Heavy Neutral Lepton final-state product generator.
Definition HNLDecayer.h:41
std::string GetHNLInterestingChannels() const
double GetHNLMass() const
void ProcessEventRecord(GHepRecord *event) const
double GetHNLLifetime() const
bool IsHNLMajorana() const
std::vector< double > GetHNLCouplings() const
A GENIE flux container specific for HNL containers. Based on the dk2nu flux paradigm and genie::flux:...
double nimpwt
Weight of parent.
double zetaPlus
maximum angular deviation from parent momentum to reach detector [deg]
double Ecm
Parent rest-frame energy of HNL [GeV].
double delay
delay HNL would have wrt SMv [ns]
int nuPdg
PDG code of SM neutrino that would have been produced.
int prodChan
Decay mode that produced HNL.
TVector3 polz
HNL polarisation vector, in HNL rest frame, in NEAR coords.
TLorentzVector parp4
parent momentum at HNL production in NEAR coords [GeV/c]
double zetaMinus
minimum angular deviation from parent momentum to reach detector [deg]
TVector3 targetPoint
point in detector HNL is forced towards in NEAR coords [m]
int nuProdChan
Decay mode that would have produced SM neutrino.
int parPdg
parent PDG code
TLorentzVector parp4User
parent momentum at HNL production in USER coords [GeV/c]
TLorentzVector p4User
HNL momentum in USER coords [GeV/c].
double accCorr
acceptance correction (collimation effect. SM v == 1)
TVector3 startPointUser
parent decay vertex in USER coords [m]
double nuEcm
Parent rest-frame energy of equivalent SM neutrino [GeV].
double XYWgt
geometric acceptance (angular size of detector in parent rest frame)
double boostCorr
boost correction wrt parent rest-frame (ELAB = ECM * boostCorr)
TVector3 startPoint
parent decay vertex in NEAR coords [m]
TLorentzVector p4
HNL momentum in NEAR coords [GeV/c].
TVector3 targetPointUser
point in detector HNL is forced towards in USER coords [m]
double acceptance
full acceptance == XYWgt * boostCorr^2 * accCorr
int lepPdg
PDG code of lepton produced with HNL on parent decay.
Calculates HNL production kinematics & production vertex. Is a concrete implementation of the FluxRec...
void SetFirstFluxEntry(int iFirst) const
void ProcessEventRecord(GHepRecord *event_rec) const
FluxContainer RetrieveFluxInfo() const
void SetGeomFile(string geomfile) const
void SetInputFluxPath(std::string finpath) const
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
Heavy Neutral Lepton decay vertex generator given production vertex and momentum ***.
void SetGeomFile(std::string geomfile) const
void ProcessEventRecord(GHepRecord *event_rec) const
long int gOptRanSeed
double gOptGeomLUnits
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
NtpMCFormat_t kDefOptNtpFormat
bool gOptUsingRootGeom
Long_t gOptRunNu
string kDefOptGeomDUnits
int gOptNev
string gOptRootGeom
string gOptEvFilePrefix
string gOptRootGeomTopVol
double gOptECoupling
TLorentzVector GeneratePosition(GHepRecord *event)
bool gOptIsMajorana
double NTP_FS1_E
double NTP_IS_E
double NTP_FS1_PX
double NTP_FS0_PZ
double foy
double gOptMassHNL
HNLDecayMode_t gOptDecayMode
double NTP_FS1_PY
double gOptMCoupling
double NTP_FS0_PY
double fdz
std::vector< HNLDecayMode_t > gOptIntChannels
double decayMod
double foz
const EventRecordVisitorI * HNLGenerator(void)
double NTP_FS0_PX
double evWeight
double fox
bool gOptIsMonoEnFlux
double gOptEnergyHNL
double NTP_IS_PX
double NTP_FS0_E
double fdy
double NTP_FS2_PZ
int SelectDecayMode(std::vector< HNLDecayMode_t > *intChannels, SimpleHNL sh)
double NTP_FS2_PY
int NTP_FS1_PDG
double fdx
int NTP_FS2_PDG
string kDefOptFluxFilePath
void GetCommandLineArgs(int argc, char **argv)
double NTP_FS2_PX
void PrintSyntax(void)
double NTP_IS_PY
string kDefOptSTune
double CoMLifetime
double NTP_FS2_E
int NTP_FS0_PDG
string kDefOptSName
double gOptTCoupling
double NTP_FS1_PZ
double NTP_IS_PZ
string kDefOptSConfig
string lunits
string geom
map< string, string > gOptFluxShortNames
GENIE flux drivers.
Typedef enums.
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)
enum genie::hnl::EHNLDecayMode HNLDecayMode_t
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
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=' *')
string P4AsString(const TLorentzVector *p)
double UnitFromString(string u)
Definition UnitUtils.cxx:18
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
@ kIStInitialState
Definition GHepStatus.h:29
bool gAbortingInErr
Definition Messenger.cxx:34
const int kPdgKP
Definition PDGCodes.h:172
@ kNFGHEP
Definition NtpMCFormat.h:30
const int kPdgHNL
Definition PDGCodes.h:224
enum genie::ENtpMCFormat NtpMCFormat_t