GENIEGenerator
Loading...
Searching...
No Matches
gUpMuFluxGen.cxx File Reference
#include <cassert>
#include <cstdlib>
#include <cctype>
#include <string>
#include <vector>
#include <sstream>
#include <map>
#include <TRotation.h>
#include <TH1D.h>
#include <TH3D.h>
#include <TTree.h>
#include <TLorentzVector.h>
#include "Framework/Algorithm/Algorithm.h"
#include "Framework/Algorithm/AlgFactory.h"
#include "Framework/EventGen/XSecAlgorithmI.h"
#include "Framework/Conventions/Constants.h"
#include "Framework/Conventions/Units.h"
#include "Framework/EventGen/GFluxI.h"
#include "Framework/EventGen/GMCJDriver.h"
#include "Framework/EventGen/GMCJMonitor.h"
#include "Framework/Numerical/RandomGen.h"
#include "Framework/Messenger/Messenger.h"
#include "Framework/ParticleData/PDGCodes.h"
#include "Framework/ParticleData/PDGLibrary.h"
#include "Framework/Utils/XSecSplineList.h"
#include "Framework/Utils/StringUtils.h"
#include "Framework/Utils/SystemUtils.h"
#include "Framework/Utils/PrintUtils.h"
#include "Framework/Utils/UnitUtils.h"
#include "Framework/Utils/KineUtils.h"
#include "Framework/Utils/AppInit.h"
#include "Framework/Utils/RunOpt.h"
#include "Framework/Utils/CmdLnArgParser.h"
Include dependency graph for gUpMuFluxGen.cxx:

Go to the source code of this file.

Functions

void GetCommandLineArgs (int argc, char **argv)
void PrintSyntax (void)
GFluxIGetFlux (void)
void GenerateUpNu (GFluxI *flux_driver)
TH3D * BuildEmuEnuCosThetaPdf (int nu_code)
TH1D * GetEmuPdf (double Enu, double costheta, const TH3D *pdf3d)
TVector3 GetDetectorVertex (double CosTheta, double Enu)
double GetCrossSection (int nu_code, double Enu, double Emu)
double ProbabilityEmu (int nu_code, double Enu, double Emu)
int main (int argc, char **argv)

Variables

Long_t gOptRunNu
string gOptFluxSim
map< int, string > gOptFluxFiles
int gOptNev = -1
double gOptDetectorSide
long int gOptRanSeed
string gOptInpXSecFile
const double a = 2e+6
const double e = 500e+9
double kDefOptDetectorSide = 1e+5

Function Documentation

◆ BuildEmuEnuCosThetaPdf()

TH3D * BuildEmuEnuCosThetaPdf ( int nu_code)

Definition at line 444 of file gUpMuFluxGen.cxx.

445{
446// Set up a 3D histogram, with axes Emu, Enu, CosTheta.
447// Bin convention is defined at the start.
448// Content of each bin is given by the probability of getting a muon
449// of energy Emu from a neutrino of energy Enu and zenith ange costheta
450
451 // Bin convention for Enu, consistent with BGLRS
452 const double Enumin = 0.1;
453 const int nEnubinsPerDecade = 10;
454 const int nEnubins = 31;
455
456 // Bin convention for Emu, consistent with BGLRS
457 const double Emumin = 0.1;
458 const int nEmubinsPerDecade = 10;
459 const int nEmubins = 31;
460
461 // Bin convention for CosTheta, consistent with BGLRS
462 const double costheta_min = -1;
463 const double costheta_max = +1;
464 const int ncostheta = 20;
465
466 // Set up an array of CosTheta bins.
467 double costhetabinwidth = ((costheta_max-costheta_min)/ncostheta);
468 double CosThetaBinEdges[ncostheta+1];
469 for (int i=0; i<=ncostheta; i++)
470 {
471 CosThetaBinEdges[i] = costheta_min + i*costhetabinwidth;
472 }
473
474 // Set up an array of Enu bins.
475 Double_t MinLogEnu = log(Enumin);
476 Double_t MaxLogEnu = log(10*Enumin);
477 Double_t LogBinWidthEnu = ((MaxLogEnu-MinLogEnu)/nEnubinsPerDecade);
478 Double_t EnuBinEdges[nEnubins+1];
479 for (int i=0; i<=nEnubins; i++)
480 {
481 EnuBinEdges[i] = exp(MinLogEnu + i*LogBinWidthEnu);
482 }
483
484 // Set up an array of Emu bins.
485 Double_t MinLogEmu = log(Emumin);
486 Double_t MaxLogEmu = log(10*Emumin);
487 Double_t LogBinWidthEmu = ((MaxLogEmu-MinLogEmu)/nEmubinsPerDecade);
488 Double_t EmuBinEdges[nEmubins+2];
489 for (int i=0; i<=nEmubins+1; i++)
490 {
491 EmuBinEdges[i] = exp(MinLogEmu + i*LogBinWidthEmu);
492 }
493
494 // Create 3D histogram. X-axis: Emu; Y-axis: Enu; Z-axis: CosTheta.
495 TH3D *h3 = new TH3D("h3","",nEmubins+1,EmuBinEdges,nEnubins,EnuBinEdges,ncostheta,CosThetaBinEdges);
496
497 // Draw histogram.
498 //h3->Draw();
499
500 // Calculate Probability for each triplet.
501 for (int i=1; i< h3->GetNbinsX(); i++)
502 {
503 double Emu = h3->GetXaxis()->GetBinCenter(i);
504 for (int j=1; j<= h3->GetNbinsY(); j++)
505 {
506 double Enu = h3->GetBinCenter(j);
507 double Int = ProbabilityEmu(nu_code,Enu,Emu);
508 for (int k=1; k<= h3->GetNbinsZ(); k++)
509 {
510 h3->SetBinContent(i,j,k,Int); // Bin Content will is Int.
511 //h3->SetBinContent(i,j,k,1); // Bin content constant (to test)
512 }
513 }
514 }
515 return h3;
516}
double ProbabilityEmu(int nu_code, double Enu, double Emu)

