LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
genie::flux::GPowerSpectrumAtmoFlux Class Reference

A driver for a power spectrum atmospheric neutrino flux. Elements extensively reused from GAtmoFlux. More...

#include "GPowerSpectrumAtmoFlux.h"

Inheritance diagram for genie::flux::GPowerSpectrumAtmoFlux:

Public Member Functions

 GPowerSpectrumAtmoFlux ()
 
 ~GPowerSpectrumAtmoFlux ()
 
const PDGCodeList & FluxParticles (void) override
 declare list of flux neutrinos that can be generated (for init. purposes) More...
 
double MaxEnergy (void) override
 declare the max flux neutrino energy that can be generated (for init. purposes) More...
 
bool GenerateNext (void) override
 generate the next flux neutrino (return false in err) More...
 
int PdgCode (void) override
 returns the flux neutrino pdg code More...
 
double Weight (void) override
 returns the flux neutrino weight (if any) More...
 
const TLorentzVector & Momentum (void) override
 returns the flux neutrino 4-momentum More...
 
const TLorentzVector & Position (void) override
 returns the flux neutrino 4-position (note: expect SI rather than physical units) More...
 
bool End (void) override
 true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples) More...
 
long int Index (void) override
 returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current entry number) More...
 
void Clear (Option_t *opt) override
 reset state variables based on opt More...
 
void GenerateWeighted (bool gen_weighted) override
 set whether to generate weighted or unweighted neutrinos More...
 
double MinEnergy (void)
 
void SetSpectralIndex (double index)
 
void SetUserCoordSystem (TRotation &rotation)
 
void SetRadii (double Rlongitudinal, double Rtransverse)
 
void SetFlavors (std::vector< int > flavors)
 
void SetMinEnergy (double Emin)
 
void SetMaxEnergy (double Emax)
 
long int NFluxNeutrinos (void) const
 
void ResetNFluxNeutrinos (void)
 
void AddFluxFile (int neutrino_pdg, string filename)
 
bool LoadFluxData (void)
 
double GetFlux (int flavour, double energy, double costh, double phi)
 
double ComputeWeight (int flavour, double energy, double costh, double phi)
 
void InitializeWeight ()
 

Private Member Functions

bool FillFluxHisto (int nu_pdg, string filename)
 
void AddAllFluxes (void)
 
TH3D * CreateNormalisedFluxHisto (TH3D *hist)
 
void ResetSelection (void)
 
void Initialize (void)
 
void CleanUp (void)
 

Private Attributes

PDGCodeList * fPdgCList
 input list of neutrino pdg-codes More...
 
double fMaxEvCut
 user-defined cut: maximum energy More...
 
double fMinEvCut
 user-defined cut: minimum energy More...
 
int fgPdgC
 current generated nu pdg-code More...
 
TLorentzVector fgP4
 current generated nu 4-momentum More...
 
TLorentzVector fgX4
 current generated nu 4-position More...
 
double fWeight
 current generated nu weight More...
 
double fSpectralIndex
 power law function More...
 
double fRl
 defining flux neutrino generation surface: longitudinal radius More...
 
double fRt
 defining flux neutrino generation surface: transverse radius More...
 
double fGlobalGenWeight
 global generation weight to apply to all events More...
 
double fAgen
 current generation area More...
 
TRotation fRotTHz2User
 coord. system rotation: THZ -> Topocentric user-defined More...
 
long int fNNeutrinos
 number of flux neutrinos thrown so far More...
 
vector< int > fFluxFlavour
 input flux file for each neutrino species More...
 
vector< string > fFluxFile
 input flux file for each neutrino species More...
 
map< int, TH3D * > fRawFluxHistoMap
 flux = f(Ev,cos8,phi) for each neutrino species More...
 
map< int, TH2D * > fRawFluxHistoMap2D
 flux = f(Ev,cos8) for each neutrino species More...
 

Detailed Description

A driver for a power spectrum atmospheric neutrino flux. Elements extensively reused from GAtmoFlux.

Author
Pierre Granger grang.nosp@m.er@a.nosp@m.pc.in.nosp@m.2p3..nosp@m.fr APC (CNRS)

December 16, 2022

Definition at line 28 of file GPowerSpectrumAtmoFlux.h.

