GENIEGenerator
Loading...
Searching...
No Matches
genie::flux::GSimpleNtpMeta Class Reference

#include <GSimpleNtpFlux.h>

Inheritance diagram for genie::flux::GSimpleNtpMeta:
[legend]
Collaboration diagram for genie::flux::GSimpleNtpMeta:
[legend]

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
Double_t maxEnergy
 maximum energy
Double_t minWgt
 minimum weight
Double_t maxWgt
 maximum weight
Double_t protons
 represented number of protons-on-target
Double_t windowBase [3]
 x,y,z position of window base point
Double_t windowDir1 [3]
 dx,dy,dz of window direction 1
Double_t windowDir2 [3]
 dx,dy,dz of window direction 2
std::vector< std::string > auxintname
 tagname of aux ints associated w/ entry
std::vector< std::string > auxdblname
 tagname of aux doubles associated w/ entry
std::vector< std::string > infiles
 list of input files
Int_t seed
 random seed used in generation
UInt_t metakey
 index key to tie to individual entries

Static Public Attributes

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

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::GSimpleNtpMeta ( )

Definition at line 909 of file GSimpleNtpFlux.cxx.

910 : TObject() //, nflavors(0), flavor(0)
911{
912 Reset();
913}

References Reset().

◆ ~GSimpleNtpMeta()

GSimpleNtpMeta::~GSimpleNtpMeta ( )
virtual

Definition at line 915 of file GSimpleNtpFlux.cxx.

916{
917 Reset();
918}

References Reset().

Member Function Documentation

◆ AddFlavor()

void GSimpleNtpMeta::AddFlavor ( Int_t nupdg)

Definition at line 946 of file GSimpleNtpFlux.cxx.

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

References pdglist.

◆ Print()

void GSimpleNtpMeta::Print ( const Option_t * opt = "") const

Definition at line 967 of file GSimpleNtpFlux.cxx.

968{
969 std::cout << *this << std::endl;
970}

◆ Reset()

void GSimpleNtpMeta::Reset ( )

Definition at line 920 of file GSimpleNtpFlux.cxx.

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

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

Referenced by GSimpleNtpMeta(), and ~GSimpleNtpMeta().

◆ operator<<

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

Definition at line 1026 of file GSimpleNtpFlux.cxx.

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

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

Member Data Documentation

◆ auxdblname

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 operator<<, and Reset().

◆ auxintname

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 operator<<, and Reset().

◆ infiles

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

list of input files

Definition at line 184 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ maxEnergy

Double_t genie::flux::GSimpleNtpMeta::maxEnergy

maximum energy

Definition at line 173 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ maxWgt

Double_t genie::flux::GSimpleNtpMeta::maxWgt

maximum weight

Definition at line 175 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ metakey

UInt_t genie::flux::GSimpleNtpMeta::metakey

index key to tie to individual entries

Definition at line 187 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ minWgt

Double_t genie::flux::GSimpleNtpMeta::minWgt

minimum weight

Definition at line 174 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ mxfileprint

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 operator<<.

◆ pdglist

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

list of neutrino flavors

Definition at line 171 of file GSimpleNtpFlux.h.

Referenced by AddFlavor(), operator<<, and Reset().

◆ protons

Double_t genie::flux::GSimpleNtpMeta::protons

represented number of protons-on-target

Definition at line 176 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ seed

Int_t genie::flux::GSimpleNtpMeta::seed

random seed used in generation

Definition at line 186 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ windowBase

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 operator<<, and Reset().

◆ windowDir1

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

dx,dy,dz of window direction 1

Definition at line 179 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().

◆ windowDir2

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

dx,dy,dz of window direction 2

Definition at line 180 of file GSimpleNtpFlux.h.

Referenced by operator<<, and Reset().


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