References ProbabilityEmu().

Referenced by main().

◆ GenerateUpNu()

void GenerateUpNu ( GFluxI * flux_driver)

Definition at line 416 of file gUpMuFluxGen.cxx.

417{
418// Generate a Neutrino. Keep generating neutrinos until an upgoing one is generated.
419
420 while (1)
421 {
422 LOG("gevgen_upmu", pINFO) << "Pulling next neurtino from the flux driver...";
423 flux_driver->GenerateNext();
424
425 int nu_code = flux_driver->PdgCode();
426
427 const TLorentzVector & p4 = flux_driver->Momentum();
428 double costheta = -p4.Pz() / p4.Vect().Mag();
429
430 bool keep = (costheta<0) && (nu_code==kPdgNuMu || nu_code==kPdgAntiNuMu);
431 if(keep) return;
432 }
433}
#define pINFO
Definition Messenger.h:62
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
virtual bool GenerateNext(void)=0
generate the next flux neutrino (return false in err)
virtual int PdgCode(void)=0
returns the flux neutrino pdg code
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
const int kPdgAntiNuMu
Definition PDGCodes.h:31
const int kPdgNuMu
Definition PDGCodes.h:30

References genie::GFluxI::GenerateNext(), genie::kPdgAntiNuMu, genie::kPdgNuMu, LOG, genie::GFluxI::Momentum(), genie::GFluxI::PdgCode(), and pINFO.

Referenced by main().

◆ GetCommandLineArgs()

void GetCommandLineArgs ( int argc,
char ** argv )

Definition at line 643 of file gUpMuFluxGen.cxx.