Constructor & Destructor Documentation

GPowerSpectrumAtmoFlux::GPowerSpectrumAtmoFlux ( )

Definition at line 27 of file GPowerSpectrumAtmoFlux.cxx.

References Initialize().

28 {
29  LOG("Flux", pNOTICE)
30  << "Instantiating the GENIE Power Spectrum atmospheric neutrino flux driver";
31  this->Initialize();
32 }
GPowerSpectrumAtmoFlux::~GPowerSpectrumAtmoFlux ( )

Definition at line 36 of file GPowerSpectrumAtmoFlux.cxx.

37 {
38 
39 }

Member Function Documentation

void genie::flux::GPowerSpectrumAtmoFlux::AddAllFluxes ( void  )
private
void GPowerSpectrumAtmoFlux::AddFluxFile ( int  neutrino_pdg,
string  filename 
)

Definition at line 402 of file GPowerSpectrumAtmoFlux.cxx.

References f, fFluxFile, and fFluxFlavour.

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

403 {
404  // Check file exists
405  std::ifstream f(filename.c_str());
406  if (!f.good()) {
407  LOG("Flux", pFATAL) << "Flux file does not exist "<<filename;
408  exit(-1);
409  }
410  if ( pdg::IsNeutrino(nu_pdg) || pdg::IsAntiNeutrino(nu_pdg) ) {
411  fFluxFlavour.push_back(nu_pdg);
412  fFluxFile.push_back(filename);
413  } else {
414  LOG ("Flux", pWARN)
415  << "Input particle code: " << nu_pdg << " not a neutrino!";
416  }
417 }
TFile f
Definition: plotHisto.C:6
vector< int > fFluxFlavour
input flux file for each neutrino species
vector< string > fFluxFile
input flux file for each neutrino species
void GPowerSpectrumAtmoFlux::CleanUp ( void  )
private

Definition at line 339 of file GPowerSpectrumAtmoFlux.cxx.

References fPdgCList.

340 {
341  LOG("Flux", pNOTICE) << "Cleaning up...";
342 
343  if (fPdgCList) delete fPdgCList;
344 }
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
void GPowerSpectrumAtmoFlux::Clear ( Option_t *  opt)
override

reset state variables based on opt

Definition at line 245 of file GPowerSpectrumAtmoFlux.cxx.

References ResetNFluxNeutrinos().

246 {
247  LOG("Flux", pWARN) << "GPowerSpectrumAtmoFlux::Clear(" << opt << ") called";
248  this->ResetNFluxNeutrinos();
249 }
double GPowerSpectrumAtmoFlux::ComputeWeight ( int  flavour,
double  energy,
double  costh,
double  phi 
)

Definition at line 492 of file GPowerSpectrumAtmoFlux.cxx.

References fAgen, fGlobalGenWeight, fSpectralIndex, and GetFlux().

Referenced by GenerateNext().

493 {
494  double flux = this->GetFlux(flavour, energy, costh, phi);
495 
496  return flux*fAgen*fGlobalGenWeight*pow(energy, fSpectralIndex);
497 }
double fAgen
current generation area
double fSpectralIndex
power law function
double fGlobalGenWeight
global generation weight to apply to all events
double energy
Definition: plottest35.C:25
double GetFlux(int flavour, double energy, double costh, double phi)
TH3D* genie::flux::GPowerSpectrumAtmoFlux::CreateNormalisedFluxHisto ( TH3D *  hist)
private
bool GPowerSpectrumAtmoFlux::End ( void  )
override

true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)

Definition at line 224 of file GPowerSpectrumAtmoFlux.cxx.

225 {
226  return false;
227 }
bool GPowerSpectrumAtmoFlux::FillFluxHisto ( int  nu_pdg,
string  filename 
)
private

Definition at line 372 of file GPowerSpectrumAtmoFlux.cxx.

References f, fRawFluxHistoMap, fRawFluxHistoMap2D, and h2.

Referenced by LoadFluxData().

