Public Member Functions | Public Attributes | Static Public Attributes | Friends | List of all members
genie::flux::GSimpleNtpMeta Class Reference

#include "/cvmfs/nova.opensciencegrid.org/externals/genie/v3_00_06_p01/Linux64bit+3.10-2.17-e17-debug/GENIE-Generator/src/Tools/Flux/GSimpleNtpFlux.h"

Inheritance diagram for genie::flux::GSimpleNtpMeta:

Public Member Functions

 GSimpleNtpMeta ()
 
virtual ~GSimpleNtpMeta ()
 
void Reset ()
 
void AddFlavor (Int_t nupdg)
 
void Print (const Option_t *opt="") const
 

Public Attributes

std::vector< Int_t > pdglist
 list of neutrino flavors More...
 
Double_t maxEnergy
 maximum energy More...
 
Double_t minWgt
 minimum weight More...
 
Double_t maxWgt
 maximum weight More...
 
Double_t protons
 represented number of protons-on-target More...
 
Double_t windowBase [3]
 x,y,z position of window base point More...
 
Double_t windowDir1 [3]
 dx,dy,dz of window direction 1 More...
 
Double_t windowDir2 [3]
 dx,dy,dz of window direction 2 More...
 
std::vector< std::stringauxintname
 tagname of aux ints associated w/ entry More...
 
std::vector< std::stringauxdblname
 tagname of aux doubles associated w/ entry More...
 
std::vector< std::stringinfiles
 list of input files More...
 
Int_t seed
 random seed used in generation More...
 
UInt_t metakey
 index key to tie to individual entries More...
 

Static Public Attributes

static UInt_t mxfileprint
 allow user to limit # of files to print More...
 

Friends

ostream & operator<< (ostream &stream, const GSimpleNtpMeta &info)
 

Detailed Description

GSimpleNtpMeta

A small persistable C-struct -like class that holds metadata about the the SimpleNtpFlux ntple.

Definition at line 158 of file GSimpleNtpFlux.h.

Constructor & Destructor Documentation

GSimpleNtpMeta::GSimpleNtpMeta ( )

Definition at line 906 of file GSimpleNtpFlux.cxx.

References Reset().

907  : TObject() //, nflavors(0), flavor(0)
908 {
909  Reset();
910 }
GSimpleNtpMeta::~GSimpleNtpMeta ( )
virtual

Definition at line 912 of file GSimpleNtpFlux.cxx.

References Reset().

913 {
914  Reset();
915 }

Member Function Documentation

void GSimpleNtpMeta::AddFlavor ( Int_t  nupdg)

Definition at line 943 of file GSimpleNtpFlux.cxx.

References modifyFHiCL::found, MECModelEnuComparisons::i, and pdglist.

944 {
945  bool found = false;
946  for (size_t i=0; i < pdglist.size(); ++i)
947  if ( pdglist[i] == nupdg) found = true;
948  if ( ! found ) pdglist.push_back(nupdg);
949 
950  /* // OLD fashion array
951  bool found = false;
952  for (int i=0; i < nflavors; ++i) if ( flavor[i] == nupdg ) found = true;
953  if ( ! found ) {
954  Int_t* old_list = flavor;
955  flavor = new Int_t[nflavors+1];
956  for (int i=0; i < nflavors; ++i) flavor[i] = old_list[i];
957  flavor[nflavors] = nupdg;
958  nflavors++;
959  delete [] old_list;
960  }
961  */
962 }
std::vector< Int_t > pdglist
list of neutrino flavors
void GSimpleNtpMeta::Print ( const Option_t *  opt = "") const

Definition at line 964 of file GSimpleNtpFlux.cxx.

References om::cout, allTimeWatchdog::endl, and flux.

965 {
966  std::cout << *this << std::endl;
967 }
OStream cout
Definition: OStream.cxx:6
void GSimpleNtpMeta::Reset ( void  )

Definition at line 917 of file GSimpleNtpFlux.cxx.

References auxdblname, auxintname, infiles, maxEnergy, maxWgt, metakey, minWgt, pdglist, protons, seed, windowBase, windowDir1, and windowDir2.

Referenced by GSimpleNtpMeta(), genie::flux::GSimpleNtpFlux::ProcessMeta(), and ~GSimpleNtpMeta().

918 {
919 
920  pdglist.clear();
921  maxEnergy = 0.;
922  minWgt = 0.;
923  maxWgt = 0.;
924  protons = 0.;
925  windowBase[0] = 0.;
926  windowBase[1] = 0.;
927  windowBase[2] = 0.;
928  windowDir1[0] = 0.;
929  windowDir1[1] = 0.;
930  windowDir1[2] = 0.;
931  windowDir2[0] = 0.;
932  windowDir2[1] = 0.;
933  windowDir2[2] = 0.;
934 
935  auxintname.clear();
936  auxdblname.clear();
937  infiles.clear();
938 
939  seed = 0;
940  metakey = 0;
941 }
Double_t maxEnergy
maximum energy
std::vector< std::string > auxintname
tagname of aux ints associated w/ entry
Double_t windowDir1[3]
dx,dy,dz of window direction 1
Double_t protons
represented number of protons-on-target
Int_t seed
random seed used in generation
Double_t windowDir2[3]
dx,dy,dz of window direction 2
std::vector< Int_t > pdglist
list of neutrino flavors
Double_t windowBase[3]
x,y,z position of window base point
Double_t minWgt
minimum weight
std::vector< std::string > infiles
list of input files
Double_t maxWgt
maximum weight
UInt_t metakey
index key to tie to individual entries
std::vector< std::string > auxdblname
tagname of aux doubles associated w/ entry