644{
645// Get the command line arguments
646
647 LOG("gevgen_upmu", pINFO) << "Parsing command line arguments";
648
649 // Common run options. Set defaults and read.
651
652 // Parse run options for this app
653
654 CmdLnArgParser parser(argc,argv);
655
656 // help?
657 bool help = parser.OptionExists('h');
658 if(help) {
659 PrintSyntax();
660 exit(0);
661 }
662
663 //
664 // *** run number
665 //
666 if( parser.OptionExists('r') ) {
667 LOG("gevgen_upmu", pDEBUG) << "Reading MC run number";
668 gOptRunNu = parser.ArgAsLong('r');
669 } else {
670 LOG("gevgen_upmu", pDEBUG) << "Unspecified run number - Using default";
671 gOptRunNu = 1000;
672 } //-r
673
674 //
675 // *** exposure
676 //
677
678 // in number of events
679 bool have_required_statistics = false;
680 if( parser.OptionExists('n') ) {
681 LOG("gevgen_upmu", pDEBUG)
682 << "Reading number of events to generate";
683 gOptNev = parser.ArgAsInt('n');
684 have_required_statistics = true;
685 }//-n
686 if(!have_required_statistics) {
687 LOG("gevgen_upmu", pFATAL)
688 << "You must request exposure either in terms of number of events"
689 << "\nUse the -n option";
690 PrintSyntax();
691 gAbortingInErr = true;
692 exit(1);
693 }
694
695 //
696 // *** detector side length
697 //
698
699 if( parser.OptionExists('d') ) {
700 LOG("gevgen_upmu", pDEBUG)
701 << "Reading detector side length";
702 gOptDetectorSide = parser.ArgAsDouble('d');
703 } else {
704 LOG("gevgen_upmu", pDEBUG)
705 << "Unspecified detector side length - Using default";
707 }//-d
708
709 //
710 // *** flux files
711 //
712
713 // syntax:
714 // simulation:/path/file.data[neutrino_code],/path/file.data[neutrino_code],...
715 //
716 if( parser.OptionExists('f') ) {
717 LOG("gevgen_upmu", pDEBUG) << "Getting input flux files";
718 string flux = parser.ArgAsString('f');
719
720 // get flux simulation info (FLUKA,BGLRS) so as to instantiate the
721 // appropriate flux driver
722 string::size_type jsimend = flux.find_first_of(":",0);
723 if(jsimend==string::npos) {
724 LOG("gevgen_upmu", pFATAL)
725 << "You need to specify the flux file source";
726 PrintSyntax();
727 gAbortingInErr = true;
728 exit(1);
729 }
730 gOptFluxSim = flux.substr(0,jsimend);
731 for(string::size_type i=0; i<gOptFluxSim.size(); i++) {
732 gOptFluxSim[i] = toupper(gOptFluxSim[i]);
733 }
734 if((gOptFluxSim != "FLUKA") && (gOptFluxSim != "BGLRS")) {
735 LOG("gevgen_upmu", pFATAL)
736 << "The flux file source needs to be one of <FLUKA,BGLRS>";
737 PrintSyntax();
738 gAbortingInErr = true;
739 exit(1);
740 }
741 // now get the list of input files and the corresponding neutrino codes.
742 flux.erase(0,jsimend+1);
743 vector<string> fluxv = utils::str::Split(flux,",");
744 vector<string>::const_iterator fluxiter = fluxv.begin();
745 for( ; fluxiter != fluxv.end(); ++fluxiter) {
746 string filename_and_pdg = *fluxiter;
747 string::size_type open_bracket = filename_and_pdg.find("[");
748 string::size_type close_bracket = filename_and_pdg.find("]");
749 if (open_bracket ==string::npos ||
750 close_bracket==string::npos)
751 {
752 LOG("gevgen_upmu", pFATAL)
753 << "You made an error in specifying the flux info";
754 PrintSyntax();
755 gAbortingInErr = true;
756 exit(1);
757 }
758 string::size_type ibeg = 0;
759 string::size_type iend = open_bracket;
760 string::size_type jbeg = open_bracket+1;
761 string::size_type jend = close_bracket;
762 string flux_filename = filename_and_pdg.substr(ibeg,iend-ibeg);
763 string neutrino_pdg = filename_and_pdg.substr(jbeg,jend-jbeg);
764 gOptFluxFiles.insert(
765 map<int,string>::value_type(atoi(neutrino_pdg.c_str()), flux_filename));
766 }
767 if(gOptFluxFiles.size() == 0) {
768 LOG("gevgen_upmu", pFATAL)
769 << "You must specify at least one flux file!";
770 PrintSyntax();
771 gAbortingInErr = true;
772 exit(1);
773 }
774
775 } else {
776 LOG("gevgen_upmu", pFATAL) << "No flux info was specified! Use the -f option.";
777 PrintSyntax();
778 gAbortingInErr = true;
779 exit(1);
780 }
781
782 //
783 // *** random number seed
784 //
785 if( parser.OptionExists("seed") ) {
786 LOG("gevgen_upmu", pINFO) << "Reading random number seed";
787 gOptRanSeed = parser.ArgAsLong("seed");
788 } else {
789 LOG("gevgen_upmu", pINFO) << "Unspecified random number seed - Using default";
790 gOptRanSeed = -1;
791 }
792
793 //
794 // *** input cross-section file
795 //
796 if( parser.OptionExists("cross-sections") ) {
797 LOG("gevgen_upmu", pINFO) << "Reading cross-section file";
798 gOptInpXSecFile = parser.ArgAsString("cross-sections");
799 } else {
800 LOG("gevgen_upmu", pINFO) << "Unspecified cross-section file";
801 gOptInpXSecFile = "";
802 }
803
804
805 //
806 // print-out summary
807 //
808
809 PDGLibrary * pdglib = PDGLibrary::Instance();
810
811 ostringstream fluxinfo;
812 fluxinfo << "Using " << gOptFluxSim << " flux files: ";
813 map<int,string>::const_iterator file_iter = gOptFluxFiles.begin();
814 for( ; file_iter != gOptFluxFiles.end(); ++file_iter) {
815 int neutrino_code = file_iter->first;
816 string filename = file_iter->second;
817 TParticlePDG * p = pdglib->Find(neutrino_code);
818 if(p) {
819 string name = p->GetName();
820 fluxinfo << "(" << name << ") -> " << filename << " / ";
821 }
822 }
823
824 ostringstream expinfo;
825 if(gOptNev > 0) { expinfo << gOptNev << " events"; }
826
827 LOG("gevgen_atmo", pNOTICE)
828 << "\n\n"
829 << utils::print::PrintFramedMesg("gevgen_upmu job configuration");
830
831 LOG("gevgen_upmu", pNOTICE)
832 << "\n"
833 << "\n @@ Run number: " << gOptRunNu
834 << "\n @@ Random number seed: " << gOptRanSeed
835 << "\n @@ Using cross-section file: " << gOptInpXSecFile
836 << "\n @@ Flux"
837 << "\n\t" << fluxinfo.str()
838 << "\n @@ Exposure"
839 << "\n\t" << expinfo.str()
840 << "\n\n";
841}
#define pNOTICE
Definition Messenger.h:61
#define pFATAL
Definition Messenger.h:56
#define pDEBUG
Definition Messenger.h:63
Command line argument parser.
Singleton class to load & serve a TDatabasePDG.
Definition PDGLibrary.h:36
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
map< int, string > gOptFluxFiles
long int gOptRanSeed
Long_t gOptRunNu
string gOptFluxSim
int gOptNev
string gOptInpXSecFile
double kDefOptDetectorSide
void PrintSyntax(void)
double gOptDetectorSide
GENIE flux drivers.
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f=' *')
vector< string > Split(string input, string delim)
bool gAbortingInErr
Definition Messenger.cxx:34

References genie::CmdLnArgParser::ArgAsDouble(), genie::CmdLnArgParser::ArgAsInt(), genie::CmdLnArgParser::ArgAsLong(), genie::CmdLnArgParser::ArgAsString(), genie::PDGLibrary::Find(), genie::gAbortingInErr, gOptDetectorSide, gOptFluxFiles, gOptFluxSim, gOptInpXSecFile, gOptNev, gOptRanSeed, gOptRunNu, genie::PDGLibrary::Instance(), genie::RunOpt::Instance(), kDefOptDetectorSide, LOG, genie::CmdLnArgParser::OptionExists(), pDEBUG, pFATAL, pINFO, pNOTICE, genie::utils::print::PrintFramedMesg(), PrintSyntax(), genie::RunOpt::ReadFromCommandLine(), and genie::utils::str::Split().

Referenced by main().

◆ GetCrossSection()

double GetCrossSection ( int nu_code,
double Enu,
double Emu )

Definition at line 538 of file gUpMuFluxGen.cxx.