372  {
373  LOG("Flux", pNOTICE)
374  << "Loading fine grained flux for neutrino: " << nu_pdg
375  << " from file: " << filename;
376 
377  TH3D* histo = nullptr;
378 
379  TFile *f = new TFile(filename.c_str(), "READ");
380  f->GetObject("flux", histo);
381 
382  if(!histo) {
383  LOG("Flux", pERROR) << "Null flux histogram!";
384  f->Close();
385  return false;
386  }
387 
388  TH3D* h = static_cast<TH3D*>(histo->Clone());
389  f->Close();
390 
391  fRawFluxHistoMap.insert(std::make_pair(nu_pdg, h));
392 
393  TH2D* h2 = static_cast<TH2D*>(h->Project3D("yx"));
394 
395  fRawFluxHistoMap2D.insert(std::make_pair(nu_pdg, h2));
396 
397  return true;
398 }
map< int, TH3D * > fRawFluxHistoMap
flux = f(Ev,cos8,phi) for each neutrino species
TFile f
Definition: plotHisto.C:6
TH1F * h2
Definition: plot.C:44
map< int, TH2D * > fRawFluxHistoMap2D
flux = f(Ev,cos8) for each neutrino species
const PDGCodeList & GPowerSpectrumAtmoFlux::FluxParticles ( void  )
override

declare list of flux neutrinos that can be generated (for init. purposes)

Definition at line 44 of file GPowerSpectrumAtmoFlux.cxx.

References fPdgCList.

45 {
46  return *fPdgCList;
47 }
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
bool GPowerSpectrumAtmoFlux::GenerateNext ( void  )
override

generate the next flux neutrino (return false in err)

Definition at line 65 of file GPowerSpectrumAtmoFlux.cxx.

References ComputeWeight(), fgP4, fgPdgC, fgX4, fNNeutrinos, fPdgCList, fRl, fRotTHz2User, fRt, fSpectralIndex, fWeight, MaxEnergy(), MinEnergy(), ResetSelection(), weight, x, y, and z.