Friends And Related Function Documentation

ostream& operator<< ( ostream &  stream,
const GSimpleNtpMeta info 
)
friend

Definition at line 1023 of file GSimpleNtpFlux.cxx.

1025  {
1026  size_t nf = meta.pdglist.size();
1027  stream << "\nGSimpleNtpMeta " << nf << " flavors: ";
1028  for (size_t i=0; i<nf; ++i) stream << " " << meta.pdglist[i];
1029 
1030  //stream << "\nGSimpleNtpMeta " << meta.nflavors
1031  // << " flavors: ";
1032  //for (int i=0; i< meta.nflavors; ++i) stream << " " << meta.flavor[i];
1033 
1034  stream << "\n maxEnergy " << meta.maxEnergy
1035  << " min/maxWgt " << meta.minWgt << "/" << meta.maxWgt
1036  << " protons " << meta.protons
1037  << " metakey " << meta.metakey
1038  << "\n windowBase [" << meta.windowBase[0] << ","
1039  << meta.windowBase[1] << "," << meta.windowBase[2] << "]"
1040  << "\n windowDir1 [" << meta.windowDir1[0] << ","
1041  << meta.windowDir1[1] << "," << meta.windowDir1[2] << "]"
1042  << "\n windowDir2 [" << meta.windowDir2[0] << ","
1043  << meta.windowDir2[1] << "," << meta.windowDir2[2] << "]";
1044 
1045  size_t nInt = meta.auxintname.size();
1046  if ( nInt > 0 ) stream << "\n aux ints: ";
1047  for (size_t ijInt=0; ijInt < nInt; ++ijInt)
1048  stream << " " << meta.auxintname[ijInt];
1049 
1050  size_t nDbl = meta.auxdblname.size();
1051  if ( nDbl > 0 ) stream << "\n aux doubles: ";
1052  for (size_t ijDbl=0; ijDbl < nDbl; ++ijDbl)
1053  stream << " " << meta.auxdblname[ijDbl];
1054 
1055  size_t nfiles = meta.infiles.size();
1056  stream << "\n " << nfiles << " input files: ";
1057  UInt_t nprint = TMath::Min(UInt_t(nfiles),
1059  for (UInt_t ifiles=0; ifiles < nprint; ++ifiles)
1060  stream << "\n " << meta.infiles[ifiles];
1061 
1062  stream << "\n input seed: " << meta.seed;
1063 
1064  return stream;
1065  }
Class for parsing *.fcl files for their metadata information.
static UInt_t mxfileprint
allow user to limit # of files to print

Member Data Documentation

std::vector<std::string> genie::flux::GSimpleNtpMeta::auxdblname

tagname of aux doubles associated w/ entry

Definition at line 183 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::operator<<(), and Reset().

std::vector<std::string> genie::flux::GSimpleNtpMeta::auxintname

tagname of aux ints associated w/ entry

Definition at line 182 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::operator<<(), and Reset().

std::vector<std::string> genie::flux::GSimpleNtpMeta::infiles

list of input files

Definition at line 184 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::operator<<(), Reset(), and warp_meta().

Double_t genie::flux::GSimpleNtpMeta::maxEnergy

maximum energy

Definition at line 173 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), main(), genie::flux::operator<<(), Reset(), and warp_gsimple().

Double_t genie::flux::GSimpleNtpMeta::maxWgt

maximum weight

Definition at line 175 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), main(), genie::flux::operator<<(), Reset(), and warp_gsimple().

UInt_t genie::flux::GSimpleNtpMeta::metakey

index key to tie to individual entries

Definition at line 187 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::GSimpleNtpFlux::GenerateNext_weighted(), main(), genie::flux::operator<<(), Reset(), and warp_gsimple().

Double_t genie::flux::GSimpleNtpMeta::minWgt

minimum weight

Definition at line 174 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), main(), genie::flux::operator<<(), Reset(), and warp_gsimple().

UInt_t genie::flux::GSimpleNtpMeta::mxfileprint
static

allow user to limit # of files to print

Definition at line 189 of file GSimpleNtpFlux.h.

Referenced by genie::flux::operator<<().

std::vector<Int_t> genie::flux::GSimpleNtpMeta::pdglist

list of neutrino flavors

Definition at line 171 of file GSimpleNtpFlux.h.

Referenced by AddFlavor(), fill_simple(), main(), genie::flux::operator<<(), Reset(), and update_flavors().

Double_t genie::flux::GSimpleNtpMeta::protons

represented number of protons-on-target

Definition at line 176 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), main(), genie::flux::operator<<(), and Reset().

Int_t genie::flux::GSimpleNtpMeta::seed

random seed used in generation

Definition at line 186 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::operator<<(), and Reset().

Double_t genie::flux::GSimpleNtpMeta::windowBase[3]

x,y,z position of window base point

Definition at line 178 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::GSimpleNtpFlux::GetFluxWindow(), genie::flux::operator<<(), and Reset().

Double_t genie::flux::GSimpleNtpMeta::windowDir1[3]

dx,dy,dz of window direction 1

Definition at line 179 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::GSimpleNtpFlux::GetFluxWindow(), genie::flux::operator<<(), and Reset().

Double_t genie::flux::GSimpleNtpMeta::windowDir2[3]

dx,dy,dz of window direction 2

Definition at line 180 of file GSimpleNtpFlux.h.

Referenced by fill_simple(), genie::flux::GSimpleNtpFlux::GetFluxWindow(), genie::flux::operator<<(), and Reset().


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