539{
540// Get the interaction cross section for a neutrino of energy Enu
541// to produce a muon of energy Emu.
542
543 double dxsec_dxdy = 0;
544
545 if ( Emu >= Enu ) {
546 dxsec_dxdy = 0;
547 }
548 else {
549 // get the algorithm factory & config pool
551
552 // instantiate some algorithms
553 AlgId id("genie::QPMDISPXSec","Default");
554 Algorithm * alg = algf->AdoptAlgorithm(id);
555 XSecAlgorithmI * xsecalg = dynamic_cast<XSecAlgorithmI*>(alg);
556
559
560 // Integrate over x [0,1] and y [0,1-Emu/Enu], 100 steps for each.
561 double ymax = 1 - Emu/Enu;
562
563 double dy = ymax/10;
564 double dx = 0.1;
565 double x = 0;
566 double y = 0;
567
568 for (int i = 0; i<=10; i++){
569 x = i*dx;
570 for (int j = 0; j<=10; j++){
571 y = j * dy;
572 double W = 0;
573 double Q2 = 0;
575 vp->KinePtr()->Setx(x);
576 vp->KinePtr()->Sety(y);
577 vp->KinePtr()->SetQ2(Q2);
578 vp->KinePtr()->SetW(W);
579 vn->KinePtr()->Setx(x);
580 vn->KinePtr()->Sety(y);
581 vn->KinePtr()->SetQ2(Q2);
582 vn->KinePtr()->SetW(W);
583
584 dxsec_dxdy += 0.5*(xsecalg->XSec(vp,kPSxyfE) + xsecalg->XSec(vn,kPSxyfE)) / units::cm2;
585 }
586 }
587
588 delete vp;
589 delete vn;
590 }
591
592 return dxsec_dxdy;
593}
The GENIE Algorithm Factory.
Definition AlgFactory.h:39
static AlgFactory * Instance()
Algorithm * AdoptAlgorithm(const AlgId &algid) const
Algorithm ID (algorithm name + configuration set name)
Definition AlgId.h:34
Algorithm abstract base class.
Definition Algorithm.h:54
Summary information for an interaction.
Definition Interaction.h:56
static Interaction * DISCC(int tgt, int nuc, int probe, double E=0)
Kinematics * KinePtr(void) const
Definition Interaction.h:76
void Setx(double x, bool selected=false)
void SetQ2(double Q2, bool selected=false)
void Sety(double y, bool selected=false)
void SetW(double W, bool selected=false)
Cross Section Calculation Interface.
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
static constexpr double cm2
Definition Units.h:69
double W(const Interaction *const i)
void XYtoWQ2(double Ev, double M, double &W, double &Q2, double x, double y)
double Q2(const Interaction *const i)
const int kPdgProton
Definition PDGCodes.h:81
const int kPdgNeutron
Definition PDGCodes.h:83
const int kPdgTgtFreeN
Definition PDGCodes.h:200
const int kPdgTgtFreeP
Definition PDGCodes.h:199

References genie::AlgFactory::AdoptAlgorithm(), genie::units::cm2, genie::Interaction::DISCC(), genie::AlgFactory::Instance(), genie::Interaction::KinePtr(), genie::constants::kNucleonMass, genie::kPdgNeutron, genie::kPdgProton, genie::kPdgTgtFreeN, genie::kPdgTgtFreeP, genie::kPSxyfE, genie::Kinematics::SetQ2(), genie::Kinematics::SetW(), genie::Kinematics::Setx(), genie::Kinematics::Sety(), genie::XSecAlgorithmI::XSec(), and genie::utils::kinematics::XYtoWQ2().

Referenced by main(), and ProbabilityEmu().

◆ GetDetectorVertex()

TVector3 GetDetectorVertex ( double CosTheta,
double Enu )

Definition at line 288 of file gUpMuFluxGen.cxx.