66 {
67  // Reset previously generated neutrino code / 4-p / 4-x
68  this->ResetSelection();
69 
70  // Get a RandomGen instance
71  RandomGen * rnd = RandomGen::Instance();
72 
73  // Generate (Ev, costheta, phi)
74  double Ev = 0.;
75  double costheta = 0.;
76  double phi = 0;
77  double weight = 0;
78  int nu_pdg = 0;
79 
80  // generate events according to a power law spectrum,
81  // then weight events by inverse power law
82  // (note: cannot use index alpha=1)
83  double alpha = fSpectralIndex;
84 
85  double emin = TMath::Power(this->MinEnergy(),1.0-alpha);
86  double emax = TMath::Power(this->MaxEnergy(),1.0-alpha);
87  Ev = TMath::Power(emin+(emax-emin)*rnd->RndFlux().Rndm(),1.0/(1.0-alpha));
88  costheta = -1+2*rnd->RndFlux().Rndm();
89  phi = 2.*kPi* rnd->RndFlux().Rndm();
90 
91  unsigned int nnu = fPdgCList->size();
92  unsigned int inu = rnd->RndFlux().Integer(nnu);
93  nu_pdg = (*fPdgCList)[inu];
94  weight = this->ComputeWeight(nu_pdg, Ev, costheta, phi);
95 
96  // Compute etc trigonometric numbers
97  double sintheta = TMath::Sqrt(1-costheta*costheta);
98  double cosphi = TMath::Cos(phi);
99  double sinphi = TMath::Sin(phi);
100 
101  // Set the neutrino pdg code
102  fgPdgC = nu_pdg;
103 
104  // Set the neutrino weight
105  fWeight = weight;
106 
107  // Compute the neutrino momentum
108  // The `-1' means it is directed towards the detector.
109  double pz = -1.* Ev * costheta;
110  double py = -1.* Ev * sintheta * sinphi;
111  double px = -1.* Ev * sintheta * cosphi;
112 
113  // Default vertex is at the origin
114  double z = 0.0;
115  double y = 0.0;
116  double x = 0.0;
117 
118  // Shift the neutrino position onto the flux generation surface.
119  // The position is computed at the surface of a sphere with R=fRl
120  // at the topocentric horizontal (THZ) coordinate system.
121  if( fRl>0.0 ){
122  z += fRl * costheta;
123  y += fRl * sintheta * sinphi;
124  x += fRl * sintheta * cosphi;
125  }
126 
127  // Apply user-defined rotation from THZ -> user-defined topocentric
128  // coordinate system.
129  if( !fRotTHz2User.IsIdentity() )
130  {
131  TVector3 tx3(x, y, z );
132  TVector3 tp3(px,py,pz);
133 
134  tx3 = fRotTHz2User * tx3;
135  tp3 = fRotTHz2User * tp3;
136 
137  x = tx3.X();
138  y = tx3.Y();
139  z = tx3.Z();
140  px = tp3.X();
141  py = tp3.Y();
142  pz = tp3.Z();
143  }
144 
145  // If the position is left as is, then all generated neutrinos
146  // would point towards the origin.
147  // Displace the position randomly on the surface that is
148  // perpendicular to the selected point P(xo,yo,zo) on the sphere
149  if( fRt>0.0 ){
150  TVector3 vec(x,y,z); // vector towards selected point
151  TVector3 dvec1 = vec.Orthogonal(); // orthogonal vector
152  TVector3 dvec2 = dvec1; // second orthogonal vector
153  dvec2.Rotate(-kPi/2.0,vec); // rotate second vector by 90deg,
154  // now forming a new orthogonal cartesian coordinate system
155  double psi = 2.*kPi* rnd->RndFlux().Rndm(); // rndm angle [0,2pi]
156  double random = rnd->RndFlux().Rndm(); // rndm number [0,1]
157  dvec1.SetMag(TMath::Sqrt(random)*fRt*TMath::Cos(psi));
158  dvec2.SetMag(TMath::Sqrt(random)*fRt*TMath::Sin(psi));
159  x += dvec1.X() + dvec2.X();
160  y += dvec1.Y() + dvec2.Y();
161  z += dvec1.Z() + dvec2.Z();
162  }
163 
164  // Set the neutrino momentum and position 4-vectors with values
165  // calculated at previous steps.
166  fgP4.SetPxPyPzE(px, py, pz, Ev);
167  fgX4.SetXYZT (x, y, z, 0.);
168 
169  // Increment flux neutrino counter used for sample normalization purposes.
170  fNNeutrinos++;
171 
172  // Report and exit
173  LOG("Flux", pINFO)
174  << "Generated neutrino: "
175  << "\n pdg-code: " << fgPdgC
176  << "\n p4: " << utils::print::P4AsShortString(&fgP4)
177  << "\n x4: " << utils::print::X4AsString(&fgX4);
178 
179  return true;
180 }
Float_t x
Definition: compare.C:6
long int fNNeutrinos
number of flux neutrinos thrown so far
TLorentzVector fgP4
current generated nu 4-momentum
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:276
double fSpectralIndex
power law function
double fRl
defining flux neutrino generation surface: longitudinal radius
double fRt
defining flux neutrino generation surface: transverse radius
TRotation fRotTHz2User
coord. system rotation: THZ -> Topocentric user-defined
double ComputeWeight(int flavour, double energy, double costh, double phi)
int fgPdgC
current generated nu pdg-code
double fWeight
current generated nu weight
double weight
Definition: plottest35.C:25
TLorentzVector fgX4
current generated nu 4-position
double MaxEnergy(void) override
declare the max flux neutrino energy that can be generated (for init. purposes)
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
void GPowerSpectrumAtmoFlux::GenerateWeighted ( bool  gen_weighted)
override

set whether to generate weighted or unweighted neutrinos

Definition at line 253 of file GPowerSpectrumAtmoFlux.cxx.

254 {
255 // Dummy clear method needed to conform to GFluxI interface
256 //
257  LOG("Flux", pERROR) <<
258  "The neutrinos are always generated weighted with this implementation!! GenerateWeighted method not implemented.";
259 }
double GPowerSpectrumAtmoFlux::GetFlux ( int  flavour,
double  energy,
double  costh,
double  phi 
)

Definition at line 468 of file GPowerSpectrumAtmoFlux.cxx.

References fRawFluxHistoMap, fRawFluxHistoMap2D, and h2.

Referenced by ComputeWeight().

469 {
470  TH3D* flux_hist = nullptr;
472  if(it != fRawFluxHistoMap.end())
473  {
474  flux_hist = it->second;
475  }
476 
477  LOG("Flux", pERROR) << "flux_hist: " << flux_hist;
478 
479  if(!flux_hist) return 0.0;
480 
481  if(flux_hist->GetZaxis()->GetNbins() == 1){ //no binning in phi, bilinear interpolation only so using the 2D hist
482  TH2D* h2 = fRawFluxHistoMap2D[it->first];
483  return h2->Interpolate(energy, costh);
484  }
485  else {
486  return flux_hist->Interpolate(energy, costh, phi);
487  }
488 }
intermediate_table::iterator iterator
map< int, TH3D * > fRawFluxHistoMap
flux = f(Ev,cos8,phi) for each neutrino species
double energy
Definition: plottest35.C:25
TH1F * h2
Definition: plot.C:44
map< int, TH2D * > fRawFluxHistoMap2D
flux = f(Ev,cos8) for each neutrino species
long int GPowerSpectrumAtmoFlux::Index ( void  )
override

returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current entry number)

Definition at line 231 of file GPowerSpectrumAtmoFlux.cxx.

232 {
233  return -1;
234 }
void GPowerSpectrumAtmoFlux::Initialize ( void  )
private

Definition at line 284 of file GPowerSpectrumAtmoFlux.cxx.

References fAgen, fMaxEvCut, fMinEvCut, fNNeutrinos, fPdgCList, fRl, fRotTHz2User, fRt, fSpectralIndex, InitializeWeight(), and ResetSelection().

285 {
286  LOG("Flux", pNOTICE) << "Initializing power spectrum atmospheric flux driver";
287 
288  bool allow_dup = false;
289  fPdgCList = new PDGCodeList(allow_dup);
290 
291  fSpectralIndex = 2.0;
292 
293  fMinEvCut = 0.01;
294  fMaxEvCut = 9999999999;
295 
296  // Default radii
297  fRl = 0.0;
298  fRt = 0.0;
299  fAgen = 0.0;
300 
301  // Default detector coord system: Topocentric Horizontal Coordinate system
302  fRotTHz2User.SetToIdentity();
303 
304  // Reset `current' selected flux neutrino
305  this->ResetSelection();
306 
307  // Reset number of neutrinos thrown so far
308  fNNeutrinos = 0;
309 
310  // Init the global gen weight
311  this->InitializeWeight();
312 }
long int fNNeutrinos
number of flux neutrinos thrown so far
double fAgen
current generation area
double fSpectralIndex
power law function
double fRl
defining flux neutrino generation surface: longitudinal radius
double fRt
defining flux neutrino generation surface: transverse radius
TRotation fRotTHz2User
coord. system rotation: THZ -> Topocentric user-defined
double fMinEvCut
user-defined cut: minimum energy
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
double fMaxEvCut
user-defined cut: maximum energy
void GPowerSpectrumAtmoFlux::InitializeWeight ( )

Definition at line 501 of file GPowerSpectrumAtmoFlux.cxx.

References fGlobalGenWeight, fMaxEvCut, fMinEvCut, and fSpectralIndex.

Referenced by Initialize(), SetMaxEnergy(), SetMinEnergy(), and SetSpectralIndex().

502 {
503  double IE = 1.;
504 
505  if(fSpectralIndex!=1){
506  IE = (pow(fMaxEvCut,(1.-fSpectralIndex))-pow(fMinEvCut,(1.-fSpectralIndex)))/(1.-fSpectralIndex);
507  }
508 
509  double ITheta = 4*kPi;
510 
511  fGlobalGenWeight=IE*ITheta;
512 }
double fSpectralIndex
power law function
double fGlobalGenWeight
global generation weight to apply to all events
double fMinEvCut
user-defined cut: minimum energy
double fMaxEvCut
user-defined cut: maximum energy
bool GPowerSpectrumAtmoFlux::LoadFluxData ( void  )

Definition at line 421 of file GPowerSpectrumAtmoFlux.cxx.

References fFluxFile, fFluxFlavour, FillFluxHisto(), fPdgCList, fRawFluxHistoMap, and n.

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