289{
290// Find the point at which the muon crosses the detector. Returns 0,0,0 if the muon misses.
291// Detector is a cube of side length l (=gOptDetectorSide).
292
293 // Get a RandomGen instance
295
296 // Generate random phi.
297 double phi = 2.*kPi* rnd->RndFlux().Rndm();
298
299 // Set distance at which incoming muon is considered
300 double rad = 0.87*gOptDetectorSide;
301
302 // Set transverse radius of a circle
303 double rad_trans = 0.87*gOptDetectorSide;
304
305 // Get necessary trig
306 double sintheta = TMath::Sqrt(1-costheta*costheta);
307 double cosphi = TMath::Cos(phi);
308 double sinphi = TMath::Sin(phi);
309
310 // Compute the muon position at distance Rad.
311 double z = rad * costheta;
312 double y = rad * sintheta * cosphi;
313 double x = rad * sintheta * sinphi;
314
315 // Displace muon randomly on a circular surface of radius rad_trans,
316 // perpendicular to a sphere radius rad at that position [x,y,z].
317 TVector3 vec(x,y,z); // vector towards selected point
318 TVector3 dvec1 = vec.Orthogonal(); // orthogonal vector
319 TVector3 dvec2 = dvec1; // second orthogonal vector
320 dvec2.Rotate(-kPi/2.0,vec); // rotate second vector by 90deg -> Cartesian coords
321
322 // Select a random point on the transverse surface, within radius rad_trans
323 double psi = 2 * kPi * rnd->RndFlux().Rndm(); // rndm angle [0,2pi]
324 double random = rnd->RndFlux().Rndm(); // rndm number [0,1]
325 dvec1.SetMag(TMath::Sqrt(random)*rad_trans*TMath::Cos(psi));
326 dvec2.SetMag(TMath::Sqrt(random)*rad_trans*TMath::Sin(psi));
327 x += dvec1.X() + dvec2.X(); // x-coord of a point the muon passes through
328 y += dvec1.Y() + dvec2.Y(); // y-coord of a point the muon passes through
329 z += dvec1.Z() + dvec2.Z(); // z-coord of a point the muon passes through
330
331 // Find out if the muon passes through any side of the detector.
332
333 // Get momentum vector
334 double pz = Enu * costheta;
335 double py = Enu * sintheta * sinphi;
336 double px = Enu * sintheta * cosphi;
337 TVector3 p3(-px,-py,-pz);
338
339 // Set up other vectors needed for line-box intersection.
340 TVector3 x3(x,y,z);
341 TVector3 temp3(x3);
342 TVector3 Hit3(0,0,0);
343
344 // Find out if the line of the muon intersects the detector.
345 bool HitFound = false;
346 double l = gOptDetectorSide;
347 TVector3 unit(0,0,0);
348 unit = p3.Unit();
349 double unitx = unit.X();
350 double unity = unit.Y();
351 double unitz = unit.Z();
352
353 // Check to see if muon intersects with z=-l/2 surface.
354 if (x3.Z() < -l/2 && unitz > 0 && !HitFound){
355 while (temp3.Z() < -l/2){
356 temp3 += unit;
357 }
358 if ( -l/2<temp3.X() && l/2>temp3.X() && -l/2<temp3.Y() && l/2>temp3.Y() ){
359 Hit3.SetXYZ(temp3.X(),temp3.Y(),-l/2);
360 HitFound = true;
361 }
362 else temp3 = x3;
363 }
364
365 // Check to see if muon intersects with x=l/2 surface.
366 if (x3.X() > l/2 && unitx < 0 && !HitFound){
367 while (temp3.X() > l/2){
368 temp3 += p3.Unit();
369 }
370 if ( -l/2<temp3.Y() && l/2>temp3.Y() && -l/2<temp3.Z() && l/2>temp3.Z() ){
371 Hit3.SetXYZ(l/2,temp3.Y(),temp3.Z());
372 HitFound = true;
373 }
374 else temp3 = x3;
375 }
376
377 // Check to see if muon intersects with y=l/2 surface.
378 if (x3.Y() > l/2 && unity < 0 && !HitFound){
379 while (temp3.Y() > l/2){
380 temp3 += p3.Unit();
381 }
382 if ( -l/2<temp3.X() && l/2>temp3.X() && -l/2<temp3.Z() && l/2>temp3.Z() ){
383 Hit3.SetXYZ(temp3.X(),l/2,temp3.Z());
384 HitFound = true;
385 }
386 else temp3 = x3;
387 }
388
389 // Check to see if muon intersects with x=-l/2 surface.
390 if (x3.X() < -l/2 && unitx > 0 && !HitFound){
391 while (temp3.X() < -l/2){
392 temp3 += p3.Unit();
393 }
394 if ( -l/2<temp3.Y() && l/2>temp3.Y() && -l/2<temp3.Z() && l/2>temp3.Z() ){
395 Hit3.SetXYZ(-l/2,temp3.Y(),temp3.Z());
396 HitFound = true;
397 }
398 else temp3 = x3;
399 }
400
401 // Check to see if muon intersects with y=-l/2 surface.
402 if (x3.Y() < -l/2 && unity > 0 && !HitFound){
403 while (temp3.Y() < -l/2){
404 temp3 += p3.Unit();
405 }
406 if ( -l/2<temp3.X() && l/2>temp3.X() && -l/2<temp3.Z() && l/2>temp3.Z() ){
407 Hit3.SetXYZ(temp3.X(),-l/2,temp3.Z());
408 HitFound = true;
409 }
410 else temp3 = x3;
411 }
412
413 return Hit3;
414}
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
TRandom3 & RndFlux(void) const
rnd number generator used by flux drivers
Definition RandomGen.h:71

References gOptDetectorSide, genie::RandomGen::Instance(), genie::constants::kPi, and genie::RandomGen::RndFlux().

Referenced by main().

◆ GetEmuPdf()

TH1D * GetEmuPdf ( double Enu,
double costheta,
const TH3D * pdf3d )

Definition at line 518 of file gUpMuFluxGen.cxx.

519{
520 // Build 1-D Emu pdf
521 int Emu_nbins = pdf3d->GetXaxis()->GetNbins();
522 const double * Emu_bins = pdf3d->GetXaxis()->GetXbins()->GetArray();
523 TH1D * pdf_Emu = new TH1D("pdf_Emu","",Emu_nbins,Emu_bins);
524 pdf_Emu->SetDirectory(0);
525
526 // Find appropriate bins for Enu and CosTheta.
527 int Enu_bin = pdf3d->GetYaxis()->FindBin(Enu);
528 int costheta_bin = pdf3d->GetZaxis()->FindBin(costheta);
529
530 // For those Enu and costheta bins, build Emu pdf
531 for (int Emu_bin = 1; Emu_bin < pdf_Emu->GetNbinsX(); Emu_bin++) {
532 pdf_Emu->SetBinContent(Emu_bin, pdf3d->GetBinContent(Emu_bin,Enu_bin,costheta_bin));
533 }
534
535 return pdf_Emu;
536}

Referenced by main().

◆ GetFlux()

GFluxI * GetFlux ( void )

Definition at line 595 of file gUpMuFluxGen.cxx.

596{
597 GFluxI * flux_driver = 0;
598
599#ifdef __GENIE_FLUX_DRIVERS_ENABLED__
600
601 // Instantiate appropriate concrete flux driver
602 GAtmoFlux * atmo_flux_driver = 0;
603 if(gOptFluxSim == "FLUKA") {
604 GFLUKAAtmoFlux * fluka_flux = new GFLUKAAtmoFlux;
605 atmo_flux_driver = dynamic_cast<GAtmoFlux *>(fluka_flux);
606 } else
607 if(gOptFluxSim == "BGLRS") {
608 GBGLRSAtmoFlux * bartol_flux = new GBGLRSAtmoFlux;
609 atmo_flux_driver = dynamic_cast<GAtmoFlux *>(bartol_flux);
610 } else {
611 LOG("gevgen_upmu", pFATAL) << "Uknonwn flux simulation: " << gOptFluxSim;
612 gAbortingInErr = true;
613 exit(1);
614 }
615
616
617 atmo_flux_driver->GenerateWeighted(true);
618
619 // Configure GAtmoFlux options (common to all concrete atmospheric flux drivers)
620 // set flux files:
621 map<int,string>::const_iterator file_iter = gOptFluxFiles.begin();
622 for( ; file_iter != gOptFluxFiles.end(); ++file_iter) {
623 int neutrino_code = file_iter->first;
624 string filename = file_iter->second;
625 atmo_flux_driver->AddFluxFile(neutrino_code, filename);
626 }
627 atmo_flux_driver->LoadFluxData();
628 // configure flux generation surface:
629 atmo_flux_driver->SetRadii(1, 1);
630 // Cast to GFluxI, the generic flux driver interface
631 flux_driver = dynamic_cast<GFluxI *>(atmo_flux_driver);
632
633#else
634 LOG("gevgen_upmu", pFATAL) << "You need to enable the GENIE flux drivers first!";
635 LOG("gevgen_upmu", pFATAL) << "Use --enable-flux-drivers at the configuration step.";
636 gAbortingInErr = true;
637 exit(1);
638#endif
639
640 return flux_driver;
641}
GENIE Interface for user-defined flux classes.
Definition GFluxI.h:29
A base class for the FLUKA, BGLRS and ATMNC atmo. nu. flux drivers. The driver depends on data files ...
Definition GAtmoFlux.h:60
void AddFluxFile(int neutrino_pdg, string filename)
virtual void GenerateWeighted(bool gen_weighted)
set whether to generate weighted or unweighted neutrinos
void SetRadii(double Rlongitudinal, double Rtransverse)
A flux driver for the Bartol Atmospheric Neutrino Flux.
A flux driver for the FLUKA 3-D Atmospheric Neutrino Flux.

References genie::flux::GAtmoFlux::AddFluxFile(), genie::gAbortingInErr, genie::flux::GAtmoFlux::GenerateWeighted(), gOptFluxFiles, gOptFluxSim, genie::flux::GAtmoFlux::LoadFluxData(), LOG, pFATAL, and genie::flux::GAtmoFlux::SetRadii().

Referenced by main().

◆ main()

int main ( int argc,
char ** argv )

Definition at line 172 of file gUpMuFluxGen.cxx.

173{
174 // Parse command line arguments
175 GetCommandLineArgs(argc,argv);
176
177 if ( ! RunOpt::Instance()->Tune() ) {
178 LOG("gmkspl", pFATAL) << " No TuneId in RunOption";
179 exit(-1);
180 }
182
183 // Init
184 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
187
188
189 // Get requested flux driver
190 GFluxI * flux_driver = GetFlux();
191
192 // Create output tree to store generated up-going muons
193 TTree * ntupmuflux = new TTree("ntupmuflux","GENIE Upgoing Muon Event Tree");
194 // Tree branches
195 int brIev = 0; // Event number
196 int brNuCode = 0; // Neutrino PDG code
197 double brEmu = 0; // Muon energy (GeV)
198 double brEnu = 0; // Neutrino energy (GeV)
199 double brCosTheta = 0; // Neutrino cos(zenith angle), muon assumed to be collinear
200 double brWghtFlxNu = 0; // Weight associated with the current flux neutrino
201 double brWghtEmuPdf = 0; // Weight associated with the Emu pdf for current Enu, costheta bins
202 double brVx = 0; // Muon x (mm) - intersection with box surrounding the detector volume
203 double brVy = 0; // Muon y (mm) - intersection with box surrounding the detector volume
204 double brVz = 0; // Muon z (mm) - intersection with box surrounding the detector volume
205 double brXSec = 0; //
206 ntupmuflux->Branch("iev", &brIev, "iev/I" );
207 ntupmuflux->Branch("nu_code", &brNuCode, "nu_code/I" );
208 ntupmuflux->Branch("Emu", &brEmu, "Emu/D" );
209 ntupmuflux->Branch("Enu", &brEnu, "Enu/D" );
210 ntupmuflux->Branch("costheta", &brCosTheta, "costheta/D" );
211 ntupmuflux->Branch("wght_fluxnu", &brWghtFlxNu, "wght_fluxnu/D" );
212 ntupmuflux->Branch("wght_emupdf", &brWghtEmuPdf, "wght_emupdf/D" );
213 ntupmuflux->Branch("vx", &brVx, "vx/D" );
214 ntupmuflux->Branch("vy", &brVy, "vy/D" );
215 ntupmuflux->Branch("vz", &brVz, "vz/D" );
216 ntupmuflux->Branch("xsec", &brXSec, "xsec/D" );
217
218 // Build 3-D pdfs describing the the probability of a muon neutrino (or anti-neutrino)
219 // of energy Enu and zenith angle costheta producing a mu- (or mu+) of energy E_mu
220 TH3D * pdf3d_numu = BuildEmuEnuCosThetaPdf (kPdgNuMu );
221 TH3D * pdf3d_numubar = BuildEmuEnuCosThetaPdf (kPdgAntiNuMu);
222
223 // Up-going muon event loop
224 for(brIev = 0; brIev < gOptNev; brIev++) {
225
226 // Loop until upgoing nu is generated.
227 GenerateUpNu(flux_driver);
228
229 // Get neutrino code, Enu, costheta and weight for the generated neutrino
230 brNuCode = flux_driver->PdgCode();
231 brEnu = flux_driver->Momentum().E();
232 brCosTheta = -1. * flux_driver->Momentum().Pz() / flux_driver->Momentum().Vect().Mag();
233 brWghtFlxNu = flux_driver->Weight();
234
235 LOG("gevgen_upmu", pNOTICE)
236 << "Generated flux neutrino: code = " << brNuCode
237 << ", Ev = " << brEnu << " GeV"
238 << ", cos(theta) = " << brCosTheta
239 << ", weight = " << brWghtFlxNu;
240
241 // Get Emu pdf my slicing the 3-D Enu,Emu,costheta pdf.
242 TH3D * pdf3d = (brNuCode == kPdgNuMu) ? pdf3d_numu : pdf3d_numubar;
243 TH1D * pdfEmu = GetEmuPdf(brEnu,brCosTheta,pdf3d);
244
245 // Get a random Emu from the pdf, and get the weight for that Emu.
246 brEmu = pdfEmu->GetRandom();
247 brWghtEmuPdf = pdfEmu->Integral("width");
248 LOG("gevgen_upmu", pNOTICE)
249 << "Selected muon has energy Emu = " << brEmu
250 << " and Emu pdg weight = " << brWghtEmuPdf;
251
252 // Randomly select the point at which the neutrino crosses the detector.
253 // assumes a simple cube detector of side length gOptDetectorSide.
254 TVector3 Vertex = GetDetectorVertex(brCosTheta,brEnu);
255 brVx = Vertex.X();
256 brVy = Vertex.Y();
257 brVz = Vertex.Z();
258
259 LOG("gevgen_upmu", pNOTICE)
260 << "Generated muon position: (" << brVx << ", " << brVy << ", " << brVz << ") m";
261
262 // Save the cross section for the interaction generated.
263 brXSec = GetCrossSection(brNuCode, brEnu, brEmu);
264
265 // save all relevant values to the ntuple
266 ntupmuflux->Fill();
267
268 // Clean-up
269 delete pdfEmu;
270 }
271
272 // Save the muon ntuple and calculate 3-D pdfs
273
274 ostringstream outfname;
275 outfname << "./genie." << gOptRunNu << ".upmu.root";
276 TFile outf(outfname.str().c_str(), "recreate");
277 ntupmuflux -> Write("ntupmu");
278 pdf3d_numu -> Write("pdf3d_numu");
279 pdf3d_numubar -> Write("pdf3d_numubar");
280 outf.Close();
281
282 // Clean-up
283 delete flux_driver;
284
285 return 0;
286}
virtual double Weight(void)=0
returns the flux neutrino weight (if any)
void BuildTune()
build tune and inform XSecSplineList
Definition RunOpt.cxx:92
void GenerateUpNu(GFluxI *flux_driver)
TVector3 GetDetectorVertex(double CosTheta, double Enu)
TH3D * BuildEmuEnuCosThetaPdf(int nu_code)
TH1D * GetEmuPdf(double Enu, double costheta, const TH3D *pdf3d)
void GetCommandLineArgs(int argc, char **argv)
GFluxI * GetFlux(void)
double GetCrossSection(int nu_code, double Enu, double Emu)
hadnt Write("hadnt")
void XSecTable(string inpfile, bool require_table)
Definition AppInit.cxx:38
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99