422 {
423  LOG("Flux", pNOTICE)
424  << "Loading atmospheric neutrino flux simulation data";
425 
426  fPdgCList->clear();
427 
428  bool loading_status = true;
429 
430  for( unsigned int n=0; n<fFluxFlavour.size(); n++ ){
431  int nu_pdg = fFluxFlavour.at(n);
432  string filename = fFluxFile.at(n);
433  string pname = PDGLibrary::Instance()->Find(nu_pdg)->GetName();
434 
435  LOG("Flux", pNOTICE) << "Loading data for: " << pname;
436 
437  bool loaded = this->FillFluxHisto(nu_pdg, filename);
438 
439  loading_status = loading_status && loaded;
440 
441  if (!loaded) {
442  LOG("Flux", pERROR)
443  << "Error loading atmospheric neutrino flux simulation data from " << filename;
444  break;
445  }
446  }
447 
448  if(loading_status) {
449  map<int,TH3D*>::iterator hist_iter = fRawFluxHistoMap.begin();
450  for ( ; hist_iter != fRawFluxHistoMap.end(); ++hist_iter) {
451  int nu_pdg = hist_iter->first;
452  fPdgCList->push_back(nu_pdg);
453  }
454 
455  LOG("Flux", pNOTICE)
456  << "Atmospheric neutrino flux simulation data loaded!";
457  return true;
458  }
459 
460  LOG("Flux", pERROR)
461  << "Error loading atmospheric neutrino flux simulation data";
462  return false;
463 }
intermediate_table::iterator iterator
map< int, TH3D * > fRawFluxHistoMap
flux = f(Ev,cos8,phi) for each neutrino species
bool FillFluxHisto(int nu_pdg, string filename)
vector< int > fFluxFlavour
input flux file for each neutrino species
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
vector< string > fFluxFile
input flux file for each neutrino species
Char_t n[5]
double GPowerSpectrumAtmoFlux::MaxEnergy ( void  )
override

declare the max flux neutrino energy that can be generated (for init. purposes)

Definition at line 51 of file GPowerSpectrumAtmoFlux.cxx.

References fMaxEvCut.

Referenced by GenerateNext().

52 {
53  return fMaxEvCut;
54 }
double fMaxEvCut
user-defined cut: maximum energy
double GPowerSpectrumAtmoFlux::MinEnergy ( void  )

Definition at line 58 of file GPowerSpectrumAtmoFlux.cxx.

References fMinEvCut.

Referenced by GenerateNext().

59 {
60  return fMinEvCut;
61 }
double fMinEvCut
user-defined cut: minimum energy
const TLorentzVector & GPowerSpectrumAtmoFlux::Momentum ( void  )
override

returns the flux neutrino 4-momentum

Definition at line 210 of file GPowerSpectrumAtmoFlux.cxx.

References fgP4.

211 {
212  return fgP4;
213 }
TLorentzVector fgP4
current generated nu 4-momentum
long int GPowerSpectrumAtmoFlux::NFluxNeutrinos ( void  ) const

Definition at line 365 of file GPowerSpectrumAtmoFlux.cxx.

References fNNeutrinos.

366 {
367  return fNNeutrinos;
368 }
long int fNNeutrinos
number of flux neutrinos thrown so far
int GPowerSpectrumAtmoFlux::PdgCode ( void  )
override

returns the flux neutrino pdg code

Definition at line 196 of file GPowerSpectrumAtmoFlux.cxx.

References fgPdgC.

197 {
198  return fgPdgC;
199 }
int fgPdgC
current generated nu pdg-code
const TLorentzVector & GPowerSpectrumAtmoFlux::Position ( void  )
override

returns the flux neutrino 4-position (note: expect SI rather than physical units)

Definition at line 217 of file GPowerSpectrumAtmoFlux.cxx.

References fgX4.

218 {
219  return fgX4;
220 }
TLorentzVector fgX4
current generated nu 4-position
void GPowerSpectrumAtmoFlux::ResetNFluxNeutrinos ( void  )

Definition at line 238 of file GPowerSpectrumAtmoFlux.cxx.

References fNNeutrinos.

Referenced by Clear().

239 {
240  fNNeutrinos = 0;
241 }
long int fNNeutrinos
number of flux neutrinos thrown so far
void GPowerSpectrumAtmoFlux::ResetSelection ( void  )
private

Definition at line 184 of file GPowerSpectrumAtmoFlux.cxx.

References fgP4, fgPdgC, fgX4, and fWeight.