References BuildEmuEnuCosThetaPdf(), genie::RunOpt::BuildTune(), GenerateUpNu(), GetCommandLineArgs(), GetCrossSection(), GetDetectorVertex(), GetEmuPdf(), GetFlux(), gOptInpXSecFile, gOptNev, gOptRanSeed, gOptRunNu, genie::RunOpt::Instance(), genie::kPdgAntiNuMu, genie::kPdgNuMu, LOG, genie::utils::app_init::MesgThresholds(), genie::GFluxI::Momentum(), genie::GFluxI::PdgCode(), pFATAL, pNOTICE, genie::utils::app_init::RandGen(), genie::GFluxI::Weight(), Write(), and genie::utils::app_init::XSecTable().

◆ PrintSyntax()

void PrintSyntax ( void )

Definition at line 843 of file gUpMuFluxGen.cxx.

844{
845 LOG("gevgen_upmu", pFATAL)
846 << "\n **Syntax**"
847 << "\n gevgen_upmu [-h]"
848 << "\n [-r run#]"
849 << "\n -f simulation:flux_file[neutrino_code],..."
850 << "\n -n n_of_events,"
851 << "\n [-d detector side length (mm)]"
852 << "\n [--seed random_number_seed]"
853 << "\n --cross-sections xml_file"
854 << "\n [--message-thresholds xml_file]"
855 << "\n"
856 << " Please also read the detailed documentation at http://www.genie-mc.org"
857 << "\n";
858}

References LOG, and pFATAL.

Referenced by GetCommandLineArgs().

◆ ProbabilityEmu()

double ProbabilityEmu ( int nu_code,
double Enu,
double Emu )

Definition at line 435 of file gUpMuFluxGen.cxx.

436{
437// Calculate the probability of an incoming neutrino of energy Enu
438// generating a muon of energy Emu.
439 double dxsec_dxdy = GetCrossSection(nu_code,Enu,Emu);
440 double Int = e * constants::kNA * dxsec_dxdy / (a * (1 + Emu/e));
441 return Int;
442}
const double e
const double a

References a, e, GetCrossSection(), and genie::constants::kNA.

Referenced by BuildEmuEnuCosThetaPdf().

Variable Documentation

◆ a

const double a = 2e+6

Definition at line 164 of file gUpMuFluxGen.cxx.

Referenced by genie::Spline::Add(), genie::AGKYLowW2019::AverageChMult(), genie::utils::nuclear::BindEnergy(), genie::QPMDISStrucFuncBase::Calculate(), genie::QPMDMDISStrucFuncBase::Calculate(), genie::alvarezruso::AREikonalSolution::Cc(), genie::utils::kinematics::CohQ2Lim(), genie::RSHelicityAmplModelCC::Compute(), genie::RSHelicityAmplModelNCn::Compute(), genie::RSHelicityAmplModelNCp::Compute(), genie::mueloss::BetheBlochModel::dE_dx(), genie::Pythia8Decayer2023::Decay(), genie::alvarezruso::AlvarezRusoCOHPiPDXSec::DeltaSelfEnergyConstant(), genie::utils::nuclear::DensityGaus(), genie::Spline::Divide(), genie::mueloss::gsl::BezrukovBugaevIntegrand::DoEval(), genie::mueloss::gsl::PetrukhinShestakovIntegrand::DoEval(), genie::utils::gsl::dXSec_Log_Wrapper::DoEval(), genie::PhotonCOHPXSec::F2_Q(), genie::utils::nuclear::FmI1(), genie::utils::nuclear::FmI2(), genie::BBA03ELFormFactorsModel::Gen(), genie::LeptoHadPythia8::getRandomZ(), genie::utils::kinematics::electromagnetic::InelYLim_X(), genie::utils::kinematics::InelYLim_X(), isclose(), genie::Born::Lambda(), genie::NievesQELCCPXSec::LmunuAnumu(), genie::Spline::Multiply(), genie::operator<<(), genie::operator==(), genie::operator>>(), genie::P33PaschosLalakulichPXSec::PPiStar(), ProbabilityEmu(), genie::FGMBodekRitchie::ProbDistro(), genie::Born::PXSecCCRNC(), genie::Born::PXSecCCVNC(), genie::Born::PXSecNCVnu(), genie::Born::PXSecNCVnubar(), genie::SmithMonizUtils::Q2lim1_SM(), genie::SmithMonizUtils::Q2lim2_SM(), osetUtils::quadraticFunction(), genie::NievesQELCCPXSec::relLindhardIm(), genie::alvarezruso::integrationtools::RG202D(), genie::alvarezruso::integrationtools::RG482D(), genie::alvarezruso::integrationtools::RGN2D(), genie::BYStrucFunc::ScalingVar(), genie::DMBYStrucFunc::ScalingVar(), setA(), genie::alvarezruso::integrationtools::SG20R(), genie::alvarezruso::integrationtools::SG48R(), genie::alvarezruso::integrationtools::SGNR(), genie::utils::kinematics::TransformMatched(), genie::SmithMonizUtils::vlim1_SM(), genie::SmithMonizUtils::vlim2_SM(), genie::SmithMonizUtils::vQES_SM_lim(), and genie::BSKLNBaseRESPXSec2014::XSec().

◆ e

const double e = 500e+9

Definition at line 165 of file gUpMuFluxGen.cxx.

Referenced by genie::OutgoingDarkGenerator::AddToEventRecord(), genie::PrimaryLeptonGenerator::AddToEventRecord(), genie::utils::nuclear::BindEnergy(), genie::CollinsSpillerFragm::BuildFunction(), genie::PetersonFragm::BuildFunction(), genie::DMELEventGenerator::ComputeMaxXSec(), genie::QELEventGenerator::ComputeMaxXSec(), genie::QELEventGeneratorSM::ComputeMaxXSec(), genie::QELEventGeneratorSM::ComputeMaxXSec(), genie::SPPEventGenerator::ComputeMaxXSec(), genie::SpectralFunc::Convert2Graph(), genie::GEVGDriver::CreateXSecSumSpline(), genie::Pythia8Decayer2023::Decay(), genie::HEDISPXSec::ds_dxdy_mass(), FindhAFate(), genie::Target::ForceHitNucOnMassShell(), genie::LeptoHadPythia8::Hadronize(), genie::Pythia8Hadro2019::Hadronize(), genie::AGCharmPythia8Hadro2023::HadronizeRemnant(), genie::HEDISStrucFunc::HEDISStrucFunc(), genie::evtlib::EvtLibPXSec::Integral(), genie::NormXSec::Integral(), genie::NewQELXSec::Integrate(), genie::BardinIMDRadCorPXSec::Li2(), genie::NewQELXSec::LoadConfig(), genie::SmithMonizQELCCXSec::LoadConfig(), genie::SPPXSec::LoadConfig(), genie::utils::math::LongLorentzVector::LongLorentzVector(), main(), genie::KLVOxygenIBDPXSec::MakeAntiNuESpline(), genie::KLVOxygenIBDPXSec::MakeNuESpline(), genie::operator==(), genie::SPPEventGenerator::Vertex::operator==(), INukeNucleonCorr::OutputFiles(), ProbabilityEmu(), SaveGraphsToRootFile(), genie::DarkSectorDecayer::SetSpaceTime(), testGetTotalFluxInEnergyRange(), anonymous_namespace{Validation.cpp}::validation_plot_energy_xs_nu_nubar(), genie::GLRESWdecPythia8::Wdecay(), genie::PhotonCOHWdecPythia8::Wdecay(), genie::PhotonRESWdecPythia8::Wdecay(), genie::BardinIMDRadCorPXSec::XSec(), and genie::NormXSec::XSec().

◆ gOptDetectorSide

double gOptDetectorSide

Definition at line 159 of file gUpMuFluxGen.cxx.

Referenced by GetCommandLineArgs(), and GetDetectorVertex().

◆ gOptFluxFiles

map<int,string> gOptFluxFiles

Definition at line 157 of file gUpMuFluxGen.cxx.

◆ gOptFluxSim

string gOptFluxSim

Definition at line 156 of file gUpMuFluxGen.cxx.

◆ gOptInpXSecFile

string gOptInpXSecFile

Definition at line 161 of file gUpMuFluxGen.cxx.

◆ gOptNev

int gOptNev = -1

Definition at line 158 of file gUpMuFluxGen.cxx.

◆ gOptRanSeed

long int gOptRanSeed

Definition at line 160 of file gUpMuFluxGen.cxx.

◆ gOptRunNu

Long_t gOptRunNu

Definition at line 155 of file gUpMuFluxGen.cxx.

◆ kDefOptDetectorSide

double kDefOptDetectorSide = 1e+5

Definition at line 169 of file gUpMuFluxGen.cxx.

Referenced by GetCommandLineArgs().