Referenced by GenerateNext(), and Initialize().

185 {
186 // initializing running neutrino pdg-code, 4-position, 4-momentum
187 
188  fgPdgC = 0;
189  fgP4.SetPxPyPzE (0.,0.,0.,0.);
190  fgX4.SetXYZT (0.,0.,0.,0.);
191  fWeight = 0;
192 }
TLorentzVector fgP4
current generated nu 4-momentum
int fgPdgC
current generated nu pdg-code
double fWeight
current generated nu weight
TLorentzVector fgX4
current generated nu 4-position
void GPowerSpectrumAtmoFlux::SetFlavors ( std::vector< int >  flavors)

Definition at line 329 of file GPowerSpectrumAtmoFlux.cxx.

References fPdgCList.

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

330 {
331  fPdgCList->clear();
332  for(int flavor : flavors){
333  fPdgCList->push_back(flavor);
334  }
335 }
PDGCodeList * fPdgCList
input list of neutrino pdg-codes
void GPowerSpectrumAtmoFlux::SetMaxEnergy ( double  Emax)

Definition at line 358 of file GPowerSpectrumAtmoFlux.cxx.

References fMaxEvCut, and InitializeWeight().

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

358  {
359  fMaxEvCut = Emax;
360  this->InitializeWeight(); //Recompute the global weight
361 }
double fMaxEvCut
user-defined cut: maximum energy
void GPowerSpectrumAtmoFlux::SetMinEnergy ( double  Emin)

Definition at line 348 of file GPowerSpectrumAtmoFlux.cxx.

References fMinEvCut, and InitializeWeight().

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

349 {
350  if(Emin > 0){
351  fMinEvCut = Emin;
352  this->InitializeWeight(); //Recompute the global weight
353  }
354 }
double fMinEvCut
user-defined cut: minimum energy
void GPowerSpectrumAtmoFlux::SetRadii ( double  Rlongitudinal,
double  Rtransverse 
)

Definition at line 316 of file GPowerSpectrumAtmoFlux.cxx.

References fAgen, fRl, and fRt.

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

317 {
318  LOG ("Flux", pNOTICE) << "Setting R[longitudinal] = " << Rlongitudinal;
319  LOG ("Flux", pNOTICE) << "Setting R[transverse] = " << Rtransverse;
320 
321  fRl = Rlongitudinal;
322  fRt = Rtransverse;
323 
324  fAgen = kPi*fRt*fRt;
325 }
double fAgen
current generation area
double fRl
defining flux neutrino generation surface: longitudinal radius
double fRt
defining flux neutrino generation surface: transverse radius
void GPowerSpectrumAtmoFlux::SetSpectralIndex ( double  index)

Definition at line 263 of file GPowerSpectrumAtmoFlux.cxx.

References fSpectralIndex, and InitializeWeight().

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

264 {
265  if( index != 1.0 ){
266  fSpectralIndex = index;
267  this->InitializeWeight(); //Recompute the global weight
268  }
269  else {
270  LOG("Flux", pWARN) << "Warning: cannot use a spectral index of unity";
271  }
272  LOG("Flux", pNOTICE) << "Using Spectral Index = " << index;
273 }
double fSpectralIndex
power law function
void GPowerSpectrumAtmoFlux::SetUserCoordSystem ( TRotation &  rotation)

Definition at line 277 of file GPowerSpectrumAtmoFlux.cxx.

References fRotTHz2User.

Referenced by evgb::GENIEHelper::InitializeFluxDriver().

278 {
279  fRotTHz2User = rotation;
280 }
TRotation fRotTHz2User
coord. system rotation: THZ -> Topocentric user-defined
double GPowerSpectrumAtmoFlux::Weight ( void  )
override

returns the flux neutrino weight (if any)

Definition at line 203 of file GPowerSpectrumAtmoFlux.cxx.

References fWeight.

204 {
205  return fWeight;
206 }
double fWeight
current generated nu weight

Member Data Documentation

double genie::flux::GPowerSpectrumAtmoFlux::fAgen
private

current generation area

Definition at line 74 of file GPowerSpectrumAtmoFlux.h.

Referenced by ComputeWeight(), Initialize(), and SetRadii().

vector<string> genie::flux::GPowerSpectrumAtmoFlux::fFluxFile
private

input flux file for each neutrino species

Definition at line 78 of file GPowerSpectrumAtmoFlux.h.

Referenced by AddFluxFile(), and LoadFluxData().

vector<int> genie::flux::GPowerSpectrumAtmoFlux::fFluxFlavour
private

input flux file for each neutrino species

Definition at line 77 of file GPowerSpectrumAtmoFlux.h.

Referenced by AddFluxFile(), and LoadFluxData().

double genie::flux::GPowerSpectrumAtmoFlux::fGlobalGenWeight
private

global generation weight to apply to all events

Definition at line 73 of file GPowerSpectrumAtmoFlux.h.

Referenced by ComputeWeight(), and InitializeWeight().

TLorentzVector genie::flux::GPowerSpectrumAtmoFlux::fgP4
private

current generated nu 4-momentum

Definition at line 67 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Momentum(), and ResetSelection().

int genie::flux::GPowerSpectrumAtmoFlux::fgPdgC
private

current generated nu pdg-code

Definition at line 66 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), PdgCode(), and ResetSelection().

TLorentzVector genie::flux::GPowerSpectrumAtmoFlux::fgX4
private

current generated nu 4-position

Definition at line 68 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Position(), and ResetSelection().

double genie::flux::GPowerSpectrumAtmoFlux::fMaxEvCut
private

user-defined cut: maximum energy

Definition at line 64 of file GPowerSpectrumAtmoFlux.h.

Referenced by Initialize(), InitializeWeight(), MaxEnergy(), and SetMaxEnergy().

double genie::flux::GPowerSpectrumAtmoFlux::fMinEvCut
private

user-defined cut: minimum energy

Definition at line 65 of file GPowerSpectrumAtmoFlux.h.

Referenced by Initialize(), InitializeWeight(), MinEnergy(), and SetMinEnergy().

long int genie::flux::GPowerSpectrumAtmoFlux::fNNeutrinos
private

number of flux neutrinos thrown so far

Definition at line 76 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Initialize(), NFluxNeutrinos(), and ResetNFluxNeutrinos().

PDGCodeList* genie::flux::GPowerSpectrumAtmoFlux::fPdgCList
private

input list of neutrino pdg-codes

Definition at line 63 of file GPowerSpectrumAtmoFlux.h.

Referenced by CleanUp(), FluxParticles(), GenerateNext(), Initialize(), LoadFluxData(), and SetFlavors().

map<int, TH3D*> genie::flux::GPowerSpectrumAtmoFlux::fRawFluxHistoMap
private

flux = f(Ev,cos8,phi) for each neutrino species

Definition at line 79 of file GPowerSpectrumAtmoFlux.h.

Referenced by FillFluxHisto(), GetFlux(), and LoadFluxData().

map<int, TH2D*> genie::flux::GPowerSpectrumAtmoFlux::fRawFluxHistoMap2D
private

flux = f(Ev,cos8) for each neutrino species

Definition at line 80 of file GPowerSpectrumAtmoFlux.h.

Referenced by FillFluxHisto(), and GetFlux().

double genie::flux::GPowerSpectrumAtmoFlux::fRl
private

defining flux neutrino generation surface: longitudinal radius

Definition at line 71 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Initialize(), and SetRadii().

TRotation genie::flux::GPowerSpectrumAtmoFlux::fRotTHz2User
private

coord. system rotation: THZ -> Topocentric user-defined

Definition at line 75 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Initialize(), and SetUserCoordSystem().

double genie::flux::GPowerSpectrumAtmoFlux::fRt
private

defining flux neutrino generation surface: transverse radius

Definition at line 72 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), Initialize(), and SetRadii().

double genie::flux::GPowerSpectrumAtmoFlux::fSpectralIndex
private

power law function

Definition at line 70 of file GPowerSpectrumAtmoFlux.h.

Referenced by ComputeWeight(), GenerateNext(), Initialize(), InitializeWeight(), and SetSpectralIndex().

double genie::flux::GPowerSpectrumAtmoFlux::fWeight
private

current generated nu weight

Definition at line 69 of file GPowerSpectrumAtmoFlux.h.

Referenced by GenerateNext(), ResetSelection(), and Weight().


The documentation for this class was generated from the following files: